iris.util

Miscellaneous utility functions.

In this module:

iris.util.approx_equal(a, b, max_absolute_error=1e-10, max_relative_error=1e-10)[source]

Returns whether two numbers are almost equal, allowing for the finite precision of floating point numbers.

↑ top ↑

iris.util.array_equal(array1, array2, withnans=False)[source]

Returns whether two arrays have the same shape and elements.

Args:

  • array1, array2 (arraylike):

    args to be compared, after normalising with np.asarray().

Kwargs:

  • withnans (bool):

    When unset (default), the result is False if either input contains NaN points. This is the normal floating-point arithmetic result. When set, return True if inputs contain the same value in all elements, _including_ any NaN values.

This provides much the same functionality as numpy.array_equal(), but with additional support for arrays of strings and NaN-tolerant operation.

↑ top ↑

iris.util.as_compatible_shape(src_cube, target_cube)[source]

Return a cube with added length one dimensions to match the dimensionality and dimension ordering of target_cube.

This function can be used to add the dimensions that have been collapsed, aggregated or sliced out, promoting scalar coordinates to length one dimension coordinates where necessary. It operates by matching coordinate metadata to infer the dimensions that need modifying, so the provided cubes must have coordinates with the same metadata (see iris.common.CoordMetadata).

Note

This function will load and copy the data payload of src_cube.

Deprecated since version 3.0.0: Instead use Resolve. For example, rather than calling as_compatible_shape(src_cube, target_cube) replace with Resolve(src_cube, target_cube)(target_cube.core_data()).

Args:

  • src_cube:

    An instance of iris.cube.Cube with missing dimensions.

  • target_cube:

    An instance of iris.cube.Cube with the desired dimensionality.

Returns

A instance of iris.cube.Cube with the same dimensionality as target_cube but with the data and coordinates from src_cube suitably reshaped to fit.

↑ top ↑

iris.util.between(lh, rh, lh_inclusive=True, rh_inclusive=True)[source]

Provides a convenient way of defining a 3 element inequality such as a < number < b.

Arguments:

  • lh

    The left hand element of the inequality

  • rh

    The right hand element of the inequality

Keywords:

  • lh_inclusive - boolean

    Affects the left hand comparison operator to use in the inequality. True for <= false for <. Defaults to True.

  • rh_inclusive - boolean

    Same as lh_inclusive but for right hand operator.

For example:

between_3_and_6 = between(3, 6)
for i in range(10):
   print(i, between_3_and_6(i))


between_3_and_6 = between(3, 6, rh_inclusive=False)
for i in range(10):
   print(i, between_3_and_6(i))

↑ top ↑

iris.util.broadcast_to_shape(array, shape, dim_map)[source]

Broadcast an array to a given shape.

Each dimension of the array must correspond to a dimension in the given shape. Striding is used to repeat the array until it matches the desired shape, returning repeated views on the original array. If you need to write to the resulting array, make a copy first.

Args:

  • array (numpy.ndarray-like)

    An array to broadcast.

  • shape (list, tuple etc.):

    The shape the array should be broadcast to.

  • dim_map (list, tuple etc.):

    A mapping of the dimensions of array to their corresponding element in shape. dim_map must be the same length as the number of dimensions in array. Each element of dim_map corresponds to a dimension of array and its value provides the index in shape which the dimension of array corresponds to, so the first element of dim_map gives the index of shape that corresponds to the first dimension of array etc.

Examples:

Broadcasting an array of shape (2, 3) to the shape (5, 2, 6, 3) where the first dimension of the array corresponds to the second element of the desired shape and the second dimension of the array corresponds to the fourth element of the desired shape:

a = np.array([[1, 2, 3], [4, 5, 6]])
b = broadcast_to_shape(a, (5, 2, 6, 3), (1, 3))

Broadcasting an array of shape (48, 96) to the shape (96, 48, 12):

# a is an array of shape (48, 96)
result = broadcast_to_shape(a, (96, 48, 12), (1, 0))

↑ top ↑

iris.util.clip_string(the_str, clip_length=70, rider='...')[source]

Returns a clipped version of the string based on the specified clip length and whether or not any graceful clip points can be found.

If the string to be clipped is shorter than the specified clip length, the original string is returned.

If the string is longer than the clip length, a graceful point (a space character) after the clip length is searched for. If a graceful point is found the string is clipped at this point and the rider is added. If no graceful point can be found, then the string is clipped exactly where the user requested and the rider is added.

Args:

  • the_str

    The string to be clipped

  • clip_length

    The length in characters that the input string should be clipped to. Defaults to a preconfigured value if not specified.

  • rider

    A series of characters appended at the end of the returned string to show it has been clipped. Defaults to a preconfigured value if not specified.

Returns

The string clipped to the required length with a rider appended. If the clip length was greater than the original string, the original string is returned unaltered.

↑ top ↑

iris.util.column_slices_generator(full_slice, ndims)[source]

Given a full slice full of tuples, return a dictionary mapping old data dimensions to new and a generator which gives the successive slices needed to index correctly (across columns).

This routine deals with the special functionality for tuple based indexing e.g. [0, (3, 5), :, (1, 6, 8)] by first providing a slice which takes the non tuple slices out first i.e. [0, :, :, :] then subsequently iterates through each of the tuples taking out the appropriate slices i.e. [(3, 5), :, :] followed by [:, :, (1, 6, 8)]

This method was developed as numpy does not support the direct approach of [(3, 5), : , (1, 6, 8)] for column based indexing.

↑ top ↑

iris.util.create_temp_filename(suffix='')[source]

Return a temporary file name.

Parameters

extension. (* suffix - Optional filename) –

↑ top ↑

iris.util.delta(ndarray, dimension, circular=False)[source]

Calculates the difference between values along a given dimension.

Args:

  • ndarray:

    The array over which to do the difference.

  • dimension:

    The dimension over which to do the difference on ndarray.

  • circular:

    If not False then return n results in the requested dimension with the delta between the last and first element included in the result otherwise the result will be of length n-1 (where n is the length of ndarray in the given dimension’s direction)

    If circular is numeric then the value of circular will be added to the last element of the given dimension if the last element is negative, otherwise the value of circular will be subtracted from the last element.

    The example below illustrates the process:

    original array              -180, -90,  0,    90
    delta (with circular=360):    90,  90, 90, -270+360
    

Note

The difference algorithm implemented is forward difference:

>>> import numpy as np
>>> import iris.util
>>> original = np.array([-180, -90, 0, 90])
>>> iris.util.delta(original, 0)
array([90, 90, 90])
>>> iris.util.delta(original, 0, circular=360)
array([90, 90, 90, 90])

↑ top ↑

iris.util.demote_dim_coord_to_aux_coord(cube, name_or_coord)[source]

Demotes a dimension coordinate on the cube to an auxiliary coordinate.

The DimCoord is demoted to an auxiliary coordinate on the cube. The dimension of the cube that was associated with the DimCoord becomes anonymous. The class of the coordinate is left as DimCoord, it is not recast as an AuxCoord instance.

Args:

For example,

>>> print(cube)
air_temperature / (K)               (time: 240; latitude: 37; longitude: 49)
    Dimension coordinates:
        time                             x              -              -
        latitude                         -              x              -
        longitude                        -              -              x
    Auxiliary coordinates:
        forecast_period                  x              -              -
        year                             x              -              -
>>> demote_dim_coord_to_aux_coord(cube, "time")
>>> print(cube)
air_temperature / (K)               (-- : 240; latitude: 37; longitude: 49)
    Dimension coordinates:
        latitude                        -              x              -
        longitude                       -              -              x
    Auxiliary coordinates:
        forecast_period                 x              -              -
        time                            x              -              -
        year                            x              -              -

↑ top ↑

iris.util.describe_diff(cube_a, cube_b, output_file=None)[source]

Prints the differences that prevent compatibility between two cubes, as defined by iris.cube.Cube.is_compatible().

Args:

  • cube_a:

    An instance of iris.cube.Cube or iris.cube.CubeMetadata.

  • cube_b:

    An instance of iris.cube.Cube or iris.cube.CubeMetadata.

  • output_file:

    A file or file-like object to receive output. Defaults to sys.stdout.

Note

Compatibility does not guarantee that two cubes can be merged. Instead, this function is designed to provide a verbose description of the differences in metadata between two cubes. Determining whether two cubes will merge requires additional logic that is beyond the scope of this function.

↑ top ↑

iris.util.equalise_attributes(cubes)[source]

Delete cube attributes that are not identical over all cubes in a group.

This function simply deletes any attributes which are not the same for all the given cubes. The cubes will then have identical attributes. The given cubes are modified in-place.

Args:

  • cubes (iterable of iris.cube.Cube):

    A collection of cubes to compare and adjust.

↑ top ↑

iris.util.file_is_newer_than(result_path, source_paths)[source]

