Convention interface#

Each supported convention implements the Convention interface.

open_dataset(path, **kwargs)#

Open a dataset and determine the correct Convention implementation for it. If a valid Convention implementation can not be found, an error is raised.

Parameters:
  • path – The path to the dataset to open

  • kwargs – These are passed straight through to xarray.open_dataset().

Returns:

xarray.Dataset – The opened dataset

Example

import emsarray
dataset = emsarray.open_dataset("./tests/datasets/ugrid_mesh2d.nc")
get_dataset_convention(dataset)#

Find the most appropriate Convention subclass for this dataset.

Parameters:

dataset – The dataset to introspect

Returns:

Convention – A Convention subclass appropriate for this dataset, or None if nothing appropriate can be found.

Example

>>> ds = emsarray.tutorial.open_dataset('austen.nc')
>>> get_dataset_convention(ds)
emsarray.conventions.ugrid.UGrid
class Convention(dataset)#

Each supported geometry convention represents data differently. The Convention class abstracts these differences away, allowing developers to query and manipulate datasets without worrying about the details. See Supported dataset conventions for a list of implemented conventions.

All conventions have the concept of a cell at a geographic location, vertically stacked layers of cells, and multiple timesteps of data. A convention may support additional grids, such as face edges and vertices. Refer to Grids for more information.

A cell can be addressed using a linear index or a native index. A linear index is always an int, while the native index type will depend on the specific convention. You can convert between a linear and a native index using ravel_index() and wind_index(). Refer to Indexing for more information.

dataset#

The xarray.Dataset instance for this Convention

classmethod check_validity(dataset)#

Checks that the dataset is OK to use. Called during __init__, and raises an exception if the dataset has problems.

abstract classmethod check_dataset(dataset)#

Check if a dataset uses this convention.

This may check for variables of the correct dimensions, the presence of specific attributes, or the layout of the dataset dimensions.

Specific subclasses must implement this function. It will be called by the convention autodetector when guessing the correct convention for a dataset.

If the dataset matches, the return value indicates how specific the match is. When autodetecting the correct convention implementation the convention with the highest specicifity will be used. Many conventions extend the CF grid conventions, so the CF Grid convention classes will match many datasets. However this match is very generic. A more specific implementation such as SHOC may be supported. The SHOC convention implementation should return a higher specicifity than the CF grid convention.

Parameters:

dataset (xarray.Dataset) – The dataset instance to inspect.

Returns:

int, optional – If this convention implementation can handle this dataset some integer greater than zero is returned. The higher the number, the more specific the support. If the dataset does not match this convention, None is returned. Values on the Specificity enum are used by emsarray itself to indicated specificity. New convention implementations are free to use these values, or use any integer value.

Example

>>> import xarray
... from emsarray.conventions.shoc import ShocStandard
... from emsarray.conventions.ugrid import UGrid
... dataset = xarray.open_dataset("./tests/datasets/shoc_standard.nc")
>>> ShocStandard.check_dataset(dataset)
True
>>> UGrid.check_dataset(dataset)
False
bind()#

Bind this Convention instance as the default convention for the xarray.Dataset. This convention instance will be assigned to dataset.ems.

You can use a Convention instance without binding it to a Dataset, binding is only necessary if you need to use the dataset.ems accessor.

Note

If you use emsarray.open_dataset() or dataset.ems to autodetect the dataset convention you do not need to call this method. Convention.bind() is only useful if you manually construct a Convention.

Examples

# Open a dataset built using the GRASS convention
dataset = xarray.open_dataset("grass-dataset.nc")

# Construct a Grass instance for the dataset and bind it
convention = Grass(dataset)
convention.bind()

# dataset.ems is now the bound convention
assert dataset.ems is convention

If the dataset already has a bound convention, an error is raised. To bind a new convention to a dataset, make a copy of the dataset first:

new_dataset = dataset.copy()
convention = Grass(new_dataset)
convention.bind()
property time_coordinate#

The time coordinate for this dataset.

Returns:

xarray.DataArray – The variable for the time coordinate for this dataset.

Raises:

exceptions.NoSuchCoordinateError – If no time coordinate was found

See also

get_time_name

property depth_coordinate#

The depth coordinate for this dataset.

Returns:

xarray.DataArray – The variable for the depth coordinate for this dataset.

Raises:

exceptions.NoSuchCoordinateError – If no depth coordinate was found

See also

get_depth_name

get_time_name()#

Get the name of the time variable in this dataset.

Returns:

Hashable – The name of the time coordinate.

Raises:

exceptions.NoSuchCoordinateError – If no time coordinate was found

Notes

The CF Conventions state that a time variable is defined by having a units attribute formatted according to the UDUNITS package [1].

xarray will find all time variables and convert them to numpy datetimes.

References

get_depth_name()#

Get the name of the layer depth coordinate variable.

For datasets with multiple depth variables, this should be the one that represents the centre of the layer, not the bounds.

Note that this is the name of the coordinate variable, not the name of the dimension, for datasets where these differ.

Returns:

Hashable – The name of the depth coordinate.

Raises:

exceptions.NoSuchCoordinateError – If no time coordinate was found

Notes

The CF Conventions state that a depth variable is identifiable by units of pressure; or the presence of the positive attribute with value of up or down [2].

In practice, many datasets do not follow this convention. In addition to checking for the positive attribute, all coordinates are checked for a standard_name: "depth", coordinate_type: "Z", or axiz: "Z".

References

get_all_depth_names()#

Get the names of all depth layers. Some datasets include both a depth layer centre, and the depth layer ‘edges’.

Note that this is the names of the coordinate variables, not the names of the dimensions, for datasets where these differ.

get_depths()#

Get the depth of each vertical layer in this dataset.

Deprecated since version 0.5.0: This method is replaced by Convention.depth_coordinate.values.

Returns:

numpy.ndarray – An array of depths, one per vertical layer in the dataset.

get_times()#

Get all timesteps in this dataset.

Deprecated since version 0.5.0: This method is replaced by Convention.time_coordinate.values.

Returns:

numpy.ndarray – An array of datetimes. The datetimes will be whatever native format the dataset uses, likely numpy.datetime64.

abstract ravel_index(index)#

Convert a convention native index to a linear index.

Each conventnion has a different native index type, read the specific conventnion documentation for more information.

Parameters:

index (Index) – The native index to ravel.

Returns:

int – The linear index for that native index.

Example

If the dataset used the CF Grid conventions, across a (latitude, longitude) grid of size (30, 40):

>>> temp = dataset.data_vars['temp']
>>> temp.dims
('t', 'z', 'y', 'x')
>>> temp.shape
(10, 20, 30, 40)
>>> dataset.ems.ravel_index((3, 4))
124

Cell polygons are indexed in the same order as the linear indices for cells. To find the polygon for the cell with the native index (3, 4):

>>> index = (3, 4)
>>> linear_index = dataset.ems.ravel_index(index)
>>> polygon = dataset.ems.polygons[linear_index]

See also

Convention.wind_index()

The inverse operation

abstract wind_index(linear_index, *, grid_kind=None)#

Convert a linear index to a conventnion native index.

Each conventnion has a different native index type, read the specific conventnion documentation for more information.

Parameters:
  • linear_index (int) – The linear index to wind.

  • grid_kind (GridKind, optional) – Used to indicate what kind of index is being wound, for conventions with multiple grids. Optional, if not provided the default grid kind will be used.

Returns:

Index – The convention native index for that same cell

Example

If the dataset used the CF Grid conventions, across a (latitude, longitude) grid of size (30, 40):

>>> temp = dataset.data_vars['temp']
>>> temp.dims
('t', 'z', 'y', 'x')
>>> temp.shape
(10, 20, 30, 40)
>>> dataset.ems.wind_index(124)
(3, 4)

See also

Convention.ravel_index()

The inverse operation

unravel_index(linear_index, grid_kind=None)#

An alias for Convention.wind_index().

Deprecated since version 0.6.0: Use Convention.wind_index() instead

abstract property grid_kinds#

All of the grid kinds this dataset includes.

abstract property default_grid_kind#

The default grid kind for this dataset. For most datasets this should be the face grid.

abstract property grid_size#

The linear size of each grid kind.

abstract get_grid_kind(data_array)#

Determines the relevant grid kind for this data array.

If the data array is not indexable using the native index types a ValueError is raised.

Parameters:

data_array (xarray.DataArray) – The data array to inspect

Returns:

GridKind

Raises:

ValueError – If the data array passed in is not indexable using any native index type a ValueError is raised. Depth coordinates or time coordinates are examples of data arrays that will not be indexable and will raise an error.

Example

For a UGRID dataset with temperature data defined at the cell centres and current defined as flux through the cell edges:

>>> dataset.data_vars['temp'].dims
('time', 'depth', 'face')
>>> dataset.data_vars['u1'].dims
('time', 'depth', 'edge')
>>> dataset.ems.get_grid_kind(dataset.data_vars['temp'])
UGridKind.face
>>> dataset.ems.get_grid_kind(dataset.data_vars['u1'])
UGridKind.edge
get_grid_kind_and_size(data_array)#

Determines the relevant index kind and the extent of the linear index space for this data array.

Deprecated since version 0.6.0: This method is replaced by Convention.get_grid_kind() and Convention.grid_size.

If the data array is not indexable using the native index types a ValueError is raised.

Parameters:

data_array (xarray.DataArray) – The data array to introspect

Returns:

tuple of GridKind and int

Raises:

ValueError – If the data array passed in is not indexable using any native index type a ValueError is raised. Depth coordinates or time coordinates are examples of data arrays that will not be indexable and will raise an error.

Example

For a UGRID dataset with temperature data defined at the cell centres and current defined as flux through the cell edges:

>>> dataset.ems.topology.face_count
4
>>> dataset.ems.topology.edge_count
9
>>> dataset.data_vars['temp'].dims
('time', 'depth', 'face')
>>> dataset.data_vars['u1'].dims
('time', 'depth', 'edge')
>>> dataset.ems.get_grid_kind_and_size(dataset.data_vars['temp'])
(UGridKind.face, 4)
>>> dataset.ems.get_grid_kind_and_size(dataset.data_vars['u1'])
(UGridKind.edge, 9)
abstract ravel(data_array, *, linear_dimension=None)#

Flatten the surface dimensions of a DataArray, returning a flatter numpy.ndarray indexed in the same order as the linear index.

For DataArrays with extra dimensions such as time or depth, only the surface dimensions are flattened. Other dimensions are left as is.

For datasets with multiple grids, with data defined on edges or vertices for example, this will flatten those data arrays in the correct linear order to be indexed by the relevant index type.

Parameters:
  • data_array (xarray.DataArray) – One of the data variables from this dataset.

  • linear_dimension (Hashable, optional) – The name of the new dimension to flatten the surface dimensions to. Defaults to ‘index’.

Returns:

xarray.DataArray – A new data array, where all the surface dimensions have been flattened in to one linear array. The values for each cell, in the same order as the linear index for this dataset. Any other dimensions, such as depth or time, will be retained.

See also

utils.ravel_dimensions

A function that ravels some given dimensions in a dataset.

Convention.wind

The inverse operation.

abstract wind(data_array, *, grid_kind=None, axis=None, linear_dimension=None)#

Wind a flattened DataArray so that it has the same shape as data variables in this dataset.

By using grid_size and wind() together it is possible to construct new data variables for a dataset of any arbitrary shape.

Parameters:
  • data_array (xarray.DataArray) – One of the data variables from this dataset.

  • grid_kind (GridKind) – The kind of grid this data array represents, for those conventions with multiple grid kinds. Optional, defaults to the default grid kind.

  • axis (int, optional) – The axis number that should be wound. Optional, defaults to the last axis.

  • linear_dimension (Hashable, optional) – The axis number that should be wound. Optional, defaults to the last dimension.

Returns:

xarray.DataArray – A new data array where the linear data have been wound to match the shape of the convention. Any other dimensions, such as depth or time, will be retained.

Examples

The following will construct a data array of the correct shape for any convention supported by emsarray:

import emsarray
import numpy
import xarray

dataset = emsarray.tutorial.open_dataset('fraser')
face_size = dataset.ems.grid_size[dataset.ems.default_grid_kind]
flat_array = xarray.DataArray(
    data=numpy.arange(face_size),
    dims=['index'],
)
data_array = dataset.ems.wind(flat_array)

This will construct a boolean array indicating which cells of a dataset intersect a target geometry:

import emsarray
import numpy
import shapely
import xarray

dataset = emsarray.tutorial.open_dataset('gbr4')
target = shapely.Polygon([
    [152.8088379, -22.7863108],
    [153.9184570, -22.2280904],
    [153.4680176, -20.9614396],
    [151.8255615, -20.4012720],
    [151.4135742, -21.8309067],
    [152.0068359, -22.4313402],
    [152.8088379, -22.7863108],
])

hits = dataset.ems.strtree.query(target, predicate='intersects')
grid_size = dataset.ems.grid_size[dataset.ems.default_grid_kind]
intersecting_cells = xarray.DataArray(
    data=numpy.zeros(grid_size, dtype=bool),
    dims=['index'],
)
intersecting_cells.values[hits] = True
intersecting_cells = dataset.ems.wind(intersecting_cells)

See also

utils.wind_dimension

Reshape a particular dimension in a data array.

Convention.ravel

The inverse operation.

make_linear(data_array)#

A deprecated alias for Convention.ravel()

Deprecated since version 0.6.0: This method is replaced by Convention.ravel()

property data_crs#

The coordinate reference system that coordinates in this dataset are defined in. Used by Convention.make_poly_collection() and Convention.make_quiver(). Defaults to cartopy.crs.PlateCarree.

plot_on_figure(figure, scalar=None, vector=None, title=None, **kwargs)#

Plot values for a DataArray on a matplotlib Figure.

The data array can either be passed in directly, or the name of a data array on this Convention.dataset instance. The data array does not have to come from the same dataset, as long as the dimensions are the same.

This method will only plot a single time step and depth layer. Callers are responsible for selecting a single slice before calling this method.

Parameters:

See also

plot.plot_on_figure()

The underlying implementation

plot(*args, **kwargs)#

Plot a data array and automatically display it.

This method is most useful when working in Jupyter notebooks which display figures automatically. This method is a wrapper around plot_on_figure() that creates and shows a Figure for you. All arguments are passed on to plot_on_figure(), refer to that function for details.

See also

plot_on_figure()

animate_on_figure(figure, scalar=None, vector=None, coordinate=None, title=None, **kwargs)#

Make an animated plot of a data array.

For real world examples, refer to the examples/animation.ipynb notebook.

Parameters:
Returns:

matplotlib.animation.Animation – The data variable plotted as an animation on the figure. Call Animation.save() to save animation to a file, or display it in a Notebook using Animation.to_jshtml().

make_poly_collection(data_array=None, **kwargs)#

Make a PolyCollection from the geometry of this Dataset. This can be used to make custom matplotlib plots from your data.

If a DataArray is passed in, the values of that are assigned to the PolyCollection array parameter.

Parameters:
  • data_array (Hashable or xarray.DataArray, optional) – A data array, or the name of a data variable in this dataset. Optional. If given, the data array is ravelled and passed to PolyCollection.set_array(). The data is used to colour the patches. Refer to the matplotlib documentation for more information on styling.

  • **kwargs – Any keyword arguments are passed to the PolyCollection constructor.

Returns:

PolyCollection – A PolyCollection constructed using the geometry of this dataset.

Example

import cartopy.crs as ccrs
import matplotlib.pyplot as plt
import emsarray

figure = plt.figure(figsize=(10, 8))
axes = plt.subplot(projection=ccrs.PlateCarree())
axes.set_aspect(aspect='equal', adjustable='datalim')

ds = emsarray.open_dataset("./tests/datasets/ugrid_mesh2d.nc")
ds = ds.isel(record=0, Mesh2_layers=-1)
patches = ds.ems.make_poly_collection('temp')
axes.add_collection(patches)
figure.colorbar(patches, ax=axes, location='right', label='meters')

axes.set_title("Depth")
axes.autoscale()
figure.show()
make_quiver(axes, u=None, v=None, **kwargs)#

Make a matplotlib.quiver.Quiver instance to plot vector data.

Parameters:
  • axes (matplotlib.axes.Axes) – The axes to make this quiver on.

  • u, v (xarray.DataArray or str, optional) – The DataArrays or the names of DataArrays in this dataset that make up the u and v components of the vector. If omitted, a Quiver will be constructed with all components set to 0.

  • **kwargs – Any keyword arguments are passed on to the Quiver constructor.

Returns:

matplotlib.quiver.Quiver – A quiver instance that can be added to a plot

abstract property polygons#

A numpy.ndarray of shapely.Polygon instances representing the cells in this dataset.

The order of the polygons in the list corresponds to the linear index of this dataset. Not all valid cell indices have a polygon, these holes are represented as None in the list. If you want a list of just polygons, apply the mask:

dataset = emsarray.open_dataset("...")
only_polygons = dataset.ems.polygons[dataset.ems.mask]

See also

ravel_index(), mask

property face_centres#

A numpy ndarray of face centres, which are (x, y) pairs. The first dimension will be the same length and in the same order as Convention.polygons, while the second dimension will always be of size 2.

property mask#

A boolean numpy.ndarray indicating which cells have valid polygons. This can be used to select only items from linear arrays that have a corresponding polygon.

dataset = emsarray.open_dataset("...")
mask = dataset.ems.mask
plottable_polygons = dataset.ems.polygons[mask]
plottable_values = dataset.ems.ravel("eta")[mask]
property geometry#

A shapely.Polygon or shapely.MultiPolygon that represents the geometry of the entire dataset.

This is equivalent to the union of all polygons in the dataset, although specific conventions may have a simpler way of constructing this.

property bounds#

Returns minimum bounding region (minx, miny, maxx, maxy) of the entire dataset.

This is equivalent to the bounds of the dataset geometry, although specific conventons may have a simpler way of constructing this.

property strtree#

A shapely.strtree.STRtree spatial index of all cells in this dataset. This allows for fast spatial lookups, querying which cells lie at a point, or which cells intersect a geometry.

Querying the STRtree will return the linear indices of any matching cells. Use polygons to find the geometries associated with each index. Use wind_index() to transform this back to a native index, or ravel() to linearise a variable.

Examples

Find the indices of all cells that intersect a line:

dataset = emsarray.tutorial.open_dataset('austen')
geometry = shapely.linestring([
    [146.9311523, -15.7076628],
    [149.3261719, -19.0413488],
    [152.2485352, -21.6778483],
    [154.2480469, -24.2469646],
    [155.1049805, -27.5082714],
    [154.7753906, -31.2973280],
    [153.3911133, -34.1072564],
    [152.0947266, -36.0846213],
    [150.6005859, -38.7712164],
])
hits = dataset.ems.strtree.query(geometry, predicate='intersects')
property spatial_index#

A shapely.strtree.STRtree spatial index of all cells in this dataset. This allows for fast spatial lookups, querying which cells lie at a point, or which cells intersect a geometry.

Deprecated since version 0.6.0: Use Convention.strtree instead.

This existed as a wrapper around a Shapely STRtree, which changed its interface in Shapely 2.0. Shapely 1.8.x is no longer supported by emsarray so this compatibility wrapper is deprecated. Use Convention.strtree directly instead.

Querying this spatial index will return a list of (polygon, SpatialIndexItem) tuples corresponding to each matching cell. SpatialIndexItem instances have the cells linear index, native index, and polygon.

The query results from the STRtree contain all geometries with overlapping bounding boxes. Query results need to be refined further by comparing the cell geometry to the query geometry. Refer to the Shapely 1.8.x STRtree docs for examples.

See also

SpatialIndexItem : The dataclass returned from querying the STRtree.

Shapely 1.8.x STRtree docs

get_index_for_point(point)#

Find the index for a Point in the dataset.

Parameters:

point (shapely.Point) – The geographic point to query

Returns:

SpatialIndexItem, optional – The SpatialIndexItem for the point queried. This indicates the polygon that intersected the point and the index of that polygon in the dataset.

If the point does not intersect the dataset, None is returned.

Notes

In the case where the point intersects multiple cells the cell with the lowest linear index is returned. This can happen if the point is exactly one of the cell vertices, or falls on a cell edge, or if the geometry of the dataset contains overlapping polygons.

abstract selector_for_index(index)#

Convert a convention native index into a selector that can be passed to Dataset.isel.

Parameters:

index (Index) – A convention native index

Returns:

selector – A dict suitable for passing to xarray.Dataset.isel() that will select values at this index.

select_index(index)#

Return a new dataset that contains values only from a single index. This is much like doing a xarray.Dataset.isel() on an index, but works with convention native index types.

An index is associated with a grid kind. The returned dataset will only contain variables that were defined on this grid, with the single indexed point selected. For example, if the index of a face is passed in, the returned dataset will not contain any variables defined on an edge.

Parameters:

index (Index) – The index to select. The index must be for the default grid kind for this dataset.

Returns:

xarray.Dataset – A new dataset that is subset to the one index.

Notes

The returned dataset will most likely not have sufficient coordinate data to be used with a particular Convention any more. The dataset.ems accessor will raise an error if accessed on the new dataset.

select_point(point)#

Return a new dataset that contains values for a single point. This is a shortcut for get_index_for_point() and select_index().

If the point is not inside the dataset a ValueError is raised.

Parameters:

point (shapely.Point) – The point to select

Returns:

xarray.Dataset – A dataset of values at the point

abstract get_all_geometry_names()#

Return a list of the names of all geometry variables used by this convention.

drop_geometry()#

Return a new xarray.Dataset with all geometry variables dropped. Useful when significantly transforming the dataset, such as extracting point data.

select_variables(variables)#

Select only a subset of the variables in this dataset, dropping all others.

This will keep all coordinate variables and all geometry variables.

Parameters:

variables (iterable of Hashable) – The names of all data variables to select.

Returns:

xarray.DataArray – A new dataset with the same geometry and coordinates, but only the selected data variables.

abstract make_clip_mask(clip_geometry, *, buffer=0)#

Make a new Dataset that can be used to clip this dataset to only the cells that intersect some geometry.

This dataset can be saved to a file to be reused to cut multiple datasets with identical shapes, such as a series of files representing multiple time series of a model.

The mask can be applied to this dataset (or other datasets identical in shape) using apply_clip_mask().

Parameters:
  • clip_geometry (shapely.BaseGeometry) – The desired area to cut out. This can be any shapely geometry type, but will most likely be a polygon

  • buffer (int, optional) – If set to a positive integer, a buffer of that many cells will be added around the clip region. This is useful if you need to clip to a particular area, but also would like to do some interpolation on the output cells.

Returns:

xarray.Dataset – The mask

abstract apply_clip_mask(clip_mask, work_dir)#

Apply a clip mask to this dataset, and return a new dataset. Call make_clip_mask() to create a clip mask from a clip geometry.

The clip_mask can be saved and loaded to disk if the mask needs to be reused across multiple datasets, such as multiple time series from one model.

Depending on the implementation, the input dataset may be sliced in to multiple files during cutting, and the returned Dataset might be a multi-file Dataset built from these temporary files. The caller must either load the dataset in to memory using load() or compute(), or save the dataset to disk somewhere outside of the working directory before the working directory is cleaned up.

Parameters:
Returns:

xarray.Dataset – A new Dataset clipped using the mask

clip(clip_geomery, work_dir, *, buffer=0)#

Generates a clip mask and applies it in one step.

See the documentation for make_clip_mask() and apply_clip_mask() for more details.

Parameters:
  • clip_geometry (shapely.BaseGeometry) – The desired area to cut out. This can be any shapely geometry type, but will most likely be a polygon

  • work_dir (str or pathlib.Path) – A directory where temporary files can be written to. Callers must create and manage this temporary directory, perhaps using tempfile.TemporaryDirectory.

  • buffer (int, optional) – If set to a positive integer, a buffer of that many cells will be added around the clip region. This is useful if you need to clip to a particular area, but also would like to do some interpolation on the output cells.

Returns:

xarray.Dataset – A new Dataset clipped using the mask

to_netcdf(path, **kwargs)#

Save this dataset to a netCDF file, and also fix up the time units to make the EMS compatible.

ocean_floor()#

An alias for emsarray.operations.depth.ocean_floor()

normalize_depth_variables(*, positive_down=None, deep_to_shallow=None)#

An alias for emsarray.operations.depth.normalize_depth_variables()

class DimensionConvention(dataset)#

A Convention subclass where different grid kinds are always defined on unique sets of dimension. This covers most conventions.

This subclass adds the abstract methods and properties:

Default implementations are provided for:

grid_dimensions#

The dimensions associated with a particular grid kind.

This is a mapping between grid kinds and an ordered list of dimension names. Each dimension in the dataset must be associated with at most one grid kind. Each grid kind must be associated with at least one dimension. The dimensions must be in the order expected in a dataset, if order is significant.

This property may introspect the dataset to determine which dimensions are used. The property should be cached.

abstract unpack_index(index)#

Convert a native index in to a grid kind and dimension indices.

Parameters:

index (Index) – A native index

Returns:

  • grid_kind (GridKind) – The grid kind

  • indices (sequence of int) – The dimension indices

See also

pack_index

abstract pack_index(grid_kind, indices)#

Convert a grid kind and dimension indices in to a native index.

Parameters:
  • grid_kind (GridKind) – The grid kind

  • indices (sequence of int) – The dimension indices

Returns:

index (Index) – The corresponding native index

See also

unpack_index

class SpatialIndexItem(linear_index, index, polygon)#

Information about an item in the STRtree spatial index for a dataset.

linear_index#

The linear index of this cell

index#

The native index of this cell

polygon#

The geographic shape of this cell

GridKind = ~GridKind#

Some type that can enumerate the different grid types present in a dataset. This can be an enum.Enum listing each different kind of grid.

Index values will be included in the feature properties of exported geometry from emsarray.operations.geometry. If the index type includes the grid kind, the grid kind needs to be JSON serializable. The easiest way to achieve this is to make your GridKind type subclass str:

class MyGridKind(str, enum.Enum):
    face = 'face'
    edge = 'edge'
    node = 'node'

For cases where the convention only supports a single grid, a singleton enum can be used.

More esoteric cases involving datasets with a potentially unbounded numbers of grids can use a type that supports this instead.

Index = ~Index#

An index to a specific point on a grid in this convention. For conventions with multiple grids (e.g. cells, edges, and nodes), this should be a tuple whos first element is GridKind. For conventions with a single grid, GridKind is not required.

class Specificity(value)#

How specific a match is when autodetecting a convention. Matches with higher specificity will be prioritised.

General conventions such as CF Grid are low specificity, as many conventions extend and build on CF Grid conventions.

The SHOC conventions extend the CF grid conventions, so a SHOC file will be detected as both CF Grid and SHOC. ShocStandard should return a higher specificity so that the correct convention implementation is used.

LOW = 10#
MEDIUM = 20#
HIGH = 30#