cf.Field.regrids¶
-
Field.
regrids
(dst, method=None, src_cyclic=None, dst_cyclic=None, use_src_mask=True, use_dst_mask=False, fracfield=False, src_axes=None, dst_axes=None, axis_order=None, ignore_degenerate=True, return_operator=False, check_coordinates=False, min_weight=None, weights_file=None, src_z=None, dst_z=None, z=None, ln_z=None, verbose=None, return_esmpy_regrid_operator=False, inplace=False, i=False, _compute_field_mass=None)[source]¶ Regrid the field to a new latitude and longitude grid.
Regridding is the process of interpolating the field data values while preserving the qualities of the original data, and the metadata of the unaffected axes. The metadata for the regridded axes are taken from the dst parameter.
The 2-d or 3-d regridding takes place on a sphere, with the grid being defined by latitude and longitude spherical polar coordinates, and any available vertical coordinates. In the 3-d case, the regridding may be done assuming linear or log linear weights in the vertical.
Latitude and longitude coordinates
The source and destination grids of the regridding must both be defined by latitude and longitude coordinates, which may be 1-d dimension coordinates or 2-d auxiliary coordinates. These are automatically detected from the field being regridded and the specification of the destination grid given by the dst parameter.
When a grid is defined by 2-d latitude and longitude coordinates, it is necessary for their X and Y dimensions to be defined. This is either automatically inferred from the existence of 1-d dimension coordinates, or else must be specified with src_axes or dst_axes parameters.
Curvilinear Grids
Grids in projection coordinate systems can be regridded as long as two dimensional latitude and longitude coordinates are present.
Tripolar Grids
Connections across the bipole fold are not currently supported, but are not necessary in some cases, for example if the points on either side are together without a gap (as is the case for NEMO model outputs).
UGRID meshes
Data defined on UGRID face or node cells may be regridded to any other latitude-longitude grid, including other UGRID meshes and DSG feature types.
DSG feature types
Data on any latitude-longitude grid (including tripolar and UGRID meshes) may be regridded to any DSG feature type.
Cyclicity of the X axis
The cyclicity of the X (longitude) axes of the source and destination grids (i.e. whether or not the first and last cells of the axis are adjacent) are taken into account. By default, the cyclicity is inferred from the grids’ defining coordinates, but may be also be provided with the src_cyclic and dst_cyclic parameters.
Masked cells
By default, the data mask of the source data is taken into account during the regridding process, but the destination grid mask is not. This behaviour may be changed with the use_src_mask and use_dst_mask parameters.
In general the source data may be arbitrarily masked, meaning that the mask for the regridding axes may vary along the non-regridding axes. The exceptions to this are for second-order conservative, patch recovery regridding, and nearest source to destination methods, for which the mask of the regridding axes must be the same across all non-regridding axes. In these special cases an exception will be raised if the source data mask does not meet this requirement.
Implementation
The interpolation is carried out using regridding weights calculated by the
esmpy
package, a Python interface to the Earth System Modeling Framework regridding utility (ESMF, https://earthsystemmodeling.org/regrid). Outside ofesmpy
, these weights are then modified for masked cells (if required) and the regridded data are created as the dot product of the weights with the source data. (Note that whilst theesmpy
package is also able to create the regridded data from its weights, this feature can’t be integrated with thedask
framework that underpins the field’s data.)Logging
Whether
esmpy
logging is enabled or not is determined bycf.regrid_logging
. If it is enabled then logging takes place after every call. By default logging is disabled.New in version 1.0.4.
See also
- Parameters
- dst:
Field
,Domain
,RegridOperator
or sequence ofCoordinate
The definition of the destination grid on which to regrid the field’s data. One of:
Field
: The grid is defined by the latitude and longitude coordinates of the field construct’s domain.Domain
: The grid is defined by the latitude and longitude coordinates of the domain construct.Sequence of
Coordinate
: The grid is defined by two 1-d dimension coordinate constructs, or two 2-d auxiliary coordinate constructs, that define the spherical latitude and longitude coordinates (in any order) of the destination grid.In the 2-d case, both coordinate constructs must have their axes in the same order, which must be specified by the dst_axes parameter.
RegridOperator
: The grid is defined by a regrid operator that has been returned by a previous call with the return_operator parameter set to True.Unlike the other options, for which the regrid weights need to be calculated, the regrid operator already contains the weights. Therefore, for cases where multiple fields with the same source grids need to be regridded to the same destination grid, using a regrid operator can give performance improvements by avoiding having to calculate the weights for each source field. Note that for the other types of dst parameter, the calculation of the regrid weights is not a lazy operation.
Note
When dst is a
RegridOperator
, the source grid of the regrid operator is immediately checked for compatibility with the grid of the source field. By default only the computationally cheap tests are performed (checking that the coordinate system, cyclicity, grid shape, regridding dimensionality, mesh location, and feature type are the same), with the grid coordinates not being checked. The coordinates check will be carried out, however, if the check_coordinates parameter is True.
- method:
str
orNone
, optional Specify the regridding interpolation method. This parameter must be set unless dst is a
RegridOperator
, when the method is ignored.The method parameter may be one of the following:
'linear'
: Bilinear interpolation.'bilinear'
: Deprecated alias for'linear'
.'conservative_1st'
: First order conservative interpolation. Preserves the integral of the source field across the regridding. Weight calculation is based on the ratio of source cell area overlapped with the corresponding destination cell area.'conservative'
: Alias for'conservative_1st'
'conservative_2nd'
: Second-order conservative interpolation. Preserves the integral of the source field across the regridding. Weight calculation is based on the ratio of source cell area overlapped with the corresponding destination cell area. The second-order conservative calculation also includes the gradient across the source cell, so in general it gives a smoother, more accurate representation of the source field. This is particularly true when going from a coarse to finer grid.'patch'
Patch recovery interpolation. Patch rendezvous method of taking the least squares fit of the surrounding surface patches. This is a higher order method that may produce interpolation weights that may be slightly less than 0 or slightly greater than 1. This method typically results in better approximations to values and derivatives when compared to bilinear interpolation.'nearest_stod'
: Nearest neighbour source to destination interpolation for which each destination point is mapped to the closest source point. A source point can be mapped to multiple destination points. Useful for regridding categorical data.'nearest_dtos'
: Nearest neighbour destination to source interpolation for which each source point is mapped to the closest destination point. A destination point can be mapped to multiple source points. Some destination points may not be mapped. Useful for regridding of categorical data.None
: This is the default and can only be used when dst is aRegridOperator
.
- src_cyclic:
None
orbool
, optional Specifies whether or not the source grid longitude axis is cyclic (i.e. the first and last cells of the axis are adjacent). If
None
(the default) then the cyclicity will be inferred from the source grid coordinates, defaulting toFalse
if it can not be determined.- dst_cyclic:
None
orbool
, optional Specifies whether or not the destination grid longitude axis is cyclic (i.e. the first and last cells of the axis are adjacent). If
None
(the default) then the cyclicity will be inferred from the destination grid coordinates, defaulting toFalse
if it can not be determined.Ignored if dst is a
RegridOperator
.- use_src_mask:
bool
, optional By default the mask of the source field is taken into account during the regridding process. The only possible exception to this is when the nearest source to destination regridding method (
'nearest_stod'
) is being used. In this case, if use_src_mask is False then each destination point is mapped to the closest source point, whether or not it is masked (see the method parameter for details).Ignored if dst is a
RegridOperator
.- use_dst_mask:
bool
, optional If dst is a
Field
and use_dst_mask is False (the default) then the mask of data on the destination grid is not taken into account when performing regridding. If use_dst_mask is True then any masked cells in the dst field construct are transferred to the result. If dst has more dimensions than are being regridded, then the mask of the destination grid is taken as the subspace defined by index0
of all the non-regridding dimensions.Ignored if dst is not a
Field
.- src_axes:
dict
, optional When the source grid’s X and Y dimensions can not be inferred from the existence of 1-d dimension coordinates, then they must be identified with the src_axes dictionary, with keys
'X'
and'Y'
.The dictionary values identify a unique domain axis by passing the given axis description to a call of the field construct’s
domain_axis
method. For example, for a value of'ncdim%x'
, the domain axis construct returned byf.domain_axis('ncdim%x')
is selected.- Parameter example:
{'X': 'ncdim%x', 'Y': 'ncdim%y'}
- Parameter example:
{'X': 1, 'Y': 0}
- dst_axes:
dict
, optional When the destination grid’s X and Y dimensions can not be inferred from the existence of 1-d dimension coordinates, then they must be identified with the dst_axes dictionary, with keys
'X'
and'Y'
.If dst is a
Field
orDomain
, then the dictionary values identify a unique domain axis by passing the given axis description to a call of the destination field or domain construct’sdomain_axis
method. For example, for a value of'ncdim%x'
, the domain axis construct returned byf.domain_axis('ncdim%x')
is selected.If dst is a sequence of
Coordinate
, then the dictionary values identify a unique domain axis by its position in the 2-d coordinates’ data arrays, i.e. the dictionary values must be0
and1
:Ignored if dst is a
RegridOperator
.- Parameter example:
{'X': 'ncdim%x', 'Y': 'ncdim%y'}
- Parameter example:
{'X': 1, 'Y': 0}
- ignore_degenerate:
bool
, optional For conservative regridding methods, if True (the default) then degenerate cells (those for which enough vertices collapse to leave a cell as either a line or a point) are skipped, not producing a result. Otherwise an error will be produced if degenerate cells are found, that will be present in the
esmpy
log files.For all other regridding methods, degenerate cells are always skipped, regardless of the value of ignore_degenerate.
Ignored if dst is a
RegridOperator
.- return_operator:
bool
, optional If True then do not perform the regridding, rather return the
RegridOperator
instance that defines the regridding operation, and which can be used in subsequent calls. See the dst parameter for details.New in version 3.10.0.
- check_coordinates:
bool
, optional If True and dst is a `RegridOperator`then the source grid coordinates defined by the operator are checked for compatibility against those of the source field. By default this check is not carried out. See the dst parameter for details.
Ignored unless dst is a
RegridOperator
.New in version 3.14.0.
- min_weight: float, optional
A very small non-negative number. By default min_weight is
2.5 * np.finfo("float64").eps
, i.e.5.551115123125783e-16
. It is used during linear and first-order conservative regridding when adjusting the weights matrix to account for the data mask. It is ignored for all other regrid methods, or if data being regridded has no missing values.In some cases (described below) for which weights might only be non-zero as a result of rounding errors, the min_weight parameter controls whether or a not cell in the regridded field is masked.
The default value has been chosen empirically as the smallest value that produces the same masks as
esmpy
for the use cases defined in the cf test suite.Define
w_ji
as the multiplicative weight that defines how much ofVs_i
(the value in source grid celli
) contributes toVd_j
(the value in destination grid cellj
).Linear regridding
Destination grid cell
j
will only be masked if a) it is masked in the destination grid definition; or b)w_ji >= min_weight
for those masked source grid cellsi
for whichw_ji > 0
.Conservative first-order regridding
Destination grid cell
j
will only be masked if a) it is masked in the destination grid definition; or b) the sum ofw_ji
for all non-masked source grid cellsi
is strictly less than min_weight.New in version 3.14.0.
- weights_file:
str
orNone
, optional Provide a netCDF file that contains, or will contain, the regridding weights. If
None
(the default) then the weights are computed in memory for regridding between the source and destination grids, and no file is created.If set to a file path that does not exist then the weights will be computed and also written to that file.
If set to a file path that already exists then the weights will be read from this file, instead of being computed.
Note
No checks are performed on an existing file to ensure that the weights are appropriate for the source field and the values of the keyword parameters. Inappropriate weights will produce incorrect results.
However, when regridding using weights from a file, ensuring that the source field has the same shape over the regridding axes, and the parameter settings are the same as those used when the weights file was created, will ensure correct results.
A netCDF regridding weights file created directly by ESMF has the same structure and variable names (
S
,row
, andcol
for the weights, destination/row indices, and source/col indices respectively), so may be provided as a weights_file, noting that no checks will be applied to it.Performance
The computation of the weights can be much more costly than the regridding itself, in which case reading pre-calculated weights can improve performance.
Ignored if dst is a
RegridOperator
.Ignored if dst is a
RegridOperator
.New in version 3.15.2.
- src_z: optional
If
None
, the default, then the regridding is 2-d in the latitude-longitude plane.If not
None
then 3-d spherical regridding is enabled by identifying the source grid vertical coordinates from which to derive the vertical component of the regridding weights. The vertical coordinate construct may be 1-d or 3-d and is defined by the unique construct returned byf.coordinate(src_z)
Ignored if dst is a
RegridOperator
.New in version 3.16.2.
- dst_z: optional
If
None
, the default, then the regridding is 2-d in the latitude-longitude plane.If not
None
then 3-d spherical regridding is enabled by identifying the destination grid vertical coordinates from which to derive the vertical component of the regridding weights. The vertical coordinate construct may be 1-d or 3-d.Ignored if dst is a
RegridOperator
.New in version 3.16.2.
- z: optional
The z parameter is a convenience that may be used to replace both src_z and dst_z when they would contain identical values. If not
None
then 3-d spherical regridding is enabled. See src_z and dst_z for details.Ignored if dst is a
RegridOperator
.- Example:
z='Z'
is equivalent tosrc_z='Z', dst_z='Z'
.
New in version 3.16.2.
- ln_z:
bool
orNone
, optional If True when z, src_z or dst_z are also set, calculate the vertical component of the regridding weights using the natural logarithm of the vertical coordinate values. This option should be used if the quantity being regridded varies approximately linearly with logarithm of the vertical coordinates. If False, then the weights are calculated using unaltered vertical values. If
None
, the default, then an exception is raised if any of z, src_z or dst_z have also been set.Ignored if dst is a
RegridOperator
.New in version 3.16.2.
- verbose:
int
orstr
orNone
, optional If an integer from
-1
to3
, or an equivalent string equal ignoring case to one of:'DISABLE'
(0
)'WARNING'
(1
)'INFO'
(2
)'DETAIL'
(3
)'DEBUG'
(-1
)
set for the duration of the method call only as the minimum cut-off for the verboseness level of displayed output (log) messages, regardless of the globally-configured
cf.log_level
. Note that increasing numerical value corresponds to increasing verbosity, with the exception of-1
as a special case of maximal and extreme verbosity.Otherwise, if
None
(the default value), output messages will be shown according to the value of thecf.log_level
setting.Overall, the higher a non-negative integer or equivalent string that is set (up to a maximum of
3
/'DETAIL'
) for increasing verbosity, the more description that is printed to convey information about the operation.New in version 3.16.0.
- inplace:
bool
, optional If True then do the operation in-place and return
None
.- return_esmpy_regrid_operator:
bool
, optional If True then do not perform the regridding, rather return the
esmpy.Regrid
instance that defines the regridding operation.New in version 3.16.2.
- axis_order: sequence, optional
Deprecated at version 3.14.0.
- fracfield:
bool
, optional Deprecated at version 3.14.0.
- _compute_field_mass:
dict
, optional Deprecated at version 3.14.0.
- i: deprecated at version 3.0.0
Use the inplace parameter instead.
- dst:
- Returns
Field
orNone
orRegridOperator
The regridded field construct; or
None
if the operation was in-place; or the regridding operator if return_operator is True; or theesmpy.Regrid
operator object if return_esmpy_regrid_operator is True.
Examples
>>> src, dst = cf.example_fields(1, 0) >>> print(src) Field: air_temperature (ncvar%ta) --------------------------------- Data : air_temperature(atmosphere_hybrid_height_coordinate(1), grid_latitude(10), grid_longitude(9)) K Cell methods : grid_latitude(10): grid_longitude(9): mean where land (interval: 0.1 degrees) time(1): maximum Field ancils : air_temperature standard_error(grid_latitude(10), grid_longitude(9)) = [[0.76, ..., 0.32]] K Dimension coords: atmosphere_hybrid_height_coordinate(1) = [1.5] : grid_latitude(10) = [2.2, ..., -1.76] degrees : grid_longitude(9) = [-4.7, ..., -1.18] degrees : time(1) = [2019-01-01 00:00:00] Auxiliary coords: latitude(grid_latitude(10), grid_longitude(9)) = [[53.941, ..., 50.225]] degrees_N : longitude(grid_longitude(9), grid_latitude(10)) = [[2.004, ..., 8.156]] degrees_E : long_name=Grid latitude name(grid_latitude(10)) = [--, ..., kappa] Cell measures : measure:area(grid_longitude(9), grid_latitude(10)) = [[2391.9657, ..., 2392.6009]] km2 Coord references: grid_mapping_name:rotated_latitude_longitude : standard_name:atmosphere_hybrid_height_coordinate Domain ancils : ncvar%a(atmosphere_hybrid_height_coordinate(1)) = [10.0] m : ncvar%b(atmosphere_hybrid_height_coordinate(1)) = [20.0] : surface_altitude(grid_latitude(10), grid_longitude(9)) = [[0.0, ..., 270.0]] m >>> print(dst) Field: specific_humidity (ncvar%q) ---------------------------------- Data : specific_humidity(latitude(5), longitude(8)) 1 Cell methods : area: mean Dimension coords: latitude(5) = [-75.0, ..., 75.0] degrees_north : longitude(8) = [22.5, ..., 337.5] degrees_east : time(1) = [2019-01-01 00:00:00] >>> x = src.regrids(dst, method='linear') >>> print(x) Field: air_temperature (ncvar%ta) --------------------------------- Data : air_temperature(atmosphere_hybrid_height_coordinate(1), latitude(5), longitude(8)) K Cell methods : latitude(5): longitude(8): mean where land (interval: 0.1 degrees) time(1): maximum Dimension coords: atmosphere_hybrid_height_coordinate(1) = [1.5] : latitude(5) = [-75.0, ..., 75.0] degrees_north : longitude(8) = [22.5, ..., 337.5] degrees_east : time(1) = [2019-01-01 00:00:00] Coord references: standard_name:atmosphere_hybrid_height_coordinate Domain ancils : ncvar%a(atmosphere_hybrid_height_coordinate(1)) = [10.0] m : ncvar%b(atmosphere_hybrid_height_coordinate(1)) = [20.0] : surface_altitude(latitude(5), longitude(8)) = [[--, ..., --]] m
>>> r = src.regrids(dst, method='linear', return_operator=True) >>> y = src.regrids(r) >>> y.equals(x) True