Return whether the ‘result’ file has a later modification time than all of the ‘source’ files.

If a stored result depends entirely on known ‘sources’, it need only be re-built when one of them changes. This function can be used to test that by comparing file timestamps.

Args:

  • result_path (string):

    The filepath of a file containing some derived result data.

  • source_paths (string or iterable of strings):

    The path(s) to the original datafiles used to make the result. May include wildcards and ‘~’ expansions (like Iris load paths), but not URIs.

Returns

True if all the sources are older than the result, else False.

If any of the file paths describes no existing files, an exception will be raised.

Note

There are obvious caveats to using file timestamps for this, as correct usage depends on how the sources might change. For example, a file could be replaced by one of the same name, but an older timestamp.

If wildcards and ‘~’ expansions are used, this introduces even more uncertainty, as then you cannot even be sure that the resulting list of file names is the same as the originals. For example, some files may have been deleted or others added.

Note

The result file may often be a pickle file. In that case, it also depends on the relevant module sources, so extra caution is required. Ideally, an additional check on iris.__version__ is advised.

↑ top ↑

iris.util.find_discontiguities(cube, rel_tol=1e-05, abs_tol=1e-08)[source]

Searches coord for discontiguities in the bounds array, returned as a boolean array (True where discontiguities are present).

Args:

  • cube (iris.cube.Cube):

    The cube to be checked for discontinuities in its ‘x’ and ‘y’ coordinates.

Kwargs:

  • rel_tol (float):

    The relative equality tolerance to apply in coordinate bounds checking.

  • abs_tol (float):

    The absolute value tolerance to apply in coordinate bounds checking.

Returns:

  • result (numpy.ndarray of bool) :

    true/false map of which cells in the cube XY grid have discontiguities in the coordinate points array.

    This can be used as the input array for iris.util.mask_cube().

Examples:

# Find any unknown discontiguities in your cube's x and y arrays:
discontiguities = iris.util.find_discontiguities(cube)

# Pass the resultant boolean array to `iris.util.mask_cube`
# with a cube slice; this will use the boolean array to mask
# any discontiguous data points before plotting:
masked_cube_slice = iris.util.mask_cube(cube[0], discontiguities)

# Plot the masked cube slice:
iplt.pcolormesh(masked_cube_slice)

↑ top ↑

iris.util.format_array(arr)[source]

Returns the given array as a string, using the python builtin str function on a piecewise basis.

Useful for xml representation of arrays.

For customisations, use the numpy.core.arrayprint directly.

↑ top ↑

iris.util.guess_coord_axis(coord)[source]

Returns a “best guess” axis name of the coordinate.

Heuristic categorisation of the coordinate into either label ‘T’, ‘Z’, ‘Y’, ‘X’ or None.

Args:

Returns

‘T’, ‘Z’, ‘Y’, ‘X’, or None.

↑ top ↑

iris.util.is_regular(coord)[source]

Determine if the given coord is regular.

↑ top ↑

iris.util.mask_cube(cube, points_to_mask)[source]

Masks any cells in the data array which correspond to cells marked True in the points_to_mask array.

Args:

  • cube (iris.cube.Cube):

    A 2-dimensional instance of iris.cube.Cube.

  • points_to_mask (numpy.ndarray of bool):

    A 2d boolean array of Truth values representing points to mask in the x and y arrays of the cube.

Returns:

  • result (iris.cube.Cube):

    A cube whose data array is masked at points specified by input array.

↑ top ↑

iris.util.monotonic(array, strict=False, return_direction=False)[source]

Return whether the given 1d array is monotonic.

Note that, the array must not contain missing data.

Kwargs:

  • strict (boolean)

    Flag to enable strict monotonic checking

  • return_direction (boolean)

    Flag to change return behaviour to return (monotonic_status, direction). Direction will be 1 for positive or -1 for negative. The direction is meaningless if the array is not monotonic.

Returns:

  • monotonic_status (boolean)

    Whether the array was monotonic.

    If the return_direction flag was given then the returned value will be:

    (monotonic_status, direction)

↑ top ↑

iris.util.new_axis(src_cube, scalar_coord=None)[source]

Create a new axis as the leading dimension of the cube, promoting a scalar coordinate if specified.

Args:

  • src_cube (iris.cube.Cube)

    Source cube on which to generate a new axis.

Kwargs:

  • scalar_coord (iris.coord.Coord or ‘string’)

    Scalar coordinate to promote to a dimension coordinate.

Returns

A new iris.cube.Cube instance with one extra leading dimension (length 1).

