iris.analysis.cartography

Various utilities and numeric transformations relevant to cartography.

In this module:

iris.analysis.cartography.area_weights(cube, normalize=False)[source]

Returns an array of area weights, with the same dimensions as the cube.

This is a 2D lat/lon area weights array, repeated over the non lat/lon dimensions.

Args:

Kwargs:

  • normalize (False/True):

    If False, weights are grid cell areas. If True, weights are grid cell areas divided by the total grid area.

The cube must have coordinates ‘latitude’ and ‘longitude’ with bounds.

Area weights are calculated for each lat/lon cell as:

\[r^2 (lon_1 - lon_0) (\sin(lat_1) - \sin(lat_0))\]

Currently, only supports a spherical datum. Uses earth radius from the cube, if present and spherical. Defaults to iris.analysis.cartography.DEFAULT_SPHERICAL_EARTH_RADIUS.

↑ top ↑

iris.analysis.cartography.cosine_latitude_weights(cube)[source]

Returns an array of latitude weights, with the same dimensions as the cube. The weights are the cosine of latitude.

These are n-dimensional latitude weights repeated over the dimensions not covered by the latitude coordinate.

The cube must have a coordinate with ‘latitude’ in the name. Out of range values (greater than 90 degrees or less than -90 degrees) will be clipped to the valid range.

Weights are calculated for each latitude as:

\[w_l = \cos \phi_l\]

Examples:

Compute weights suitable for averaging type operations:

from iris.analysis.cartography import cosine_latitude_weights
cube = iris.load_cube(iris.sample_data_path('air_temp.pp'))
weights = cosine_latitude_weights(cube)

Compute weights suitable for EOF analysis (or other covariance type analyses):

import numpy as np
from iris.analysis.cartography import cosine_latitude_weights
cube = iris.load_cube(iris.sample_data_path('air_temp.pp'))
weights = np.sqrt(cosine_latitude_weights(cube))

↑ top ↑

iris.analysis.cartography.get_xy_contiguous_bounded_grids(cube)[source]

Return 2d arrays for x and y bounds.

Returns array of shape (n+1, m+1).

Example:

xs, ys = get_xy_contiguous_bounded_grids(cube)

↑ top ↑

iris.analysis.cartography.get_xy_grids(cube)[source]

Return 2D X and Y points for a given cube.

Parameters

points. (* cube - The cube for which to generate 2D X and Y) –

Example:

x, y = get_xy_grids(cube)

↑ top ↑

iris.analysis.cartography.gridcell_angles(x, y=None, cell_angle_boundpoints='mid-lhs, mid-rhs')[source]

Calculate gridcell orientations for an arbitrary 2-dimensional grid.

The input grid is defined by two 2-dimensional coordinate arrays with the same dimensions (ny, nx), specifying the geolocations of a 2D mesh.

Input values may be coordinate points (ny, nx) or bounds (ny, nx, 4). However, if points, the edges in the X direction are assumed to be connected by wraparound.

Input can be either two arrays, two coordinates, or a single cube containing two suitable coordinates identified with the ‘x’ and ‘y’ axes.

Args:

The inputs (x [,y]) can be any of the following :

  • x (Cube):

    a grid cube with 2D X and Y coordinates, identified by ‘axis’. The coordinates must be 2-dimensional with the same shape. The two dimensions represent grid dimensions in the order Y, then X.

  • x, y (Coord):

    X and Y coordinates, specifying grid locations on the globe. The coordinates must be 2-dimensional with the same shape. The two dimensions represent grid dimensions in the order Y, then X. If there is no coordinate system, they are assumed to be true longitudes and latitudes. Units must convertible to ‘degrees’.

  • x, y (2-dimensional arrays of same shape (ny, nx)):

    longitude and latitude cell center locations, in degrees. The two dimensions represent grid dimensions in the order Y, then X.

  • x, y (3-dimensional arrays of same shape (ny, nx, 4)):

    longitude and latitude cell bounds, in degrees. The first two dimensions are grid dimensions in the order Y, then X. The last index maps cell corners anticlockwise from bottom-left.

