UXarray logo

UXarray Dataset & DataArray

In this section, you’ll learn:

  • What are UxDataset and UxDataArray?

  • How do they relate to their Xarray counterparts?


Overview

As stated before in this chapter, UxDataset and UxDataArray are the UXarray classes that are inherited from Xarray counterparts and enable to perform analysis and visualization on unstructured mesh datasets.

An instance of UxDataset can be thought of as a collection of data variables, i.e. UxDataArray objects, that reside on a particular unstructured grid as defined in the class property (UxDataset.uxgrid or UxDataArray.uxgrid).

In other words, each of these data structures is linked to its own Grid object through uxgrid to enable grid-awareness, and for UxDataArrays that belong to a specific UxDataset, uxgrid of all of them points to the same Grid object.

Unstructured Mesh Datasets

Since most of the unstructured grid datasets contain a single grid definition file along with one or more data files, UXarray’s core API allows to open such datasets using those files.

Let the following dataset have a single data file along with the grid definition:

grid_path = "../../meshfiles/outCSne30.grid.ug"
data_path = "../../meshfiles/outCSne30.data.nc"

Loading a UxDataset

To open a grid file with a single data file, the uxarray.open_dataset() method can be used. UXarray also provides uxarray.open_mfdataset() for opening multiple data files at once (refer to the UXarray Opening Multipled Data Files Documentation for that), but this notebook will cover only the former.

import uxarray as ux

uxds = ux.open_dataset(grid_path, data_path)
uxds
<xarray.UxDataset> Size: 43kB
Dimensions:  (n_face: 5400)
Dimensions without coordinates: n_face
Data variables:
    psi      (n_face) float64 43kB ...

Grid Accessor

Let us first see the Gridobject this dataset object is linked to:

uxds.uxgrid
<uxarray.Grid>
Original Grid Type: UGRID
Grid Dimensions:
  * n_node: 5402
  * n_face: 5400
  * n_max_face_nodes: 4
  * n_nodes_per_face: (5400,)
Grid Coordinates (Spherical):
  * node_lon: (5402,)
  * node_lat: (5402,)
Grid Coordinates (Cartesian):
Grid Connectivity Variables:
  * face_node_connectivity: (5400, 4)
Grid Descriptor Variables:
  * n_nodes_per_face: (5400,)

Accessing Data Variables (UxDataArray)

Similar to the Xarray conterparts, an UxDataset can have multiple UxDataArrays in it, which are for different data variables. A variable of interest can be accessed in the same way Xarray provides:

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

Grid Accessor

The Gridobject of this data array is the same Grid object as uxds’s above:

uxda = uxds["psi"]
uxda.uxgrid
<uxarray.Grid>
Original Grid Type: UGRID
Grid Dimensions:
  * n_node: 5402
  * n_face: 5400
  * n_max_face_nodes: 4
  * n_nodes_per_face: (5400,)
Grid Coordinates (Spherical):
  * node_lon: (5402,)
  * node_lat: (5402,)
Grid Coordinates (Cartesian):
Grid Connectivity Variables:
  * face_node_connectivity: (5400, 4)
Grid Descriptor Variables:
  * n_nodes_per_face: (5400,)

Relationship to Xarray

For users coming from an Xarray background, much of UXarray’s design is familiar thanks to the UXarray’s inheritance of the core data structures from Xarray’s counterparts and employment of design choices such as accessors. Inheritance has been chosen for the UXarray design both to provide the users with a Xarray-like user experience and to use the built-in Xarray functionality whenever possible while overwrite the others for unstructured grid arithmetic.

As a simple Xarray example that would replicate the above UXarray dataset opening and variable access with UXarray, please see the following code:

# Import
import xarray as xr

# Structured grid data file
data_path_structured = "../../meshfiles/outCSne30.structured.nc"

# Open dataset with Xarray
xrds = xr.open_dataset(data_path_structured)
xrds
<xarray.Dataset> Size: 30kB
Dimensions:  (lat: 45, lon: 80)
Coordinates:
  * lat      (lat) int64 360B -90 -86 -82 -78 -74 -70 -66 ... 66 70 74 78 82 86
  * lon      (lon) float64 640B -180.0 -175.5 -171.0 ... 166.5 171.0 175.5
Data variables:
    psi      (lat, lon) float64 29kB ...
xrda = xrds["psi"]
xrda
<xarray.DataArray 'psi' (lat: 45, lon: 80)> Size: 29kB
[3600 values with dtype=float64]
Coordinates:
  * lat      (lat) int64 360B -90 -86 -82 -78 -74 -70 -66 ... 66 70 74 78 82 86
  * lon      (lon) float64 640B -180.0 -175.5 -171.0 ... 166.5 171.0 175.5
Attributes:
    regrid_method:  nearest_s2d

See also:

For further information, please visit the UXarray Inheritance from Xarray documentation

What is next?

The next section is the final section in this chapter that will provide information about how to operate spatial selection methods on UxDataArray.