iris.experimental.regrid

Regridding functions.

In this module:

iris.experimental.regrid.regrid_area_weighted_rectilinear_src_and_grid(src_cube, grid_cube, mdtol=0)[source]

Return a new cube with data values calculated using the area weighted mean of data values from src_grid regridded onto the horizontal grid of grid_cube.

This function requires that the horizontal grids of both cubes are rectilinear (i.e. expressed in terms of two orthogonal 1D coordinates) and that these grids are in the same coordinate system. This function also requires that the coordinates describing the horizontal grids all have bounds.

Note

Elements in data array of the returned cube that lie either partially or entirely outside of the horizontal extent of the src_cube will be masked irrespective of the value of mdtol.

Args:

  • src_cube:

    An instance of iris.cube.Cube that supplies the data, metadata and coordinates.

  • grid_cube:

    An instance of iris.cube.Cube that supplies the desired horizontal grid definition.

Kwargs:

  • mdtol:

    Tolerance of missing data. The value returned in each element of the returned cube’s data array will be masked if the fraction of masked data in the overlapping cells of the source cube exceeds mdtol. This fraction is calculated based on the area of masked cells within each target cell. mdtol=0 means no missing data is tolerated while mdtol=1 will mean the resulting element will be masked if and only if all the overlapping cells of the source cube are masked. Defaults to 0.

Returns

A new iris.cube.Cube instance.

↑ top ↑

iris.experimental.regrid.regrid_weighted_curvilinear_to_rectilinear(src_cube, weights, grid_cube)[source]

Return a new cube with the data values calculated using the weighted mean of data values from src_cube and the weights from weights regridded onto the horizontal grid of grid_cube.

This function requires that the src_cube has a horizontal grid defined by a pair of X- and Y-axis coordinates which are mapped over the same cube dimensions, thus each point has an individually defined X and Y coordinate value. The actual dimensions of these coordinates are of no significance. The src_cube grid cube must have a normal horizontal grid, i.e. expressed in terms of two orthogonal 1D horizontal coordinates. Both grids must be in the same coordinate system, and the grid_cube must have horizontal coordinates that are both bounded and contiguous.

Note that, for any given target grid_cube cell, only the points from the src_cube that are bound by that cell will contribute to the cell result. The bounded extent of the src_cube will not be considered here.

A target grid_cube cell result will be calculated as, \(\sum (src\_cube.data_{ij} * weights_{ij}) / \sum weights_{ij}\), for all \(ij\) src_cube points that are bound by that cell.

Warning

  • All coordinates that span the src_cube that don’t define the horizontal curvilinear grid will be ignored.

Args:

  • src_cube:

    A iris.cube.Cube instance that defines the source variable grid to be regridded.

  • weights (array or None):

    A numpy.ndarray instance that defines the weights for the source variable grid cells. Must have the same shape as the X and Y coordinates. If weights is None, all-ones will be used.

  • grid_cube:

    A iris.cube.Cube instance that defines the target rectilinear grid.

Returns

A iris.cube.Cube instance.

↑ top ↑

This class describes the point-in-cell regridding scheme for use typically with iris.cube.Cube.regrid().

Warning

This class is now disabled.

The functionality has been moved to iris.analysis.PointInCell.

class iris.experimental.regrid.PointInCell(weights=None)[source]

Point-in-cell regridding scheme suitable for regridding over one or more orthogonal coordinates.

Warning

This class is now disabled.

The functionality has been moved to iris.analysis.PointInCell.

↑ top ↑

This class describes the linear regridding scheme which uses the scipy.interpolate.griddata to regrid unstructured data on to a grid.

The source cube and the target cube will be projected into a common projection for the scipy calculation to be performed.

class iris.experimental.regrid.ProjectedUnstructuredLinear(projection=None)[source]

Linear regridding scheme that uses scipy.interpolate.griddata on projected unstructured data.

Optional Args:

  • projection: cartopy.crs instance

    The projection that the scipy calculation is performed in. If None is given, a PlateCarree projection is used. Defaults to None.

regridder(src_cube, target_grid)[source]

Creates a linear regridder to perform regridding, using scipy.interpolate.griddata from unstructured source points to the target grid. The regridding calculation is performed in the given projection.

Typically you should use iris.cube.Cube.regrid() for regridding a cube. There are, however, some situations when constructing your own regridder is preferable. These are detailed in the user guide.

Does not support lazy regridding.

Args:

  • src_cube:

    The Cube defining the unstructured source points.

  • target_grid:

    The Cube defining the target grid.

Returns

callable(cube)

where cube is a cube with the same grid as src_cube that is to be regridded to the target_grid.

Return type

A callable with the interface

↑ top ↑

This class describes the nearest regridding scheme which uses the scipy.interpolate.griddata to regrid unstructured data on to a grid.

The source cube and the target cube will be projected into a common projection for the scipy calculation to be performed.

Note

The iris.analysis.UnstructuredNearest scheme performs essentially the same job. That calculation is more rigorously correct and may be applied to larger data regions (including global). This one however, where applicable, is substantially faster.

class iris.experimental.regrid.ProjectedUnstructuredNearest(projection=None)[source]

Nearest regridding scheme that uses scipy.interpolate.griddata on projected unstructured data.

Optional Args:

  • projection: cartopy.crs instance

    The projection that the scipy calculation is performed in. If None is given, a PlateCarree projection is used. Defaults to None.

regridder(src_cube, target_grid)[source]

Creates a nearest-neighbour regridder to perform regridding, using scipy.interpolate.griddata from unstructured source points to the target grid. The regridding calculation is performed in the given projection.

Typically you should use iris.cube.Cube.regrid() for regridding a cube. There are, however, some situations when constructing your own regridder is preferable. These are detailed in the user guide.

Does not support lazy regridding.

Args:

  • src_cube:

    The Cube defining the unstructured source points.

  • target_grid:

    The Cube defining the target grid.

Returns

callable(cube)

where cube is a cube with the same grid as src_cube that is to be regridded to the target_grid.

Return type

A callable with the interface