Quick overview#

Here are a number of quick examples of how to get started with xugrid. More detailed explanation can be found in the rest of the documentation.

We’ll start by importing a few essential packages.

import numpy as np
import xarray as xr

import xugrid as xu

Create a UgridDataArray#

There are three ways to create a UgridDataArray:

  • From an xarray Dataset containing the grid topology stored according to the UGRID conventions.

  • From a xugrid Ugrid object and an xarray DataArray containing the data.

  • From a UGRID netCDF file, via xugrid.open_dataset().

From xarray Dataset#

xugrid will automatically find the UGRID topological variables, and separate them from the main data variables.

Details on the required variables can be found in the UGRID conventions. For 1D and 2D UGRID topologies, the required variables are:

  • x-coordinates of the nodes

  • y-coordinates of the nodes

  • edge node connectivity (1D) or face node connectivity (2D)

  • a “dummy” variable storing the names of the above variables in its attributes

We’ll start by fetching a dataset:

ds = xu.data.adh_san_diego(xarray=True)
ds
<xarray.Dataset> Size: 4MB
Dimensions:                 (node: 9140, time: 49, face: 16869, nmax_face: 3)
Coordinates:
  * node                    (node) int64 73kB 0 1 2 3 4 ... 9136 9137 9138 9139
  * time                    (time) datetime64[ns] 392B 2000-01-01 ... 2000-01-02
    node_x                  (node) float64 73kB ...
    node_y                  (node) float64 73kB ...
Dimensions without coordinates: face, nmax_face
Data variables:
    elevation               (node) float64 73kB ...
    depth                   (time, node) float64 4MB ...
    mesh2d                  int32 4B ...
    face_node_connectivity  (face, nmax_face) float64 405kB ...


There are a number of topology coordinates and variables: node_x and node_y, mesh2d and face_node_connectivity. The dummy variable is mesh2d contains only a 0 for data; its attributes contain a mapping of UGRID roles to dataset variables.

We can convert this dataset to a UgridDataset which will automatically separate the variables:

uds = xu.UgridDataset(ds)
uds
<xarray.Dataset> Size: 4MB
Dimensions:    (node: 9140, time: 49)
Coordinates:
  * time       (time) datetime64[ns] 392B 2000-01-01 ... 2000-01-02
    node_x     (node) float64 73kB ...
    node_y     (node) float64 73kB ...
  * node       (node) int64 73kB 0 1 2 3 4 5 6 ... 9134 9135 9136 9137 9138 9139
Data variables:
    elevation  (node) float64 73kB ...
    depth      (time, node) float64 4MB ...


We can then grab one of the data variables as usual for xarray:

elev = uds["elevation"]
elev
<xarray.DataArray 'elevation' (node: 9140)> Size: 73kB
[9140 values with dtype=float64]
Coordinates:
    node_x   (node) float64 73kB ...
    node_y   (node) float64 73kB ...
  * node     (node) int64 73kB 0 1 2 3 4 5 6 ... 9134 9135 9136 9137 9138 9139


From Ugrid and DataArray#

Alternatively, we can build a Ugrid topology object first from vertices and connectivity numpy arrays, for example when using the topology data generated by a mesh generator (at which stage there is no data asssociated with the nodes, edges, or faces).

There are many ways to construct such arrays, typically via mesh generators or Delaunay triangulation, but we will construct two simple triangles and some data by hand here:

nodes = np.array([[0, 0], [0, 1.1], [1, 0], [1, 1]])
faces = np.array([[2, 3, 0], [3, 1, 0]])
fill_value = -1

grid = xu.Ugrid2d(nodes[:, 0], nodes[:, 1], fill_value, faces)
da = xr.DataArray(
    data=[1.0, 2.0],
    dims=[grid.face_dimension],
)
uda = xu.UgridDataArray(da, grid)
uda
<xarray.DataArray (mesh2d_nFaces: 2)> Size: 16B
array([1., 2.])
Coordinates:
  * mesh2d_nFaces  (mesh2d_nFaces) int64 16B 0 1


From netCDF file#

xugrid.open_dataset() is demonstrated in the last section of this guide. Internally, it opens the netCDF as a regular dataset, then converts it as seen in the first example.

Plotting#

elev.ugrid.plot(cmap="viridis")
quick overview
<matplotlib.collections.PolyCollection object at 0x7f98ebff12e0>

Data selection#

A UgridDataArray behaves identical to an xarray DataArray:

whole = xu.data.disk()["face_z"]

To select based on the topology, use the .ugrid attribute:

