improver.utilities.spatial module

Provides support utilities.

class DifferenceBetweenAdjacentGridSquares[source]

Bases: BasePlugin

Calculate the difference between adjacent grid squares within a cube. The difference is calculated along the x and y axis individually.

_abc_impl = <_abc_data object>
static _update_metadata(diff_cube, coord_name, cube_name)[source]

Rename cube, add attribute and cell method to describe difference.

Parameters:
  • diff_cube (Cube) –

  • coord_name (str) –

  • cube_name (str) –

Return type:

None

static calculate_difference(cube, coord_name)[source]

Calculate the difference along the axis specified by the coordinate.

Parameters:
  • cube (Cube) – Cube from which the differences will be calculated.

  • coord_name (str) – Name of coordinate along which the difference is calculated.

Return type:

ndarray

Returns:

Array after the differences have been calculated along the specified axis.

static create_difference_cube(cube, coord_name, diff_along_axis)[source]

Put the difference array into a cube with the appropriate metadata.

Parameters:
  • cube (Cube) – Cube from which the differences have been calculated.

  • coord_name (str) – The name of the coordinate over which the difference have been calculated.

  • diff_along_axis (ndarray) – Array containing the differences.

Return type:

Cube

Returns:

Cube containing the differences calculated along the specified axis.

process(cube)[source]

Calculate the difference along the x and y axes and return the result in separate cubes. The difference along each axis is calculated using numpy.diff.

Parameters:

cube (Cube) – Cube from which the differences will be calculated.

Return type:

Tuple[Cube, Cube]

Returns:

  • Cube after the differences have been calculated along the x axis.

  • Cube after the differences have been calculated along the y axis.

class GradientBetweenAdjacentGridSquares(regrid=False)[source]

Bases: BasePlugin

Calculate the gradient between adjacent grid squares within a cube. The gradient is calculated along the x and y axis individually.

__init__(regrid=False)[source]

Initialise plugin.

Parameters:

regrid (bool) – If True, the gradient cube is regridded to match the spatial dimensions of the input cube. If False, the length of the spatial dimensions of the gradient cube are one less than for the input cube.

_abc_impl = <_abc_data object>
static _create_output_cube(gradient, diff, cube, axis)[source]

Create the output gradient cube.

Parameters:
  • gradient (ndarray) – Gradient values used in the data array of the resulting cube.

  • diff (Cube) – Cube containing differences along the x or y axis

  • cube (Cube) – Cube with correct output dimensions

  • axis (str) – Short-hand reference for the x or y coordinate, as allowed by iris.util.guess_coord_axis.

Return type:

Cube

Returns:

A cube of the gradients in the coordinate direction specified.

static _gradient_from_diff(diff, axis)[source]

Calculate the gradient along the x or y axis from differences between adjacent grid squares.

Parameters:
  • diff (Cube) – Cube containing differences along the x or y axis

  • axis (str) – Short-hand reference for the x or y coordinate, as allowed by iris.util.guess_coord_axis.

Return type:

ndarray

Returns:

Array of the gradients in the coordinate direction specified.

process(cube)[source]

Calculate the gradient along the x and y axes and return the result in separate cubes. The difference along each axis is calculated using numpy.diff.

Parameters:

cube (Cube) – Cube from which the differences will be calculated.

Return type:

Tuple[Cube, Cube]

Returns:

  • Cube after the gradients have been calculated along the x axis.

  • Cube after the gradients have been calculated along the y axis.

class OccurrenceWithinVicinity(radii=None, grid_point_radii=None, land_mask_cube=None)[source]

Bases: PostProcessingPlugin

Calculate whether a phenomenon occurs within the specified radii about a point. These radii can be given in metres, or as numbers of grid points. Each radius provided will result in a distinct output, with these demarked using a radius_of_vicinity coordinate on the resulting cube. If a single radius is provided, this will be a scalar coordinate.

Radii in metres may be used with data on a equal areas projection only. Grid_point_radii will work with any projection, with caveats.

Using grid_point_radius

A grid_point_radius will work with any projection but the effective kernel shape in real space may be irregular. An exagerated example is shown in the figure below. A vicinity radius of 1 grid point leads to a 3x3 grid point area, defined here on a very coarse lat-lon (PlateCaree) grid. This has been reprojected to an equal areas grid to demonstrate the resulting trapezoidal shape of the vicinity in “real space”.

Grid point vicinity reprojected from Lat-Lon grid to equal areas.