For example:

>>> cube.shape
(360, 360)
>>> ncube = iris.util.new_axis(cube, 'time')
>>> ncube.shape
(1, 360, 360)

↑ top ↑

iris.util.points_step(points)[source]

Determine whether a NumPy array has a regular step.

↑ top ↑

iris.util.promote_aux_coord_to_dim_coord(cube, name_or_coord)[source]

Promotes an AuxCoord on the cube to a DimCoord. This AuxCoord must be associated with a single cube dimension. If the AuxCoord is associated with a dimension that already has a DimCoord, that DimCoord gets demoted to an AuxCoord.

Args:

For example,

>>> print(cube)
air_temperature / (K)               (time: 240; latitude: 37; longitude: 49)
    Dimension coordinates:
        time                             x              -              -
        latitude                         -              x              -
        longitude                        -              -              x
    Auxiliary coordinates:
        forecast_period                  x              -              -
        year                             x              -              -
>>> promote_aux_coord_to_dim_coord(cube, "year")
>>> print(cube)
air_temperature / (K)               (year: 240; latitude: 37; longitude: 49)
    Dimension coordinates:
        year                             x              -              -
        latitude                         -              x              -
        longitude                        -              -              x
    Auxiliary coordinates:
        forecast_period                  x              -              -
        time                             x              -              -

↑ top ↑

iris.util.regular_step(coord)[source]

Return the regular step from a coord or fail.

↑ top ↑

iris.util.reverse(cube_or_array, coords_or_dims)[source]

Reverse the cube or array along the given dimensions.

Args:

  • cube_or_array: iris.cube.Cube or numpy.ndarray

    The cube or array to reverse.

  • coords_or_dims: int, str, iris.coords.Coord or sequence of these

    Identify one or more dimensions to reverse. If cube_or_array is a numpy array, use int or a sequence of ints, as in the examples below. If cube_or_array is a Cube, a Coord or coordinate name (or sequence of these) may be specified instead.

>>> import numpy as np
>>> a = np.arange(24).reshape(2, 3, 4)
>>> print(a)
[[[ 0  1  2  3]
  [ 4  5  6  7]
  [ 8  9 10 11]]

 [[12 13 14 15]
  [16 17 18 19]
  [20 21 22 23]]]
>>> print(reverse(a, 1))
[[[ 8  9 10 11]
  [ 4  5  6  7]
  [ 0  1  2  3]]

 [[20 21 22 23]
  [16 17 18 19]
  [12 13 14 15]]]
>>> print(reverse(a, [1, 2]))
[[[11 10  9  8]
  [ 7  6  5  4]
  [ 3  2  1  0]]

 [[23 22 21 20]
  [19 18 17 16]
  [15 14 13 12]]]

↑ top ↑

iris.util.rolling_window(a, window=1, step=1, axis=- 1)[source]

Make an ndarray with a rolling window of the last dimension

Args:

  • aarray_like

    Array to add rolling window to

Kwargs:

  • windowint

    Size of rolling window

  • stepint

    Size of step between rolling windows

  • axisint

    Axis to take the rolling window over

Returns

Array that is a view of the original array with an added dimension of the size of the given window at axis + 1.

Examples:

>>> x = np.arange(10).reshape((2, 5))
>>> rolling_window(x, 3)
array([[[0, 1, 2], [1, 2, 3], [2, 3, 4]],
       [[5, 6, 7], [6, 7, 8], [7, 8, 9]]])

Calculate rolling mean of last dimension:

>>> np.mean(rolling_window(x, 3), -1)
array([[ 1.,  2.,  3.],
       [ 6.,  7.,  8.]])

↑ top ↑

iris.util.squeeze(cube)[source]

Removes any dimension of length one. If it has an associated DimCoord or AuxCoord, this becomes a scalar coord.

Args:

  • cube (iris.cube.Cube)

    Source cube to remove length 1 dimension(s) from.

Returns

A new iris.cube.Cube instance without any dimensions of length 1.

For example:

>>> cube.shape
(1, 360, 360)
>>> ncube = iris.util.squeeze(cube)
>>> ncube.shape
(360, 360)

↑ top ↑

iris.util.unify_time_units(cubes)[source]

Performs an in-place conversion of the time units of all time coords in the cubes in a given iterable. One common epoch is defined for each calendar found in the cubes to prevent units being defined with inconsistencies between epoch and calendar.

Each epoch is defined from the first suitable time coordinate found in the input cubes.

Arg: