Data Structures


UXarray extends upon Xarray’s core data structures (i.e. Dataset and DataArray) and provides functionality for operating on unstructured (a.k.a. non-regular) grids.

This notebook will showcase the core UXarray data structures and how to interact with them

import uxarray as ux

First, lets specify the paths to our Grid and Data files. As mentioned in the previous notebook, unstructured grids are typically separated into two files: one containing the grid topology and one containing data variables that reside on that grid.

file_dir = "../../meshfiles/"
grid_filename = file_dir + "oQU480.grid.nc"
data_filename = file_dir + "oQU480.data.nc"

Grid Data Structure

An unstructured grid file can be opened standalone using the ux.open_grid method, which is a catch-all method that parses the provided input and returns a Grid object with grid coordinate and connectivity variables represented in the UGRID conventions.

Printing our Grid instance shows all available grid variables and original grid type.

grid = ux.open_grid(grid_filename)
grid
<uxarray.Grid>
Original Grid Type: MPAS
Grid Dimensions:
  * n_node: 3947
  * n_edge: 5754
  * n_face: 1791
  * n_max_face_nodes: 6
  * n_max_face_edges: 6
  * n_max_face_faces: 6
  * n_max_node_faces: 3
  * two: 2
  * n_nodes_per_face: (1791,)
Grid Coordinates (Spherical):
  * node_lon: (3947,)
  * node_lat: (3947,)
  * edge_lon: (5754,)
  * edge_lat: (5754,)
  * face_lon: (1791,)
  * face_lat: (1791,)
Grid Coordinates (Cartesian):
  * node_x: (3947,)
  * node_y: (3947,)
  * node_z: (3947,)
  * edge_x: (5754,)
  * edge_y: (5754,)
  * edge_z: (5754,)
  * face_x: (1791,)
  * face_y: (1791,)
  * face_z: (1791,)
Grid Connectivity Variables:
  * face_node_connectivity: (1791, 6)
  * face_edge_connectivity: (1791, 6)
  * face_face_connectivity: (1791, 6)
  * edge_node_connectivity: (5754, 2)
  * edge_face_connectivity: (5754, 2)
  * node_face_connectivity: (3947, 3)
Grid Descriptor Variables:
  * face_areas: (1791,)
  * n_nodes_per_face: (1791,)
  * edge_face_distances: (5754,)
  * edge_node_distances: (5754,)

We can access our coordinate, connectivity, and other descriptor variables as attributes through this Grid instance

grid.node_lon
<xarray.DataArray 'node_lon' (n_node: 3947)> Size: 32kB
array([-174.95294523, -172.5318284 , -173.42204784, ...,  113.04705503,
       -102.95294492, -174.9529447 ])
Dimensions without coordinates: n_node
Attributes:
    standard_name:  longitude
    long name:      Longitude of the corner nodes of each face
    units:          degrees_east
grid.face_node_connectivity
<xarray.DataArray 'face_node_connectivity' (n_face: 1791, n_max_face_nodes: 6)> Size: 86kB
array([[                   2,                    3,                    4,
                           0,                    1, -9223372036854775808],
       [                   6,                    7,                    8,
                           9,                    5, -9223372036854775808],
       [                  14,                   10,                   11,
                          12,                   13, -9223372036854775808],
       ...,
       [                2943,                 2948,                 3695,
                        2955,                 2960,                 3946],
       [                2950,                 2949,                 3686,
                        2944,                 2943,                 3946],
       [                3946,                 2960,                 2959,
                        3689,                 2951,                 2950]])
Dimensions without coordinates: n_face, n_max_face_nodes
Attributes:
    cf_role:      face_node_connectivity
    long name:    Maps every face to its corner nodes.
    start_index:  0
    _FillValue:   -9223372036854775808

UxDataset & UxDataArray Data Structures

UXarray inherits from Xarray’s two core data structures Dataset & DataArray to provide a grid-informed implementation through the UxDataset & UxDataArray data structures.

The major difference between them is that UXarray’s implementation is paired with a Grid object, accessed through the .uxgrid property.

UXarray also provides an overloaded ux.open_dataset method, which takes in both a Grid and Data file path to construct a UxDataset

uxds = ux.open_dataset(grid_filename, data_filename)
uxds
<xarray.UxDataset> Size: 14kB
Dimensions:      (n_face: 1791)
Dimensions without coordinates: n_face
Data variables:
    bottomDepth  (n_face) float64 14kB ...

We can see that our UxDataset has a single data variable bottomDepth, which is mapped to each face (as specified by the n_face dimension).

We can access this variable by indexing our UxDataset to obtain a UxDataArray

uxds["bottomDepth"]
<xarray.UxDataArray 'bottomDepth' (n_face: 1791)> Size: 14kB
[1791 values with dtype=float64]
Dimensions without coordinates: n_face

You can access the Grid instance using the .uxgrid attribute, which is linked to every UxDataset and UxDataArray

uxds.uxgrid
<uxarray.Grid>
Original Grid Type: MPAS
Grid Dimensions:
  * n_node: 3947
  * n_edge: 5754
  * n_face: 1791
  * n_max_face_nodes: 6
  * n_max_face_edges: 6
  * n_max_face_faces: 6
  * n_max_node_faces: 3
  * two: 2
  * n_nodes_per_face: (1791,)
Grid Coordinates (Spherical):
  * node_lon: (3947,)
  * node_lat: (3947,)
  * edge_lon: (5754,)
  * edge_lat: (5754,)
  * face_lon: (1791,)
  * face_lat: (1791,)
Grid Coordinates (Cartesian):
  * node_x: (3947,)
  * node_y: (3947,)
  * node_z: (3947,)
  * edge_x: (5754,)
  * edge_y: (5754,)
  * edge_z: (5754,)
  * face_x: (1791,)
  * face_y: (1791,)
  * face_z: (1791,)
Grid Connectivity Variables:
  * face_node_connectivity: (1791, 6)
  * face_edge_connectivity: (1791, 6)
  * face_face_connectivity: (1791, 6)
  * edge_node_connectivity: (5754, 2)
  * edge_face_connectivity: (5754, 2)
  * node_face_connectivity: (3947, 3)
Grid Descriptor Variables:
  * face_areas: (1791,)
  * n_nodes_per_face: (1791,)
  * edge_face_distances: (5754,)
  * edge_node_distances: (5754,)

Each UxDataArray under a UxDataset is linked to the same Grid object

uxds.uxgrid == uxds["bottomDepth"].uxgrid
True