cf.Field.regrids¶
-
Field.
regrids
(dst, method, 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, i=False, _compute_field_mass=None)[source]¶ Return the field regridded onto a new latitude-longitude grid.
Regridding, also called remapping or interpolation, is the process of changing the grid underneath field data values while preserving the qualities of the original data.
The regridding method must be specified. First-order conservative interpolation conserves the global area integral of the field, but may not give approximations to the values as good as bilinear interpolation. Bilinear interpolation is available. The latter method is particular useful for cases when the latitude and longitude coordinate cell boundaries are not known nor inferrable. Higher order patch recovery is available as an alternative to bilinear interpolation. This typically results in better approximations to values and derivatives compared to the latter, but the weight matrix can be larger than the bilinear matrix, which can be an issue when regridding close to the memory limit on a machine. Nearest neighbour interpolation is also available. Nearest source to destination is particularly useful for regridding integer fields such as land use.
Metadata
The field’s domain must have well defined X and Y axes with latitude and longitude coordinate values, which may be stored as dimension coordinate objects or two dimensional auxiliary coordinate objects. If the latitude and longitude coordinates are two dimensional then the X and Y axes must be defined by dimension coordinates if present or by the netCDF dimensions. In the latter case the X and Y axes must be specified using the src_axes or dst_axes keyword. The same is true for the destination grid, if it provided as part of another field.
The cyclicity of the X axes of the source field and destination grid is taken into account. If an X axis is in fact cyclic but is not registered as such by its parent field (see
cf.Field.iscyclic
), then the cyclicity may be set with the src_cyclic or dst_cyclic parameters. In the case of two dimensional latitude and longitude dimension coordinates without bounds it will be necessary to specify src_cyclic or dst_cyclic manually if the field is global.The output field’s coordinate objects which span the X and/or Y axes are replaced with those from the destination grid. Any fields contained in coordinate reference objects will also be regridded, if possible.
Mask
The data array mask of the field is automatically taken into account, such that the regridded data array will be masked in regions where the input data array is masked. By default the mask of the destination grid is not taken into account. If the destination field data has more than two dimensions then the mask, if used, is taken from the two dimensionsal section of the data where the indices of all axes other than X and Y are zero.
Implementation
The interpolation is carried by out using the
ESMF
package - a Python interface to the Earth System Modeling Framework (ESMF) regridding utility.Logging
Whether ESMF logging is enabled or not is determined by
cf.REGRID_LOGGING
. If it is logging takes place after every call. By default logging is disabled.Latitude-Longitude Grid
The canonical grid with independent latitude and longitude coordinates.
Curvilinear Grids
Grids in projection coordinate systems can be regridded as long as two dimensional latitude and longitude coordinates are present.
Rotated Pole Grids
Rotated pole grids can be regridded as long as two dimensional latitude and longitude coordinates are present. It may be necessary to explicitly identify the grid latitude and grid longitude coordinates as being the X and Y axes and specify the src_cyclic or dst_cyclic keywords.
Tripolar Grids
Tripolar grids are logically rectangular and so may be able to be regridded. If no dimension coordinates are present it will be necessary to specify which netCDF dimensions are the X and Y axes using the src_axes or dst_axes keywords. Connections across the bipole fold are not currently supported, but are not be necessary in some cases, for example if the points on either side are together without a gap. It will also be necessary to specify src_cyclic or dst_cyclic if the grid is global.
New in version 1.0.4.
Examples 1: Regrid field
f
conservatively onto a grid contained in fieldg
:>>> h = f.regrids(g, 'conservative')
Parameters: - dst:
cf.Field
ordict
The field containing the new grid. If dst is a field list the first field in the list is used. Alternatively a dictionary can be passed containing the keywords ‘longitude’ and ‘latitude’ with either two 1D dimension coordinates or two 2D auxiliary coordinates. In the 2D case both coordinates must have their axes in the same order and this must be specified by the keyword ‘axes’ as either of the tuples
('X', 'Y')
or('Y', 'X')
.- method:
str
Specify the regridding method. The method parameter must be one of:
method Description 'bilinear'
Bilinear interpolation. 'patch'
Higher order patch recovery. 'conservative_1st'
First-order conservative regridding or 'conservative'
will be used (requires both of the fields to have contiguous, non- overlapping bounds). 'nearest_stod'
Nearest neighbor interpolation is used where each destination point is mapped to the closest source point 'nearest_dtos'
Nearest neighbor interpolation is used where each source point is mapped to the closest destination point. A given destination point may receive input from multiple source points, but no source point will map to more than one destination point. - src_cyclic:
bool
, optional Specifies whether the longitude for the source grid is periodic or not. If None then, if possible, this is determined automatically otherwise it defaults to False.
- dst_cyclic:
bool
, optional Specifies whether the longitude for the destination grid is periodic of not. If None then, if possible, this is determined automatically otherwise it defaults to False.
- use_src_mask:
bool
, optional For all methods other than ‘nearest_stod’, this must be True as it does not make sense to set it to False. For the ‘nearest_stod’ method if it is True then points in the result that are nearest to a masked source point are masked. Otherwise, if it is False, then these points are interpolated to the nearest unmasked source points.
- use_dst_mask:
bool
, optional By default the mask of the data on the destination grid is not taken into account when performing regridding. If this option is set to true then it is. If the destination field has more than two dimensions then the first 2D slice in index space is used for the mask e.g. for an field varying with (X, Y, Z, T) the mask is taken from the slice (X, Y, 0, 0).
- fracfield:
bool
, optional If the method of regridding is conservative the fraction of each destination grid cell involved in the regridding is returned instead of the regridded data if this is True. Otherwise this is ignored.
- src_axes:
dict
, optional A dictionary specifying the axes of the 2D latitude and longitude coordinates of the source field when no dimension coordinates are present. It must have keys ‘X’ and ‘Y’.
- Example:
src_axes={'X': 'ncdim%x', 'Y': 'ncdim%y'}
- dst_axes:
dict
, optional A dictionary specifying the axes of the 2D latitude and longitude coordinates of the destination field when no dimension coordinates are present. It must have keys ‘X’ and ‘Y’.
- Example:
dst_axes={'X': 'ncdim%x', 'Y': 'ncdim%y'}
- axis_order: sequence, optional
A sequence of items specifying dimension coordinates as retrieved by the
dim
method. These determine the order in which to iterate over the other axes of the field when regridding X-Y slices. The slowest moving axis will be the first one specified. Currently the regridding weights are recalculated every time the mask of an X-Y slice changes with respect to the previous one, so this option allows the user to minimise how frequently the mask changes.- ignore_degenerate:
bool
, optional True by default. Instructs ESMPy to ignore degenerate cells when checking the grids for errors. Regridding will proceed and degenerate cells will be skipped, not producing a result, when set to True. Otherwise an error will be produced if degenerate cells are found. This will be present in the ESMPy log files if
cf.REGRID_LOGGING
is set to True. As of ESMF 7.0.0 this only applies to conservative regridding. Other methods always skip degenerate cells.- i:
bool
, optional If True then update the field in place. By default a new field is created. In either case, a field is returned.
- _compute_field_mass:
dict
, optional If this is a dictionary then the field masses of the source and destination fields are computed and returned within the dictionary. The keys of the dictionary indicates the lat/long slice of the field and the corresponding value is a tuple containing the source field’s mass and the destination field’s mass. The calculation is only done if conservative regridding is being performed. This is for debugging purposes.
Returns: - out:
cf.Field
The regridded field.
Examples 2: Regrid f to the grid of g using bilinear regridding and forcing the source field f to be treated as cyclic.
>>> h = f.regrids(g, src_cyclic=True, method='bilinear')
Regrid f to the grid of g using the mask of g.
>>> h = f.regrids(g, 'conservative_1st', use_dst_mask=True)
Regrid f to 2D auxiliary coordinates lat and lon, which have their dimensions ordered ‘Y’ first then ‘X’.
>>> lat <CF AuxiliaryCoordinate: latitude(110, 106) degrees_north> >>> lon <CF AuxiliaryCoordinate: longitude(110, 106) degrees_east> >>> h = f.regrids({'longitude': lon, 'latitude': lat, 'axes': ('Y', 'X')}, 'conservative')
Regrid field, f, on tripolar grid to latitude-longitude grid of field, g.
>>> h = f.regrids(g, 'bilinear, src_axes={'X': 'ncdim%x', 'Y': 'ncdim%y'}, ... src_cyclic=True)
Regrid f to the grid of g iterating over the ‘Z’ axis last and the ‘T’ axis next to last to minimise the number of times the mask is changed.
>>> h = f.regrids(g, 'nearest_dtos', axis_order='ZT')
Examples 3: Regrid f to a constructed 2 degree spherical grid with Voronoi bounds conservatively.
>>> lon = cf.DimensionCoordinate() >>> lat = cf.DimensionCoordinate()
>>> lon.standard_name = 'longitude' >>> lat.standard_name = 'latitude'
>>> import numpy as np >>> lon.insert_data(cf.Data(np.arange(0, 360, 2), 'degrees_east')) >>> lat.insert_data(cf.Data(np.arange(-90, 91, 2), 'degrees_north'))
>>> lon.get_bounds(create=True, insert=True) <CF Bounds: (180, 2) degrees_east> >>> lat.get_bounds(create=True, insert=True, max=90, min=-90) <CF Bounds: (90, 2) degrees_north>
>>> g = f.regrids({'longitude': lon, 'latitude': lat}, method='conservative')
- dst: