ants.analysis package#

Many of the routines you will need to perform the processing of source data into an ancillary field are found in ants.analysis. You can also call iris routines to help your ancillary field processing. Some ANTS functions are relatively minimal wrappers around iris functions. Where appropriate, call the ANTS functions as these provide additional benefits such as updating meta-data on the cubes.

Some of the most common routines you will want to use are:

  • Calculating the area weighted mean of a field using mean()

  • Filling a field to ensure it has values on every land/ocean/both point using UMSpiralSearch().

  • Merging two datasets using merge().

Note, ANTS merge is used in a different context to iris. In ANTS it means combining two data sources either

  1. to build a larger coverage than the source data or

  2. to embed a high quality local source within a dataset with greater coverage.

The meaning in iris is discussed in the iris documentation.

class ants.analysis.FillABC(source, search_mask=None, target_mask=None)[source]#

Bases: ABC

An abstract base class which provides an interface for fill operations.

A custom fill class may be created by inheriting this base class, providing a concrete implementation of the _calculate_fill_lookup() method.

__call__(destination)[source]#

Apply the cached results of the search algorithm to resolve points in the destination.

Parameters:

destination (Cube) – Destination dataset with unresolved points. That is, missing points or perhaps coastlines which do not match a provided target landsea mask. This destination need not be the same as the source used to instantiate the class, however it must be compatible (the same grid and with the same points requiring filling).

Returns:

In-place operation.

Return type:

None

Raises:
  • ValueError – When the provided destination coordinates do not match those cached for the nearest neighbour search.

  • ValueError – If the destination mask is not compatible with the cached nearest neighbours.

__init__(source, search_mask=None, target_mask=None)[source]#

Fill missing data in the source with their nearest acceptable neighbours. Nearest acceptable neighbours include all locations in the source which are unmasked. Optionally, we can provide a search_mask to further constrain the search (thereby not having to override the source mask). Missing data is those masked source locations which are unmasked in the target (target_mask). All masked source locations which are masked in the target mask, will remain masked.

Return a callable which applies pre-calculated (cached) results to the destination, calculated at the point of instantiation of this object with the provided source-target pair.

Parameters:
  • source (Cube) – Source dataset with unresolved points. That is, missing points or perhaps coastlines which do not match a provided target landsea mask.

  • search_mask (Cube, optional) – A mask used to constrain the set of points which are considered valid by the search in the input source.

  • target_mask (Cube, optional) – Target mask field. The source will inherit this mask. Often this represents a land or ocean mask.

Returns:

Callable object for applying the spiral search algorithm.

Return type:

Callable

Raises:
  • ValueError – Where no suitable neighbour has been found.

  • ValueError – Where source, search_mask or target_mask are incompatible (not defined on identical horizontal grids).

abstract _calculate_fill_lookup(unresolved_mask, to_fill_mask, target_mask)[source]#

Calculate the cached lookup arrays to be used to fill the destination cube.

Parameters:
  • unresolved_mask (ndarray) –

    A 2-dimensional boolean array representing the unresolved points in the source cube. These points are not considered valid fill candidates. Unresolved points are those which are masked in the source, or have been specified as points to ignore in the search_mask.

    • True if point is unresolved (not valid fill candidate)

    • False if point is resolved (valid fill candidate)

    The array must be of the same shape as the horizontal dimensions of the source cube.

  • to_fill_mask (ndarray) –

    A 2-dimensional boolean array representing the missing points in the source cube that require filling. Missing points are those which are masked in the source and not masked in the target_mask.

    • True if point is missing (requires filling)

    • False if point is not missing (does not require filling)

    The array must be of the same shape as the horizontal dimensions of the source cube.

  • target_mask (ndarray) – A 2-dimensional boolean array representing the target mask. The array must be of the same shape as the horizontal dimensions of the source cube.

Returns:

First item corresponds to indices of missing data point, while the second item corresponds to indices of those elements which are to populate these missing data points.

