geocat.f2py.linint2_wrapper

Module Contents

Functions

_linint1(xi, fi, xo, icycx, msg_py, shape)

_linint2(xi, yi, fi, xo, yo, icycx, msg_py, shape)

_linint2pts(xi, yi, fi, xo, yo, icycx, msg_py, shape)

linint1(→ supported_types)

Interpolates from one series to another using piecewise linear

linint2(→ supported_types)

Interpolates a regular grid to a rectilinear one using bi-linear

linint2pts(→ supported_types)

Interpolates from a rectilinear grid to an unstructured grid or

linint2_points(fi, xo, yo, icycx[, msg, meta, xi, yi])

Attributes

supported_types

geocat.f2py.linint2_wrapper.supported_types
geocat.f2py.linint2_wrapper._linint1(xi, fi, xo, icycx, msg_py, shape)
geocat.f2py.linint2_wrapper._linint2(xi, yi, fi, xo, yo, icycx, msg_py, shape)
geocat.f2py.linint2_wrapper._linint2pts(xi, yi, fi, xo, yo, icycx, msg_py, shape)
geocat.f2py.linint2_wrapper.linint1(fi: supported_types, xo: supported_types, xi: supported_types = None, icycx: numpy.number = 0, msg_py: numpy.number = None) supported_types

Interpolates from one series to another using piecewise linear interpolation across the rightmost dimension. The series may be cyclic in the X direction.

If missing values are present, then linint1 will perform the piecewise linear interpolation at all points possible, but will return missing values at coordinates which could not be used.

If any of the output coordinates xo are outside those of the input coordinates xi, the fo values at those coordinates will be set to missing (i.e. no extrapolation is performed).

Parameters
  • fi (xarray.DataArray, numpy.ndarray) –

    An array of one or more dimensions. If xi is passed in as an argument, then the size of the rightmost dimension of fi must match the rightmost dimension of xi.

    If missing values are present, then linint1 will perform the piecewise linear interpolation at all points possible, but will return missing values at coordinates which could not be used.

    Note: This variable must be supplied as a xarray.DataArray in order to copy the dimension names to the output. Otherwise, default names will be used.

  • xo (xarray.DataArray, numpy.ndarray) –

    A one-dimensional array that specifies the X coordinates of the return array. It must be strictly monotonically increasing or decreasing, but may be unequally spaced.

    If the output coordinates xo are outside those of the input coordinates xi, then the fo values at those coordinates will be set to missing (i.e. no extrapolation is performed).

  • xi (xarray.DataArray, numpy.ndarray) –

    An array that specifies the X coordinates of the fi array. Most frequently, this array is one-dimensional. It must be strictly monotonically increasing or decreasing, but can be unequally spaced. If xi is multi-dimensional, then itsimensions must be the same as fi’s dimensions. If it is one-dimensional, its length must be the same as the rightmost (fastest varying) dimension of fi.

    Note: If fi is of type xarray.DataArray and xi is left unspecified, then the rightmost coordinate dimension of fi will be used. If fi is not of type xarray.DataArray, then xi becomes a mandatory parameter. This parameter must be specified as a keyword argument.

  • icycx (bool) – An option to indicate whether the rightmost dimension of fi is cyclic. This should be set to True only if you have global data, but your longitude values don’t quite wrap all the way around the globe. For example, if your longitude values go from, say, -179.75 to 179.75, or 0.5 to 359.5, then you would set this to True.

  • msg_py (numpy.number) – A numpy scalar value that represent a missing value in fi. This argument allows a user to use a missing value scheme other than NaN or masked arrays, similar to what NCL allows.

Returns

fo – The interpolated series. The returned value will have the same dimensions as fi, except for the rightmost dimension which will have the same dimension size as the length of xo. The return type will be double if fi is double, and float otherwise.

Return type

xarray.DataArray, numpy.ndarray

Examples

Example 1: Using linint1 with xarray.DataArray input

import numpy as np
import xarray as xr
import geocat.comp
fi_np = np.random.rand(80)  # random 80-element array
# xi does not have to be equally spaced, but it is
# in this example
xi = np.arange(80)
# create target coordinate array, in this case use the same
# min/max values as xi, but with different spacing
xo = np.linspace(xi.min(), xi.max(), 100)
# create :class:`xarray.DataArray` and chunk it using the
# full shape of the original array.
# note that xi is attached as a coordinate array
fi = xr.DataArray(fi_np,
                  dims=['x'],
                  coords={'x': xi}
                 ).chunk(fi_np.shape)
