esmf_regrid.esmf_regridder module#

Provides ESMF representations of grids/UGRID meshes and a modified regridder.

class esmf_regrid.esmf_regridder.GridInfo(lons, lats, lonbounds, latbounds, crs=None, circular=False, areas=None, mask=None, center=False)[source]#

Bases: SDO

Class for handling structured grids.

This class holds information about lat-lon type grids. That is, grids defined by lists of latitude and longitude values for points/bounds (with respect to some coordinate reference system i.e. rotated pole). It contains methods for translating this information into esmpy objects. In particular, there are methods for representing as a esmpy.api.grid.Grid and as a esmpy.api.field.Field containing that Grid. This esmpy Field is designed to contain enough information for area weighted regridding and may be inappropriate for other esmpy regridding schemes.

Create a GridInfo object describing the grid.

Parameters:
  • lons (ArrayLike) – A 1D or 2D array or list describing the longitudes of the grid points.

  • lats (ArrayLike) – A 1D or 2D array or list describing the latitudes of the grid points.

  • lonbounds (ArrayLike) – A 1D or 2D array or list describing the longitude bounds of the grid. Should have length one greater than lons.

  • latbounds (ArrayLike) – A 1D or 2D array or list describing the latitude bounds of the grid. Should have length one greater than lats.

  • crs (cartopy.crs.CRS, optional) – Describes how to interpret the above arguments. If None, defaults to Geodetic.

  • circular (bool, default=False) – Describes if the final longitude bounds should be considered contiguous with the first.

  • areas (ArrayLike, optional) – Array describing the areas associated with each face. If None, then esmpy will use its own calculated areas.

  • mask (ArrayLike, optional) – Array describing which elements esmpy will ignore.

  • center (bool, default=False) – Describes if the center points of the grid cells are used in regridding calculations.

class esmf_regrid.esmf_regridder.RefinedGridInfo(lonbounds, latbounds, resolution=3, crs=None, mask=None)[source]#

Bases: GridInfo

Class for handling structured grids represented in esmpy in higher resolution.

A specialised version of GridInfo. Designed to provide higher accuracy conservative regridding for rectilinear grids, especially those with particularly large cells which may not be well represented by esmpy. This class differs from GridInfo primarily in the way it represents itself as a Field in esmpy. This Field is designed to be a higher resolution version of the given grid and should contain enough information for area weighted regridding but may be inappropriate for other esmpy regridding schemes.

Create a RefinedGridInfo object describing the grid.

Parameters:
  • lonbounds (ArrayLike) – A 1D array or list describing the longitude bounds of the grid. Must be strictly increasing (for example, if a bound goes from 170 to -170 consider transposing -170 to 190).

  • latbounds (ArrayLike) – A 1D array or list describing the latitude bounds of the grid. Must be strictly increasing.

  • resolution (int, default=400) – A number describing how many latitude slices each cell should be divided into when passing a higher resolution grid to ESMF.

  • crs (cartopy.crs.CRS, optional) – Describes how to interpret the above arguments. If None, defaults to Geodetic.

class esmf_regrid.esmf_regridder.Regridder(src, tgt, method=Method.CONSERVATIVE, precomputed_weights=None)[source]#

Bases: object

Regridder for directly interfacing with esmpy.

Create a regridder from descriptions of horizontal grids/meshes.

Weights will be calculated using esmpy and stored as a scipy.sparse.csr_matrix for use in regridding. If precomputed weights are provided, these will be used instead of calculating via esmpy.

Parameters:
  • src (MeshInfo or GridInfo) – Describes the source mesh/grid. Data supplied to this regridder should be in a numpy.ndarray whose shape is compatible with src.

  • tgt (MeshInfo or GridInfo) – Describes the target mesh/grid. Data output by this regridder will be a numpy.ndarray whose shape is compatible with tgt.

  • method (Constants.Method) – The method to be used to calculate weights.

  • precomputed_weights (scipy.sparse.spmatrix, optional) – If None, esmpy will be used to calculate regridding weights. Otherwise, esmpy will be bypassed and precomputed_weights will be used as the regridding weights.

regrid(src_array, norm_type=NormType.FRACAREA, mdtol=1)[source]#

Perform regridding on an array of data.

Parameters:
  • src_array (ArrayLike) – Array whose shape is compatible with self.src

  • norm_type (Constants.NormType) – Either Constants.NormType.FRACAREA or Constants.NormType.DSTAREA. Determines the type of normalisation applied to the weights.

  • mdtol (float, default=1) – A number between 0 and 1 describing the missing data tolerance. Depending on the value of mdtol, if a cell in the target grid is not sufficiently covered by unmasked cells of the source grid, then it will be masked. mdtol=1 means that only target cells which are not covered at all will be masked, mdtol=0 means that all target cells that are not entirely covered will be masked, and mdtol=0.5 means that all target cells that are less than half covered will be masked.

Returns:

An array whose shape is compatible with self.tgt.

Return type:

ArrayLike