subset = whole.ugrid.sel(y=slice(5.0, None))
subset.ugrid.plot()
quick overview
<matplotlib.collections.PolyCollection object at 0x7f98e89875c0>

Note

ugrid.sel() currently only supports data on the faces for 2D topologies, and data on edges for 1D topologies. More flexibility may be added.

Computation#

Computation on DataArrays is unchanged from xarray:

uda + 10.0
<xarray.DataArray (mesh2d_nFaces: 2)> Size: 16B
array([11., 12.])
Coordinates:
  * mesh2d_nFaces  (mesh2d_nFaces) int64 16B 0 1


Geopandas#

Xugrid objects provide a number of conversion functions from and to geopandas GeoDataFrames using xugrid.UgridDataset.from_geodataframe(). Note that storing large grids as GeoDataFrames can be very inefficient.

gdf = uda.ugrid.to_geodataframe(name="test")
gdf
test geometry
mesh2d_nFaces
0 1.0 POLYGON ((1.00000 0.00000, 1.00000 1.00000, 0....
1 2.0 POLYGON ((1.00000 1.00000, 0.00000 1.10000, 0....


Conversion from Geopandas is easy too:

xu.UgridDataset.from_geodataframe(gdf)
<xarray.Dataset> Size: 32B
Dimensions:        (mesh2d_nFaces: 2)
Coordinates:
  * mesh2d_nFaces  (mesh2d_nFaces) int64 16B 0 1
Data variables:
    test           (mesh2d_nFaces) float64 16B 1.0 2.0


XugridDatasets#

Like an Xarray Dataset, a UgridDataset is a dict-like container of UgridDataArrays. It is required that they share the same grid topology; but the individual DataArrays may be located on different aspects of the grid (nodes, faces, edges).

xu.data.disk()
<xarray.Dataset> Size: 19kB
Dimensions:        (mesh2d_nNodes: 217, mesh2d_nFaces: 384, mesh2d_nEdges: 600)
Coordinates:
  * mesh2d_nFaces  (mesh2d_nFaces) int64 3kB 0 1 2 3 4 5 ... 379 380 381 382 383
  * mesh2d_nNodes  (mesh2d_nNodes) int64 2kB 0 1 2 3 4 5 ... 212 213 214 215 216
  * mesh2d_nEdges  (mesh2d_nEdges) int64 5kB 0 1 2 3 4 5 ... 595 596 597 598 599
Data variables:
    node_z         (mesh2d_nNodes) float64 2kB 1.933 2.091 1.875 ... 5.688 7.491
    face_z         (mesh2d_nFaces) float64 3kB 1.737 1.918 2.269 ... 5.408 6.424
    edge_z         (mesh2d_nEdges) float64 5kB 1.989 1.875 1.8 ... 4.909 6.544


A UgridDataset may be initialized without data variables, but this requires a grid object:

new_uds = xu.UgridDataset(grids=uds.ugrid.grids)
new_uds
<xarray.Dataset> Size: 0B
Dimensions:  ()
Data variables:
    *empty*


We can then add variables one-by-one, as we might with an xarray Dataset:

new_uds["elevation"] = elev
new_uds
<xarray.Dataset> Size: 292kB
Dimensions:    (node: 9140)
Coordinates:
    node_x     (node) float64 73kB ...
    node_y     (node) float64 73kB ...
  * node       (node) int64 73kB 0 1 2 3 4 5 6 ... 9134 9135 9136 9137 9138 9139
Data variables:
    elevation  (node) float64 73kB ...


Write netCDF files#

Once again like xarray, NetCDF is the recommended file format for xugrid objects. Xugrid automatically stores the grid topology according to the UGRID conventions and merges it with the main dataset containing the data variables before writing.

uds.ugrid.to_netcdf("example-ugrid.nc")
xu.open_dataset("example-ugrid.nc")
<xarray.Dataset> Size: 4MB
Dimensions:    (node: 9140, time: 49)
Coordinates:
    node_x     (node) float64 73kB ...
    node_y     (node) float64 73kB ...
  * time       (time) datetime64[ns] 392B 2000-01-01 ... 2000-01-02
  * node       (node) int64 73kB 0 1 2 3 4 5 6 ... 9134 9135 9136 9137 9138 9139
Data variables:
    elevation  (node) float64 73kB ...
    depth      (time, node) float64 4MB ...
Attributes:
    Conventions:  CF-1.9 UGRID-1.0


Total running time of the script: (0 minutes 0.541 seconds)

Gallery generated by Sphinx-Gallery