Users must be aware of this when choosing whether to use this plugin with a non-equal areas projection. It could introduce biases into the resulting data, e.g. related to latitude in the above example.

Vicinity processing of regridded coarse grid data

When applying vicinity processing to data sourced from a coarse grid some care is needed. The particular case of concern is coarse grid data that has been regridded to a finer mesh without any kind of downscaling. This means that the data on the finer mesh remains blocky, revealing the grid scale of the underlying data. In such cases, though a vicinity radius of order the fine mesh scale could be applied, there is no additional information on this scale; vicinity processing is effectively a convolution to a coarser grid, and it cannot be used to convolve to what is effectively a finer grid.

As such, when vicinity processing this kind of regridded but coarse data, the minimum vicinity radius that should be used is one which returns a vicinity region on the scale of the underlying coarse grid.

(vicinity radius)_min = underlying grid cell width / 2

__init__(radii=None, grid_point_radii=None, land_mask_cube=None)[source]
Parameters:
  • radii (Optional[List[Union[float, int]]]) – A list of radii in metres used to define the vicinities within which to search for occurrences.

  • grid_point_radii (Optional[List[Union[float, int]]]) – Alternatively, a list of numbers of grid points that define the vicinity radii over which to search for occurrences. Only one of radii or grid_point_radii should be set.

  • land_mask_cube (Optional[Cube]) – Binary land-sea mask data. True for land-points, False for sea. Restricts in-vicinity processing to only include points of a like mask value.

Raises:
  • ValueError – If both radii and grid point radii are set.

  • ValueError – If neither radii or grid point radii are set.

  • ValueError – If a provided vicinity radius is negative.

  • ValueError – Land mask not named land_binary_mask.

_abc_impl = <_abc_data object>
process(cube)[source]

Produces the vicinity processed data. The input data is sliced to yield y-x slices to which the maximum_within_vicinity method is applied. The different vicinity radii (if multiple) are looped over and a coordinate recording the radius used is added to each resulting cube. A single cube is returned with the leading coordinates of the input cube preserved. If a single vicinity radius is provided, a new scalar radius_of_vicinity coordinate will be found on the returned cube. If multiple radii are provided, this coordinate will be a dimension coordinate following any probabilistic / realization coordinates.

Parameters:

cube (Cube) – Thresholded cube.

Return type:

Cube

Returns:

Cube containing the occurrences within a vicinity for each radius, calculated for each yx slice, which have been merged to yield a single cube.

Raises:

ValueError – Cube and land mask have differing spatial coordinates.

add_vicinity_coordinate(cube, radius, native_grid_point_radius=False, radius_is_max=False)[source]

Add a coordinate to the cube that records the vicinity radius that has been applied to the data.

Parameters:
  • cube (Cube) – Vicinity processed cube.

  • radius (Union[float, int]) – The radius as a physical distance (m) or number of grid points, the value of which is recorded in the coordinate.

  • native_grid_point_radius (bool) – True if radius is “number of grid points”, else False

  • radius_is_max (bool) – True if the specified radius represents a maximum value from the source data. A comment is associated with the coord in this case.

Return type:

None

calculate_grid_spacing(cube, units, axis='x', rtol=1e-05)[source]

Returns the grid spacing of a given spatial axis. This will be positive for axes that stride negatively.

Parameters:
  • cube (Cube) – Cube of data on equal area grid

  • units (Union[Unit, str]) – Unit in which the grid spacing is required

  • axis (str) – Axis (‘x’ or ‘y’) to use in determining grid spacing

  • rtol (float) – relative tolerance

Return type:

float

Returns:

Grid spacing in required unit

Raises:

ValueError – If points are not equally spaced

check_if_grid_is_equal_area(cube, require_equal_xy_spacing=True)[source]

Identify whether the grid is an equal area grid, by checking whether points are equally spaced along each of the x- and y-axes. By default this function also checks whether the grid spacing is the same in both spatial dimensions.

Parameters:
  • cube (Cube) – Cube with coordinates that will be checked.

  • require_equal_xy_spacing (bool) – Flag to require the grid is equally spaced in the two spatial dimensions (not strictly required for equal-area criterion).

Raises:
  • ValueError – If coordinate points are not equally spaced along either axis (from calculate_grid_spacing)

  • ValueError – If point spacing is not equal for the two spatial axes

Return type:

None

create_vicinity_coord(radius, native_grid_point_radius=False)[source]

Create a coordinate that records the vicinity radius passed in. This radius may be a distance in physical units, or if native_grid_point_radius is True, it will be a number of grid cells. If the latter an attribute comment is added to note that the radius is in grid cells.

Parameters:
  • radius (Union[float, int]) – The radius as a physical distance or number of grid points, the value of which is recorded in the coordinate.

  • native_grid_point_radius (bool) – If set to True the radius is provided a a number of grid points and the metadata created will reflect that.

Return type:

AuxCoord

distance_to_number_of_grid_cells(cube, distance, axis='x', return_int=True)[source]

Return the number of grid cells in the x and y direction based on the input distance in metres. Requires an equal-area grid on which the spacing is equal in the x- and y- directions.

Parameters:
  • cube (Cube) – Cube containing the x and y coordinates, which will be used for calculating the number of grid cells in the x and y direction, which equates to the requested distance in the x and y direction.

  • distance (float) – Distance in metres. Must be positive.

  • return_int (bool) – If true only integer number of grid_cells are returned, rounded down. If false the number of grid_cells returned will be a float.

  • axis (str) – Axis (‘x’ or ‘y’) to use in determining grid spacing

Return type:

Union[float, int]

Returns:

Number of grid cells along the specified (x or y) axis equal to the requested distance in metres.

Raises:

ValueError – If a non-positive distance is provided.

get_grid_y_x_values(cube)[source]

Extract the y and x coordinate values of each points in the cube.

The result is defined over the spatial grid, of shape (ny, nx) where ny is the length of the y-axis coordinate and nx the length of the x-axis coordinate.

Parameters:

cube (Cube) – Cube with points to extract

Return type:

Tuple[ndarray, ndarray]

Returns:

  • Array of shape (ny, nx) containing y coordinate values

  • Array of shape (ny, nx) containing x coordinate values

lat_lon_determine(cube)[source]

Test whether a diagnostic cube is on a latitude/longitude grid or uses an alternative projection.

Parameters:

cube (Cube) – A diagnostic cube to examine for coordinate system.

Return type:

Optional[CRS]

Returns:

Coordinate system of the diagnostic cube in a cartopy format unless it is already a latitude/longitude grid, in which case None is returned.

maximum_within_vicinity(grid, grid_point_radius, landmask=None)[source]

Find grid points where a phenomenon occurs within a defined radius. The occurrences within this vicinity are maximised, such that all grid points within the vicinity are recorded as having an occurrence. For non-binary fields, if the vicinity of two occurrences overlap, the maximum value within the vicinity is chosen. If a land-mask has been supplied, process land and sea points separately.

Parameters:
  • grid (Union[MaskedArray, ndarray]) – An array of values to which the process is applied.

  • grid_point_radius (int) – The radius in grid points about each point within which to determine the maximum value.

  • landmask (Optional[ndarray]) – A binary grid of the same size as grid that differentiates between land and sea points to allow the different surface types to be processed independently.

Return type:

ndarray

Returns:

Array where the occurrences have been spatially spread, so that they’re equally likely to have occurred anywhere within the vicinity defined using the specified radius.

number_of_grid_cells_to_distance(cube, grid_points)[source]

Calculate distance in metres equal to the given number of gridpoints based on the coordinates on an input cube.

Parameters:
  • cube (Cube) – Cube for which the distance is to be calculated.

  • grid_points (int) – Number of grid points to convert.

Return type:

float

Returns:

The radius in metres.

rename_vicinity_cube(cube)[source]

Rename a cube in place to indicate the cube has been vicinity processed.

Parameters:

cube (Cube) – Cube to be renamed.

transform_grid_to_lat_lon(cube)[source]

Calculate the latitudes and longitudes of each points in the cube.

The result is defined over the spatial grid, of shape (ny, nx) where ny is the length of the y-axis coordinate and nx the length of the x-axis coordinate.

Parameters:

cube (Cube) – Cube with points to transform

Return type:

Tuple[ndarray, ndarray]

Returns
  • Array of shape (ny, nx) containing grid latitude values

  • Array of shape (ny, nx) containing grid longitude values

update_name_and_vicinity_coord(cube, new_name, vicinity_radius)[source]

Updates a cube with a new probabilistic-style name and replaces or adds a radius_of_vicinity coord with the specified radius.

Parameters:
  • cube (Cube) – Cube to be updated in-place

  • new_name (str) – The new name to be applied to the Cube, the threshold coord and any related cell methods.

  • vicinity_radius (float) – The point value for the radius_of_vicinity coord. The units are assumed to be the same as the x and y spatial coords of the Cube