fo = geocat.comp.linint1(fi, xo, icycx=0)
geocat.f2py.linint2_wrapper.linint2(fi: supported_types, xo: supported_types, yo: supported_types, xi: supported_types = None, yi: supported_types = None, icycx: bool = 0, msg_py: numpy.number = None) supported_types

Interpolates a regular grid to a rectilinear one using bi-linear interpolation. The input grid may be cyclic in the x-direction. The interpolation is first performed in the x-direction, and then in the y-direction.

Parameters
  • fi (xarray.DataArray, numpy.ndarray) –

    An array of two or more dimensions. If xi is passed in as an argument, then the size of the rightmost dimension of fi must match the rightmost dimension of xi. Similarly, if yi is passed in as an argument, then the size of the second-rightmost dimension of fi must match the rightmost dimension of yi.

    If missing values are present, then linint2 will perform the bilinear interpolation at all points possible, but will return missing values at coordinates which could not be used.

    Note: This variable must be supplied as a xarray.DataArray in order to copy the dimension names to the output. Otherwise, default names will be used.

  • xo (xarray.DataArray, numpy.ndarray) –

    A one-dimensional array that specifies the X-coordinates of the return array. It must be strictly monotonically increasing, but may be unequally spaced. For geo-referenced data, xo is generally the longitude array.

    If the output coordinates xo are outside those of the input coordinates xi, then the fo values at those coordinates will be set to missing (i.e. no extrapolation is performed).

  • yo (xarray.DataArray, numpy.ndarray) –

    A one-dimensional array that specifies the Y coordinates of the return array. It must be strictly monotonically increasing, but may be unequally spaced. For geo-referenced data, yo is typically the latitude array.

    If the output coordinates yo are outside those of the input coordinates yi, then the fo values at those coordinates will be set to missing (i.e. no extrapolation is performed).

  • xi (xarray.DataArray, numpy.ndarray) –

    An array that specifies the X-coordinates of the fi array. Most frequently, this is a 1D strictly monotonically increasing array that may be unequally spaced. In some cases, xi can be a multi-dimensional array (see next paragraph). The rightmost dimension (call it nxi) must have at least two elements, and is the last (fastest varying) dimension of fi.

    If xi is a multi-dimensional array, then each nxi subsection of xi must be strictly monotonically increasing, but may be unequally spaced. All but its rightmost dimension must be the same size as all but fi’s rightmost two dimensions. For geo-referenced data, xi is generally the longitude array.

    Note: If fi is of type xarray.DataArray and xi is left unspecified, then the rightmost coordinate dimension of fi will be used. If fi is not of type xarray.DataArray, then xi becomes a mandatory parameter. This parameter must be specified as a keyword argument.

  • yi (xarray.DataArray, numpy.ndarray) –

    An array that specifies the Y-coordinates of the fi array. Most frequently, this is a 1D strictly monotonically increasing array that may be unequally spaced. In some cases, yi can be a multi-dimensional array (see next paragraph). The rightmost dimension (call it nyi) must have at least two elements, and is the second-to-last dimension of fi.

    If yi is a multi-dimensional array, then each nyi subsection of yi must be strictly monotonically increasing, but may be unequally spaced. All but its rightmost dimension must be the same size as all but fi’s rightmost two dimensions. For geo-referenced data, yi is generally the latitude array.

    Note: If fi is of type xarray.DataArray and xi is left unspecified, then the second-to-rightmost coordinate dimension of fi will be used. If fi is not of type xarray.DataArray, then xi becomes a mandatory parameter. This parameter must be specified as a keyword argument.

  • icycx (bool) – An option to indicate whether the rightmost dimension of fi is cyclic. This should be set to True only if you have global data, but your longitude values don’t quite wrap all the way around the globe. For example, if your longitude values go from, say, -179.75 to 179.75, or 0.5 to 359.5, then you would set this to True.

  • msg_py (numpy.number) – A numpy scalar value that represent a missing value in fi. This argument allows a user to use a missing value scheme other than NaN or masked arrays, similar to what NCL allows.

Returns

fo – The interpolated grid. If the meta parameter is True, then the result will include named dimensions matching the input array. The returned value will have the same dimensions as fi, except for the rightmost two dimensions which will have the same dimension sizes as the lengths of yo and xo. The return type will be double if fi is double, and float otherwise.

Return type

xarray.DataArray, numpy.ndarray

Examples

Example 1: Using linint2 with xarray.DataArray input