Return type:

tuple(ndarray, ndarray)

Raises:
  • ValueError – If the provided source contains no valid data.

  • ValueError – If the algorithm cannot identify any valid neighbours for an unresolved point.

_validate_horizontal_grids(cube_name, cube)[source]#

Check a cube’s horizontal grid is consistent with that of the source cube.

Parameters:
  • cube_name (str) – Name to be used in the error message.

  • cube (Cube) – Cube to validate. The cube’s x and y coordinates are compared against those of the source cube used to initialise the class.

Raises:

ValueError – If either the x or y coordinates do not match those of the original source cube.

_validate_mask(mask_name, mask)[source]#

Check a provided mask is consistent with the source cube.

The mask is validated according to the following criteria:

  • The mask must be 2-dimensional.

  • The mask’s horizontal grid must be consistent with that of the source cube used to instantiate the class.

Parameters:
  • mask_name (str) – Name to be used in the error message.

  • mask (Cube) – Mask to validate.

Raises:
  • ValueError – If the mask is not 2-dimensional.

  • ValueError – If the mask’s horizontal grid is inconsistent with that of the source cube.

search(source)[source]#

Search for points to fill the missing data in the source.

Perform pre-calculation steps up front, independent of which underlying search algorithm is used:

  • Reshape data

  • Identify points to be filled (to_fill_mask)

  • Identify points that should be excluded from the search (unresolved_mask)

  • Call _calculate_fill_lookup(), which wraps the specific underlying search algorithm

  • Cache the resulting lookup arrays, to be applied when the filler is called

Parameters:

source (Cube) – Source cube with unresolved points, identified by being masked or NaN.

Raises:
  • ValueError – If the provided source contains no valid data.

  • ValueError – If no valid neighbours can be identified for any cell.

class ants.analysis.FillMissingPoints(source, search_mask=None, target_mask=None)[source]#

Bases: UMSpiralSearch

Warning

This class is deprecated as of ANTS 2.0. The functionality has been moved to UMSpiralSearch.

__init__(source, search_mask=None, target_mask=None)[source]#

Fill missing data in the source with their nearest acceptable neighbours. Nearest acceptable neighbours include all locations in the source which are unmasked. Optionally, we can provide a search_mask to further constrain the search (thereby not having to override the source mask). Missing data is those masked source locations which are unmasked in the target (target_mask). All masked source locations which are masked in the target mask, will remain masked.

Return a callable which applies pre-calculated (cached) results to the destination, calculated at the point of instantiation of this object with the provided source-target pair.

Parameters:
  • source (Cube) – Source dataset with unresolved points. That is, missing points or perhaps coastlines which do not match a provided target landsea mask.

  • search_mask (Cube, optional) – A mask used to constrain the set of points which are considered valid by the search in the input source.

  • target_mask (Cube, optional) – Target mask field. The source will inherit this mask. Often this represents a land or ocean mask.

Returns:

Callable object for applying the spiral search algorithm.

Return type:

Callable

Raises:
  • ValueError – Where no suitable neighbour has been found.

  • ValueError – Where source, search_mask or target_mask are incompatible (not defined on identical horizontal grids).

class ants.analysis.KDTreeFill(source, search_mask=None, target_mask=None)[source]#

Bases: FillABC

K-D tree missing data fill and coastal adjustment algorithm.

See also

https://en.wikipedia.org/wiki/K-d_tree :

for details of the K-D tree search method.

__call__(destination)[source]#

Apply the cached results of the search algorithm to resolve points in the destination.

Parameters:

destination (Cube) – Destination dataset with unresolved points. That is, missing points or perhaps coastlines which do not match a provided target landsea mask. This destination need not be the same as the source used to instantiate the class, however it must be compatible (the same grid and with the same points requiring filling).

Returns:

In-place operation.

Return type:

None

Raises:
  • ValueError – When the provided destination coordinates do not match those cached for the nearest neighbour search.

  • ValueError – If the destination mask is not compatible with the cached nearest neighbours.

__init__(source, search_mask=None, target_mask=None)[source]#

Fill missing data in the source with their nearest acceptable neighbours. Nearest acceptable neighbours include all locations in the source which are unmasked. Optionally, we can provide a search_mask to further constrain the search (thereby not having to override the source mask). Missing data is those masked source locations which are unmasked in the target (target_mask). All masked source locations which are masked in the target mask, will remain masked.

Return a callable which applies pre-calculated (cached) results to the destination, calculated at the point of instantiation of this object with the provided source-target pair.

Parameters:
  • source (Cube) – Source dataset with unresolved points. That is, missing points or perhaps coastlines which do not match a provided target landsea mask.

  • search_mask (Cube, optional) – A mask used to constrain the set of points which are considered valid by the search in the input source.

  • target_mask (Cube, optional) – Target mask field. The source will inherit this mask. Often this represents a land or ocean mask.

Returns:

Callable object for applying the spiral search algorithm.

Return type:

Callable

Raises:
  • ValueError – Where no suitable neighbour has been found.

  • ValueError – Where source, search_mask or target_mask are incompatible (not defined on identical horizontal grids).

_calculate_fill_lookup(unresolved_mask, to_fill_mask, target_mask)[source]#

Call the K-D tree search algorithm to calculate the fill lookup.

class ants.analysis.MooreNeighbourhood(acube, window_mask=None)[source]#

Bases: object

Operations on the Moore Neighbourhood of points in a cube.

__init__(acube, window_mask=None)[source]#

Each grid box in a cube has a Moore Neighbourhood: the 8 grid boxes around it on the grid. An instance of this class can be used to calculate some of the properties of the neighbourhood of the grid boxes of a single cube.

Parameters:
  • acube (Cube) – Calculations of the Moore neighbourhood will be done on this cube.

  • window_mask (bool ndarray object, optional) – Mask array which defines the Moore neighbourhood of each in the window of each grid point. The provided mask must have an odd length for each dimension. If not provided, this is assumed to be a 3x3 with all surrounding points considered.

Raises:
  • ValueError – If the provided mask does not have an odd length for each dimension.

  • ValueError – If the cube is not 2 dimensional. This is could be relaxed with a more advanced implementation.

  • ValueError – If the axes of the cube are not contiguous, and so do not look like they are neighbours.

References

https://en.wikipedia.org/wiki/Moore_neighborhood

all_equal_value(value)[source]#

Returns True for grid boxes where all neighbouring grid boxes have the specified value.

Parameters:

value (Number) – Choosing a type that can be compared to the cube data type, the value is compared to neighbourhood values.

Returns:

True where the neighbours are all equal to value.

Return type:

2d bool type ndarray

Notes

The comparison is exact; a tolerance isn’t currently supported.

all_missing()[source]#

Returns True for grid boxes where the neighbourhood contains only the missing value.

any_equal_value(value)[source]#

Return True for grid boxes where any of the Moore neighbours have the specified value.

Parameters:

value (Number or ndarray) – The value(s) to compare neighbourhood to. If an ndarray, this should match the mask shape.

Returns:

True where the neighbours are all equal to value.

Return type:

2d bool type ndarray

Notes

The comparison is exact; a tolerance isn’t currently supported.

any_non_missing()[source]#

Returns True for grid boxes where the Moore neighbourhood contains any valid data.

neighbourhood_mean()[source]#

Return the mean of the values in the Moore neighbourhood.

If the neighbourhood contains missing data then the mean is taken over only the neighbours that are not missing.

class ants.analysis.UMSpiralSearch(source, search_mask=None, target_mask=None)[source]#

Bases: FillABC

Spiral search missing data fill and coastal adjustment algorithm.

This uses the spiral search algorithm or coast_adj_method=3 as it is known in the community.

See also

https://code.metoffice.gov.uk/doc/um/latest/papers/umdp_S11.pdf : for details of the spiral search method.

Warning

Performing a spiral search for circular fields (wraparound) is known to have performance implications.

__call__(destination)[source]#

Apply the cached results of the search algorithm to resolve points in the destination.

Parameters:

destination (Cube) – Destination dataset with unresolved points. That is, missing points or perhaps coastlines which do not match a provided target landsea mask. This destination need not be the same as the source used to instantiate the class, however it must be compatible (the same grid and with the same points requiring filling).

Returns:

In-place operation.

Return type:

None

Raises:
  • ValueError – When the provided destination coordinates do not match those cached for the nearest neighbour search.

  • ValueError – If the destination mask is not compatible with the cached nearest neighbours.

__init__(source, search_mask=None, target_mask=None)[source]#

Fill missing data in the source with their nearest acceptable neighbours. Nearest acceptable neighbours include all locations in the source which are unmasked. Optionally, we can provide a search_mask to further constrain the search (thereby not having to override the source mask). Missing data is those masked source locations which are unmasked in the target (target_mask). All masked source locations which are masked in the target mask, will remain masked.

Return a callable which applies pre-calculated (cached) results to the destination, calculated at the point of instantiation of this object with the provided source-target pair.

Parameters:
  • source (Cube) – Source dataset with unresolved points. That is, missing points or perhaps coastlines which do not match a provided target landsea mask.

  • search_mask (Cube, optional) – A mask used to constrain the set of points which are considered valid by the search in the input source.

  • target_mask (Cube, optional) – Target mask field. The source will inherit this mask. Often this represents a land or ocean mask.

Returns:

Callable object for applying the spiral search algorithm.

Return type:

Callable

Raises:
  • ValueError – Where no suitable neighbour has been found.

  • ValueError – Where source, search_mask or target_mask are incompatible (not defined on identical horizontal grids).

_calculate_fill_lookup(unresolved_mask, to_fill_mask, target_mask)[source]#

Call the spiral search algorithm to calculate the fill lookup.

ants.analysis.calc_grad(source)[source]#

Attention

The calc_grad routine has been removed from the core ants library at version 2.2. It has been moved to Apps/Orography/orography_utils.py in the ancillary-file-science repository. Attempting to use this function will result in an ImportError.

ants.analysis.find_similar_region(array, seed_point, extended_neighbourhood=False, wraparound=False)[source]#

Return a set of indices where the connecting neighbours have the same value

This function is functionaly equivelent floodfill(), except that here the fill locations are returned rather than filled.

Parameters:
  • array (ndarray) – The array to search.

  • seed_point (tuple) – The starting (y, x) index (the seed point).

  • extended_neighbourhood (bool, optional) –

    In the extended neighbourhood case, also consider the diagonals in each locations neighbourhood:

    Default neighbourhood:
    [False, True, False]
    [True, True, True ]
    [False, True, False]
    Extended Neighbourhood:
    [True, True, True]
    [True, True, True]
    [True, True, True]

  • wraparound (bool, optional) – When True, support wraparound in ‘x’, otherwise stop at the boundary.

Returns:

2D array containing the row and column indices of those locations identified as similar.

Return type:

ndarray

ants.analysis.find_small_feature_seed_points(array, min_feature_size, land_features=False, extended_neighbourhood=False, wraparound=False)[source]#

Return a set of indices for seed points for filling regions of less than a specified number of points in a provided land sea mask.

Deprecated since version 2.2: Planned for removal due to no known use cases. Please contact miao@metoffice.gov.uk if you do need this functionality.

Parameters:
  • array (ndarray) – The array to search.

  • min_feature_size (int) – Find features with fewer grid boxes than this

  • land_features (bool) – Set to True to identify land seed points, False (default) for sea seed points.

  • extended_neighbourhood (bool, optional) –

    In the extended neighbourhood case, also consider the diagonals in each locations neighbourhood:

    Default neighbourhood:
    [False, True, False]
    [True, True, True ]
    [False, True, False]
    Extended Neighbourhood:
    [True, True, True]
    [True, True, True]
    [True, True, True]

  • wraparound (bool, optional) – When True, support wraparound in ‘x’, otherwise stop at the boundary.

Returns:

2D array containing the row and column indices of those locations identified as potential fill seed points for features smaller than the specified feature_size.

Return type:

ndarray

ants.analysis.floodfill(array, seed_point, fill_value, extended_neighbourhood=False, wraparound=False)[source]#

Floodfill via an iterative algorithm.

Parameters:
  • array (ndarray) – The array to apply the floodfill.

  • seed_point (tuple) – The starting (y, x) index (the seed point).

  • fill_value (int or float) – The value which results in the ‘flooded area’.

  • extended_neighbourhood (bool, optional) –

    In the extended neighbourhood case, also consider the diagonals in each locations neighbourhood:

    Default neighbourhood:
    [False, True, False]
    [True, True, True ]
    [False, True, False]
    Extended Neighbourhood:
    [True, True, True]
    [True, True, True]
    [True, True, True]

  • wraparound (bool, optional) – When True, support wraparound in ‘x’, otherwise stop at the boundary.

ants.analysis.horizontal_grid_reorder(cube)[source]#

Return a transposed view of the cube’s data, so that the outermost dimension corresponds to ‘y’ while the next corresponds to ‘x’.

Parameters:

cube (Cube) – Cube to transpose.

Returns:

A tuple of arrays, both transposed as [y, x, …]:

  • The first array corresponds to the cube’s data.

  • The second array corresponds to the cube’s mask.

Return type:

tuple(ndarray, ndarray)

Note

Both array view and mask views are returned. This minimises memory use by ensuring the mask is a view. We can not take the mask from the masked array view as this mask is not a view.

Warning

Care should be taken when reordering when there are more than three dimensions as the order of the dimension which are not mapped to the horizontal grid might not agree but be of the same length.

ants.analysis.make_consistent_with_lsm(sources, lsm, invert_mask, method='spiral')[source]#

Make the provided source(s) consistent with the provided land sea mask.

Replaces missing data values with valid data, where missing is defined as data that is either masked or NaN. Missing land points are filled with nearest valid land values, and missing sea points are filled with nearest valid sea values. Land and sea points are defined by the provided landsea mask.

Additional coordinates on the source(s), such as time and pseudo levels are handled. However, if the dimensions are not already in the order of (<other dimensions>, y, x), they will be rearranged to this order.

Parameters:
  • sources (Iterable of Cube objects) – Source cubes.

  • lsm (Cube) – Landsea mask field.

  • invert_mask (bool) – Invert the mask (land field) or not (ocean field). The landsea mask has True values to denote land.

  • method (str, optional) – Select the search method to be used when filling missing points. The methods currently supported are “spiral” and “kdtree”.

Returns:

In-place operation.

Return type:

None

Raises:

NotImplementedError – If an unsupported search method is specified.

See also

ants.analysis.UMSpiralSearch

for more details on how the nearest valid points are determined when the “spiral” method is used.

ants.analysis.KDTreeFill

for more details on how the nearest valid points are determined when the “kdtree” method is used.

ants.analysis.mean(source, target)[source]#

Calculate the source mean \(\bar{x}\) from the area mean of the source \((x_i)\) over each target grid box.

\[\bar{x} = \sum_i {x_i}{w_i}\]

where \(w_i\) is the area of the source grid box (given an index i) that overlaps the target grid box divided by the area of the target grid box.

Parameters:
  • source (Cube) – Source cube.

  • target (Cube) – Target cube.

Returns:

The area-weighted mean of the source over the target grid.

Return type:

Cube

ants.analysis.merge(primary_cube, alternate_cube, validity_polygon=None)[source]#

Merges data from the alternative cube into the primary cube.

The primary cube data is used as a base, then cells from the alternate cube which lay outside the provided polygon, override the values of the primary at those locations. Containment is defined as any cell corner which lies within the polygon. “Within” explicitly does not include those points which exactly lay on the polygon boundary. Where multiple primary and alternate cubes are provided, then these are paired appropriately where possible. Where these datasets are not defined on the same grid, the user should consider a regrid first to then utilise merge.

Parameters:
  • primary_cube (Cube or CubeList) – The primary dataset which has a highest priority i.e. overriding values in the alternate dataset.

  • alternate_cube (Cube or CubeList) – The alternate data set which is to be merged, taking a lower priority to values contained within the primary dataset.

  • validity_polygon (Iterable or shapely.Polygon instance, optional) – Polygon defining the region of valid data within the primary dataset. Data defined outside this region will not take preference of the alternate dataset. The crs of the polygon is assumed to be the same as the primary dataset to which it describes. If an iterable is provided, then each item should correspond to the x, y point definition. If not provided, the entire primary_cube dataset is considered valid (including masked value). This means that the two datasets are stacked together with the primary_cube taking priority over alternate_cube in the case of an overlap. A runtime error will be raised if the primary_cube is wholly within the validity_polygon.

Returns:

One or more Cube objects

Return type:

CubeList

Raises:

RuntimeError – If the primary cube is wholly within a provided validity_polygon.

ants.analysis.standard_deviation(source, src_mean)[source]#

Calculate the standard deviation of the source within the target grid box.

This is calculated as the square root of the variance. The variance is defined here as the mean of the squares minus the square of the mean.

\[\sigma^2 = \sum_i {x_i}^2{w_i} - \bar{x}^2\]
Parameters:
  • source (Cube) – Source cube.

  • src_mean (Cube) – Mean of the source field on the target grid.

Returns:

The sub-grid standard deviation of the source over the target grid.

Return type:

Cube

Notes

There is no correction for biases in the standard deviation (Bessel’s correction). This could be an issue when the target grid is much lower resolution than the source grid.

Submodules#

ants.analysis.cover_mapping module#

ants.analysis.cover_mapping.fetch_lct_slices(source, um_tile_ids)[source]#

Fetch the slices corresponding to the specified JULE tile ID.

Given an iris cube, derive a set of slices for this cube corresponding to this JULE tile ID. That is, indexing the pseudo_level mapped dimension.

Parameters:
  • source (Cube) – Source cube which has land cover types as defined by a pseudo level coordinate.

  • um_tile_ids (int | list[int]) – JULES tile IDs. This corresponds to the pseudo_level value(s) desired.

Returns:

A tuple containing slice objects.

Return type:

tuple

Example

>>> c3_grass = 3
>>> c3_grass_slices = fetch_veg_slice(cube, c3_grass)
>>> c3_grass.data[c3_grass] = ...
ants.analysis.cover_mapping.normalise_fractions(source)[source]#

Normalisation of fractions, ensuring the sum of the fractions is 1.

Normalisation effectively works by filling missing data by dividing equally amongst all types where data > 0 such that the ratios between fractions remain the same. Where there are no class fractions at a given point, it will remain 0. Fractions outside the range [0, 1] are considered anomalous and pulled into that range before normalisation occurs.

Parameters:

source (iris.cube.Cube) – Source land cover type fraction, with pseudo-level coordinate representing the classes.

Warning

Mask is not altered by this function.

ants.analysis.filters module#

Deprecated since version 2.2: This module is deprecated. Filtering functions can now be found in the ancillary-file-science repository, under Apps/Orography/orography_filters.py

ants.analysis.filters.raymond(source, epsilon=None, filter_length_scale=None, isotropic=False)[source]#

Attention

The Raymond filter has been removed from the core ants library at version 2.2. It has been moved to Apps/Orography/orography_filters.py in the ancillary-file-science repository. Attempting to use this function will result in an ImportError.