Optional Args:

  • cell_angle_boundpoints (string):

    Controls which gridcell bounds locations are used to calculate angles, if the inputs are bounds or bounded coordinates. Valid values are ‘lower-left, lower-right’, which takes the angle from the lower left to the lower right corner, and ‘mid-lhs, mid-rhs’ which takes an angles between the average of the left-hand and right-hand pairs of corners. The default is ‘mid-lhs, mid-rhs’.

Returns

(2-dimensional cube)

Cube of angles of grid-x vector from true Eastward direction for each gridcell, in degrees. It also has “true” longitude and latitude coordinates, with no coordinate system. When the input has coords, then the output ones are identical if the inputs are true-latlons, otherwise they are transformed true-latlon versions. When the input has bounded coords, then the output coords have matching bounds and centrepoints (possibly transformed). When the input is 2d arrays, or has unbounded coords, then the output coords have matching points and no bounds. When the input is 3d arrays, then the output coords have matching bounds, and the centrepoints are an average of the 4 boundpoints.

Return type

angles

↑ top ↑

iris.analysis.cartography.project(cube, target_proj, nx=None, ny=None)[source]

Nearest neighbour regrid to a specified target projection.

Return a new cube that is the result of projecting a cube with 1 or 2 dimensional latitude-longitude coordinates from its coordinate system into a specified projection e.g. Robinson or Polar Stereographic. This function is intended to be used in cases where the cube’s coordinates prevent one from directly visualising the data, e.g. when the longitude and latitude are two dimensional and do not make up a regular grid.

Parameters
Kwargs:
  • nx

    Desired number of sample points in the x direction for a domain covering the globe.

  • ny

    Desired number of sample points in the y direction for a domain covering the globe.

Returns

An instance of iris.cube.Cube and a list describing the extent of the projection.

Note

This function assumes global data and will if necessary extrapolate beyond the geographical extent of the source cube using a nearest neighbour approach. nx and ny then include those points which are outside of the target projection.

Note

Masked arrays are handled by passing their masked status to the resulting nearest neighbour values. If masked, the value in the resulting cube is set to 0.

Warning

This function uses a nearest neighbour approach rather than any form of linear/non-linear interpolation to determine the data value of each cell in the resulting cube. Consequently it may have an adverse effect on the statistics of the data e.g. the mean and standard deviation will not be preserved.

Warning

If the target projection is non-rectangular, e.g. Robinson, the target grid may include points outside the boundary of the projection. The latitude/longitude of such points may be unpredictable.

↑ top ↑

iris.analysis.cartography.rotate_grid_vectors(u_cube, v_cube, grid_angles_cube=None, grid_angles_kwargs=None)[source]

Rotate distance vectors from grid-oriented to true-latlon-oriented.

Can also rotate by arbitrary angles, if they are passed in.

Note

This operation overlaps somewhat in function with iris.analysis.cartography.rotate_winds(). However, that routine only rotates vectors according to transformations between coordinate systems. This function, by contrast, can rotate vectors by arbitrary angles. Most commonly, the angles are estimated solely from grid sampling points, using gridcell_angles() : This allows operation on complex meshes defined by two-dimensional coordinates, such as most ocean grids.

Args:

  • u_cube, v_cube(cube)

    Cubes of grid-u and grid-v vector components. Units should be differentials of true-distance, e.g. ‘m/s’.

Optional args:

  • grid_angles_cube(cube)

    gridcell orientation angles. Units must be angular, i.e. can be converted to ‘radians’. If not provided, grid angles are estimated from ‘u_cube’ using the gridcell_angles() method.

  • grid_angles_kwargs(dict or None)

    Additional keyword args to be passed to the gridcell_angles() method, if it is used.

Returns

(cube)

Cubes of true-north oriented vector components. Units are same as inputs.

Note

Vector magnitudes will always be the same as the inputs.

Return type

true_u, true_v

↑ top ↑

iris.analysis.cartography.rotate_pole(lons, lats, pole_lon, pole_lat)[source]

Convert arrays of longitudes and latitudes to arrays of rotated-pole longitudes and latitudes. The values of pole_lon and pole_lat should describe the rotated pole that the arrays of longitudes and latitudes are to be rotated onto.

As the arrays of longitudes and latitudes must describe a rectilinear grid, the arrays of rotated-pole longitudes and latitudes must be of the same shape as each other.

Example:

rotated_lons, rotated_lats = rotate_pole(lons, lats,         pole_lon, pole_lat)

Note

Uses proj.4 to perform the conversion.

Parameters
  • lons (*) – An array of longitude values.

  • lats (*) – An array of latitude values.

  • pole_lon (*) – The longitude of the rotated pole that the arrays of longitudes and latitudes are to be rotated onto.

  • pole_lat (*) – The latitude of the rotated pole that the arrays of longitudes and latitudes are to be rotated onto.

Returns

An array of rotated-pole longitudes and an array of rotated-pole latitudes.

↑ top ↑

iris.analysis.cartography.rotate_winds(u_cube, v_cube, target_cs)[source]

Transform wind vectors to a different coordinate system.

The input cubes contain U and V components parallel to the local X and Y directions of the input grid at each point.

The output cubes contain the same winds, at the same locations, but relative to the grid directions of a different coordinate system. Thus in vector terms, the magnitudes will always be the same, but the angles can be different.

The outputs retain the original horizontal dimension coordinates, but also have two 2-dimensional auxiliary coordinates containing the X and Y locations in the target coordinate system.

Args:

Returns

A (u’, v’) tuple of iris.cube.Cube instances that are the u and v components in the requested target coordinate system. The units are the same as the inputs.

Note

The U and V values relate to distance, with units such as ‘m s-1’. These are not the same as coordinate vectors, which transform in a different manner.

Note

The names of the output cubes are those of the inputs, prefixed with ‘transformed_’ (e.g. ‘transformed_x_wind’).

Warning

Conversion between rotated-pole and non-rotated systems can be expressed analytically. However, this function always uses a numerical approach. In locations where this numerical approach does not preserve magnitude to an accuracy of 0.1%, the corresponding elements of the returned cubes will be masked.

↑ top ↑

iris.analysis.cartography.unrotate_pole(rotated_lons, rotated_lats, pole_lon, pole_lat)[source]

Convert arrays of rotated-pole longitudes and latitudes to unrotated arrays of longitudes and latitudes. The values of pole_lon and pole_lat should describe the location of the rotated pole that describes the arrays of rotated-pole longitudes and latitudes.

As the arrays of rotated-pole longitudes and latitudes must describe a rectilinear grid, the arrays of rotated-pole longitudes and latitudes must be of the same shape as each other.

Example:

lons, lats = unrotate_pole(rotated_lons, rotated_lats,           pole_lon, pole_lat)

Note

Uses proj.4 to perform the conversion.

Parameters
  • rotated_lons (*) – An array of rotated-pole longitude values.

  • rotated_lats (*) – An array of rotated-pole latitude values.

  • pole_lon (*) – The longitude of the rotated pole that describes the arrays of rotated-pole longitudes and latitudes.

  • pole_lat (*) – The latitude of the rotated pole that describes the arrays of rotated-pole longitudes and latitudes.

Returns

An array of unrotated longitudes and an array of unrotated latitudes.

↑ top ↑

iris.analysis.cartography.wrap_lons(lons, base, period)[source]

Wrap longitude values into the range between base and base+period.

For example:
>>> print(wrap_lons(np.array([185, 30, -200, 75]), -180, 360))
[-175.   30.  160.   75.]

↑ top ↑

DistanceDifferential(dx1, dy1, dx2, dy2)

class iris.analysis.cartography.DistanceDifferential(_cls, dx1, dy1, dx2, dy2)

Create new instance of DistanceDifferential(dx1, dy1, dx2, dy2)

count(value, /)

Return number of occurrences of value.

index(value, start=0, stop=9223372036854775807, /)

Return first index of value.

Raises ValueError if the value is not present.

dx1

Alias for field number 0

dx2

Alias for field number 2

dy1

Alias for field number 1

dy2

Alias for field number 3

↑ top ↑

PartialDifferential(dx1, dy1)

class iris.analysis.cartography.PartialDifferential(_cls, dx1, dy1)

Create new instance of PartialDifferential(dx1, dy1)

count(value, /)

Return number of occurrences of value.

index(value, start=0, stop=9223372036854775807, /)

Return first index of value.

Raises ValueError if the value is not present.

dx1

Alias for field number 0

dy1

Alias for field number 1