Custom Grid Topology#

This notebook will showcase how to construct UXarray data structures from custom grid topology information, including how to convert existing Xarray data structures to UXarray.

import numpy as np
import xarray as xr

import uxarray as ux

Minimal Grid Definition#

The UGRID conventions require a minimal set of variables for representing a 2D unstructured grid.

Variable

Shape

Description

node_lon

(n_node, )

Longitude of the nodes that make up each face.

node_lat

(n_node, )

Latitude of the nodes that make up each face.

face_node_connectivity

(n_face, n_max_face_nodes])

Indices of the nodes that surround each face.

Fixed Face Sizes#

For grids where each face has the same number of nodes, such as strictly triangular grids, each row in the face_node_connectivity array will contain the indices of nodes that surround each face, with no fill values. Below is an example of the node coordinates and connectivity required for representing a single triangle in the UGRID conventions.

node_lon = np.array([-20.0, 0.0, 20.0])
node_lat = np.array([-10.0, 10.0, -10.0])
face_node_connectivity = np.array(
    [
        [0, 1, 2],
    ]
)

These variables can be passed directly into the Grid.from_topology() class-method, which allows custom grid topology information. This is especially useful for cases where a specific grid format isn’t directly supported.

uxgrid_tri = ux.Grid.from_topology(
    node_lon=node_lon, node_lat=node_lat, face_node_connectivity=face_node_connectivity
)
uxgrid_tri.plot(title="Triangle")

Mixed Face Sizes#

For grids where each face does not have the same number of nodes, the face_node_connectivity array will have it’s final dimension (n_max_face_nodes) set to the largest element shape. For example, a grid with triangles and quadrialterals will have a final dimension of 4. Any element that has less than the maximum number of nodes will be padded with a fill value. The face_node_connectivity array below showcases a basic grid with a single triangle and quadrilateral. Observe that the first row is set to [0, 1, 2, -1], with the first three integers being the indices of the triangle corners, and the final value used as a fill value.

node_lon = np.array([-20.0, 0.0, 20.0, -20, -40])
node_lat = np.array([-10.0, 10.0, -10.0, 10, -10])
face_node_connectivity = np.array([[0, 1, 2, -1], [0, 1, 3, 4]])

The fill_value parameter must be passed in when working with a mixed topology grid.

uxgrid_tri_quad = ux.Grid.from_topology(
    node_lon=node_lon,
    node_lat=node_lat,
    face_node_connectivity=face_node_connectivity,
    fill_value=-1,
)
uxgrid_tri_quad.plot(title="Triangle & Quad")

Working with Existing Xarray Structures#

The previous examples showcase how to create a Grid from custom topology. The follow sections will walk through how to match data to the Grid.

xr.DataArray to ux.UxDataArray#

Consider the previous example where we constructed the uxgrid_tri grid, consisting of a single triangle. Below is an example xr.DataArray containing four temperature values.

xrda_temp = xr.DataArray(
    name="temp",
    data=np.array(
        [
            [100, 105, 108, 109],
        ]
    ).T,
    dims=["time", "cell"],
)
xrda_temp
<xarray.DataArray 'temp' (time: 4, cell: 1)> Size: 32B
array([[100],
       [105],
       [108],
       [109]])
Dimensions without coordinates: time, cell

The original dimension must be mapped to their UGRID equivalent.

Data Mapping

UGRID Dimension Name

Faces

n_face

Edges

n_edge

Nodes

n_node

In the example above, the cell dimension corresponds to the faces of the grid, meaning we want to translate the dimension to n_face

ugrid_dims = {"cell": "n_face"}

The UxDataArray.from_xarray() takes in a user-defined Grid in addition to the original xr.DataArray and UGRID dimension mapping and returns a UxDataArray

uxda_temp = ux.UxDataArray.from_xarray(xrda_temp, uxgrid_tri, ugrid_dims)
uxda_temp
<xarray.UxDataArray 'temp' (time: 4, n_face: 1)> Size: 32B
array([[100],
       [105],
       [108],
       [109]])
Dimensions without coordinates: time, n_face

xr.Dataset to ux.UxDataset#

More commonly, you may have an entire xr.Dataset that contains data variables.

xrda_vorticity = xr.DataArray(
    name="vorticity",
    data=np.array(
        [
            [1, 2, 3, 4],
            [5, 6, 7, 8],
            [9, 10, 11, 12],
        ]
    ).T,
    dims=["time", "vertex"],
)

In this example, we have a cell-centered temp and vertex-centered vorticity variable.

ds = xr.Dataset(data_vars={"temp": xrda_temp, "vorticity": xrda_vorticity})
ds
<xarray.Dataset> Size: 128B
Dimensions:    (time: 4, cell: 1, vertex: 3)
Dimensions without coordinates: time, cell, vertex
Data variables:
    temp       (time, cell) int64 32B 100 105 108 109
    vorticity  (time, vertex) int64 96B 1 5 9 2 6 10 3 7 11 4 8 12

We must now include an additional entry into our ugrid_dims dictionary for our vertex-centered data variable.

ugrid_dims = {"cell": "n_face", "vertex": "n_node"}
uxds = ux.UxDataset.from_xarray(ds=ds, uxgrid=uxgrid_tri, ugrid_dims=ugrid_dims)

uxds now has two UxDataArray objects, one for the cell-centered data temp and one for the vertex-centered data vorticity. There are 4 time steps 0, 1, 2, and 3. Values for time step 0 are shown below.

uxds["vorticity"][0]
<xarray.UxDataArray 'vorticity' (n_node: 3)> Size: 24B
array([1, 5, 9])
Dimensions without coordinates: n_node