Accessing Grid Information#

Unstructured grids can be represented in one of many different conventions (UGRID, SCRIP, EXODUS, etc). These conventions have different definitions and representations of the attributes and variables used to describe the unstructured grid topology. Even more, the UGRID conventions does not enforce standard variable namings for most of the attributes and variables (other than just a few required ones).

UXarray unifies all of these conventions at the data loading step by representing grids in the UGRID convention regardless of the original grid type that is read in from the file. Furthermore, it uses a set of standardized names for topology attributes and variables, while still providing the user with the original attribute names and variables that came from the grid definition file.

Overview#

This notebook will showcase the different methods available for accessing the grid topology attributes and variables stored in the UXarray.Grid object.

For more details on how to load in data, check out our previous usage example

Methods

  1. Indexing with Original Variable Names

  2. Indexing with UXarray Variable Dictionary

  3. UXarray’s Standardized UGRID Names (Most convenient)

Data#

We will be using two grid files, both of which are in the UGRID convention. However, the key difference between them is the names used to describe the attributes and variables.

Let us first read in the data:

import uxarray as ux
import xarray as xr
# Base data path
base_path = "../../test/meshfiles/"

# Path to Grid files
ugrid_01_path = base_path + "outCSne30.ug"
ugrid_02_path = base_path + "geoflow-small/grid.nc"

# Load grid files and create UXarray Grid objects
ugrid_01_ds = xr.open_dataset(ugrid_01_path)
ugrid_02_ds = xr.open_dataset(ugrid_02_path)

ugrid_01 = ux.Grid(ugrid_01_ds)
ugrid_02 = ux.Grid(ugrid_02_ds)

The output of the bottom cell showcases the slight differences in variable names:

# Extract variable names
ugrid_01_names = list(ugrid_01.ds.keys()) + \
                 list(ugrid_01.ds.coords) + \
                 list(ugrid_01.ds.dims)
ugrid_02_names = list(ugrid_02.ds.keys()) + \
                 list(ugrid_02.ds.coords) + \
                 list(ugrid_02.ds.dims)

print("\nAttribute and variable names for each grid:")
print("ugrid_01 variable names:", ugrid_01_names)
print("ugrid_02 variable names:", ugrid_02_names)
Attribute and variable names for each grid:
ugrid_01 variable names: ['Mesh2', 'Mesh2_face_nodes', 'Mesh2_node_x', 'Mesh2_node_y', 'nMesh2_face', 'nMaxMesh2_face_nodes', 'nMesh2_node']
ugrid_02 variable names: ['mesh', 'mesh_face_nodes', 'mesh_depth', 'mesh_node_x', 'mesh_node_y', 'nMeshFaces', 'nFaceNodes', 'nMeshNodes', 'meshLayers']

1. Indexing with Original Variable Names#

The simplest approach is to use the original variable name as an index into the grid dataset, Grid.ds. Since ugrid_01 and ugrid_02 have different names for most of their topology attributes and variables, the index for accessing them will be different for both grids.

x_1 = ugrid_01.ds['Mesh2_node_x']
y_1 = ugrid_01.ds['Mesh2_node_y']
face_nodes_1 = ugrid_01.ds['Mesh2_face_nodes']
n_face_nodes_1 = ugrid_01.ds['nMaxMesh2_face_nodes']

x_2 = ugrid_02.ds['mesh_node_x']
y_2 = ugrid_02.ds['mesh_node_y']
face_nodes_2 = ugrid_02.ds['mesh_face_nodes']
n_face_nodes_2 = ugrid_02.ds['nFaceNodes']

2. Indexing with UXarray Variable Dictionary#

UXarray provides a dictionary, Grid.ds_var_names, to map the original topology attribute and variable names that come from the grid file into a standardized set of names. In other words, the dictionary uses a standardized set of UGRID attribute and variable names as keys, and the original variable names that come from the grid file as values.

This allows us to use the same index into either dataset. However, this makes the indexing code much longer than the previous method.

var_names_dict = ugrid_01.ds_var_names
x_1 = ugrid_01.ds[var_names_dict['Mesh2_node_x']]
y_1 = ugrid_01.ds[var_names_dict['Mesh2_node_y']]
face_nodes_1 = ugrid_01.ds[var_names_dict['Mesh2_face_nodes']]
n_face_nodes_1 = ugrid_01.ds[var_names_dict['nMesh2_node']]

var_names_dict = ugrid_02.ds_var_names
x_2 = ugrid_02.ds[var_names_dict['Mesh2_node_x']]
y_2 = ugrid_02.ds[var_names_dict['Mesh2_node_y']]
face_nodes_2 = ugrid_02.ds[var_names_dict['Mesh2_face_nodes']]
n_face_nodes_2 = ugrid_02.ds[var_names_dict['nMesh2_node']]

Please note, for instance, we accessed the actual variable mesh_node_x of ugrid_02 via indexing the dictionary with the standardized name Mesh2_node_x, likewise in ugrid_01.

3. UXarray’s Standardized UGRID Names#

The last way of accessing grid topology attributes and variables is to use the standardized UGRID namings provided by UXarray. This method still utilizes the dictionary, ds_var_names, under the hood to return a reference to the variable or attribute that is stored withing UXarray.Grid.ds.

This eliminates the need to remember the original variable names and needing to index through the ds_var_names dictionary. Because of this, we find this as the most convenient approach.

x_1 = ugrid_01.Mesh2_node_x
y_1 = ugrid_01.Mesh2_node_y
face_nodes_1 = ugrid_01.Mesh2_face_nodes
n_face_nodes_1 = ugrid_01.nMesh2_node

x_2 = ugrid_02.Mesh2_node_x
y_2 = ugrid_02.Mesh2_node_y
face_nodes_2 = ugrid_02.Mesh2_face_nodes
n_face_nodes_2 = ugrid_02.nMesh2_node

In conclusion, there are three ways of accessing the grid attributes and variables. Even though the UXarray developers recommend using the standardized UGRID names method, users can find each various pros/cons with each of these access ways.