import numpy as np
import xarray as xr
import geocat.comp
fi_np = np.random.rand(30, 80)  # random 30x80 array
# xi and yi do not have to be equally spaced, but they are
# in this example
xi = np.arange(80)
yi = np.arange(30)
# create target coordinate arrays, in this case use the same
# min/max values as xi and yi, but with different spacing
xo = np.linspace(xi.min(), xi.max(), 100)
yo = np.linspace(yi.min(), yi.max(), 50)
# create :class:`xarray.DataArray` and chunk it using the
# full shape of the original array.
# note that xi and yi are attached as coordinate arrays
fi = xr.DataArray(fi_np,
                  dims=['lat', 'lon'],
                  coords={'lat': yi, 'lon': xi}
                 ).chunk(fi_np.shape)
fo = geocat.comp.linint2(fi, xo, yo, icycx=0)
geocat.f2py.linint2_wrapper.linint2pts(fi: supported_types, xo: supported_types, yo: supported_types, icycx: bool = False, msg_py: numpy.number = None, xi: supported_types = None, yi: supported_types = None) supported_types

Interpolates from a rectilinear grid to an unstructured grid or locations using bilinear interpolation.

The linint2pts function uses bilinear interpolation to interpolate from a rectilinear grid to an unstructured grid.

If missing values are present, then linint2pts will perform the piecewise linear interpolation at all points possible, but will return missing values at coordinates which could not be used.

If one or more of the four closest grid points to a particular (xo, yo) coordinate pair are missing, then the return value for this coordinate pair will be missing.

If the user inadvertently specifies output coordinates (xo, yo) that are outside those of the input coordinates (xi, yi), the output value at this coordinate pair will be set to missing as no extrapolation is performed.

linint2pts is different from linint2 in that xo and yo are coordinate pairs, and need not be monotonically increasing. It is also different in the dimensioning of the return array. This function could be used if the user wanted to interpolate gridded data to, say, the location of rawinsonde sites or buoy/xbt locations.

Warning: if xi contains longitudes, then the xo values must be in the same range. In addition, if the xi` values span 0 to 360, then the xo values must also be specified in this range (i.e. -180 to 180 will not work).

Parameters
  • fi (xarray.DataArray, numpy.ndarray) – An array of two or more dimensions. The two rightmost dimensions (nyi x nxi) are the dimensions to be used in the interpolation. If user-defined missing values are present (other than NaNs), the value of msg_py must be set appropriately.

  • xo (xarray.DataArray, numpy.ndarray) – A one-dimensional array that specifies the X-coordinates of the unstructured grid.

  • yo (xarray.DataArray, numpy.ndarray) – A one-dimensional array that specifies the Y-coordinates of the unstructured grid. It must be the same length as xo.

  • icycx (bool) – An option to indicate whether the rightmost dimension of fi is cyclic. Default valus is 0. This should be set to True only if you have global data, but your longitude values don’t quite wrap all the way around the globe. For example, if your longitude values go from, say, -179.75 to 179.75, or 0.5 to 359.5, then you would set this to True.

  • msg_py (numpy.number) – A numpy scalar value that represent a missing value in fi. This argument allows a user to use a missing value scheme other than NaN or masked arrays, similar to what NCL allows.

  • xi (xarray.DataArray, numpy.ndarray) – A strictly monotonically increasing array that specifies the X-coordinates of the fi array. xi might be defined as the coordinates of fi when fi is of type xarray.DataArray; in this case xi may not be explicitly given as a function argument.

  • yi (xarray.DataArray, numpy.ndarray) – A strictly monotonically increasing array that specifies the Y [latitude] coordinates of the fi array. yi might be defined as the coordinates of fi when fi is of type xarray.DataArray; in this case yi may not be explicitly given as a function argument.

Returns

fo – The returned value will have the same dimensions as fi, except for the rightmost dimension which will have the same dimension size as the length of yo and xo.

Return type

xarray.DataArray, numpy.ndarray

Examples

Example 1: Using linint2pts with xarray.DataArray input

import numpy as np
import xarray as xr
import geocat.comp
fi_np = np.random.rand(30, 80)  # random 30x80 array
# xi and yi do not have to be equally spaced, but they are
# in this example
xi = np.arange(80)
yi = np.arange(30)
# create target coordinate arrays, in this case use the same
# min/max values as xi and yi, but with different spacing
xo = np.linspace(xi.min(), xi.max(), 100)
yo = np.linspace(yi.min(), yi.max(), 50)
# create :class:`xarray.DataArray` and chunk it using the
# full shape of the original array.
# note that xi and yi are attached as coordinate arrays
fi = xr.DataArray(fi_np,
                  dims=['lat', 'lon'],
                  coords={'lat': yi, 'lon': xi}
                 ).chunk(fi_np.shape)
fo = geocat.comp.linint2pts(fi, xo, yo, 0)
geocat.f2py.linint2_wrapper.linint2_points(fi, xo, yo, icycx, msg=None, meta=False, xi=None, yi=None)