Reading Structured Grids#
UXarray supports reading structured grids and representing them as unstructured grids. This user-guide section will discuss how to load in structured grids using UXarray.
import uxarray as ux
import xarray as xr
import cartopy.crs as ccrs
import warnings
warnings.filterwarnings("ignore")
Data#
For this notebook, we will be using datasets from the Xarray tutorial.
ds_air_temp = xr.tutorial.open_dataset("air_temperature")
ds_air_temp
<xarray.Dataset> Size: 31MB Dimensions: (lat: 25, time: 2920, lon: 53) Coordinates: * lat (lat) float32 100B 75.0 72.5 70.0 67.5 65.0 ... 22.5 20.0 17.5 15.0 * lon (lon) float32 212B 200.0 202.5 205.0 207.5 ... 325.0 327.5 330.0 * time (time) datetime64[ns] 23kB 2013-01-01 ... 2014-12-31T18:00:00 Data variables: air (time, lat, lon) float64 31MB ... Attributes: Conventions: COARDS title: 4x daily NMC reanalysis (1948) description: Data is from NMC initialized reanalysis\n(4x/day). These a... platform: Model references: http://www.esrl.noaa.gov/psd/data/gridded/data.ncep.reanaly...
ds_ersstv5 = xr.tutorial.open_dataset("ersstv5")
ds_ersstv5
<xarray.Dataset> Size: 40MB Dimensions: (lat: 89, lon: 180, time: 624, nbnds: 2) Coordinates: * lat (lat) float32 356B 88.0 86.0 84.0 82.0 ... -84.0 -86.0 -88.0 * lon (lon) float32 720B 0.0 2.0 4.0 6.0 ... 352.0 354.0 356.0 358.0 * time (time) datetime64[ns] 5kB 1970-01-01 1970-02-01 ... 2021-12-01 Dimensions without coordinates: nbnds Data variables: time_bnds (time, nbnds) float64 10kB ... sst (time, lat, lon) float32 40MB ... Attributes: (12/37) climatology: Climatology is based on 1971-2000 SST, Xue, Y.... description: In situ data: ICOADS2.5 before 2007 and NCEP i... keywords_vocabulary: NASA Global Change Master Directory (GCMD) Sci... keywords: Earth Science > Oceans > Ocean Temperature > S... instrument: Conventional thermometers source_comment: SSTs were observed by conventional thermometer... ... ... creator_url_original: https://www.ncei.noaa.gov license: No constraints on data access or use comment: SSTs were observed by conventional thermometer... summary: ERSST.v5 is developed based on v4 after revisi... dataset_title: NOAA Extended Reconstructed SST V5 data_modified: 2022-06-07
Grid#
A structured grid can be converted into an unstructured grid by using the ux.Grid.from_structured()
class method. An xarray.Dataset
can be passed in, with CF-compliant longitude and latitude coordinates parsed and converted to an unstructured grid.
uxgrid = ux.Grid.from_structured(ds_air_temp)
uxgrid
<uxarray.Grid> Original Grid Type: Structured Grid Dimensions: * n_node: 1404 * n_face: 1325 * n_max_face_nodes: 4 * n_nodes_per_face: (1325,) Grid Coordinates (Spherical): * node_lon: (1404,) * node_lat: (1404,) Grid Coordinates (Cartesian): Grid Connectivity Variables: * face_node_connectivity: (1325, 4) Grid Descriptor Variables: * n_nodes_per_face: (1325,)
You can also manually pass in longitude and latitude coordinates.
uxgrid = ux.Grid.from_structured(lon=ds_air_temp.lon, lat=ds_air_temp.lat)
uxgrid
<uxarray.Grid> Original Grid Type: structured Grid Dimensions: * n_node: 1404 * n_face: 1325 * n_max_face_nodes: 4 * n_nodes_per_face: (1325,) Grid Coordinates (Spherical): * node_lon: (1404,) * node_lat: (1404,) Grid Coordinates (Cartesian): Grid Connectivity Variables: * face_node_connectivity: (1325, 4) Grid Descriptor Variables: * n_nodes_per_face: (1325,)
uxgrid.plot(
title="Structured Grid loaded using UXarray",
backend="matplotlib",
width=1000,
coastline=True,
)
Dataset#
If you have a dataset that contains data variables, you can convert the entire xarray.Dataset
into a uxarray.UxDataset
using the ux.UxDataset.from_structured()
class method.
uxds = ux.UxDataset.from_structured(ds_air_temp)
uxds
<xarray.UxDataset> Size: 31MB Dimensions: (time: 2920, n_face: 1325) Coordinates: * time (time) datetime64[ns] 23kB 2013-01-01 ... 2014-12-31T18:00:00 Dimensions without coordinates: n_face Data variables: air (time, n_face) float64 31MB 241.2 242.5 243.5 ... 296.5 296.2 295.7
uxds["air"][0].plot(
title="Structured Grid with Data loaded using UXarray",
backend="matplotlib",
width=500,
coastline=True,
)
Limitations of Structured Grids#
Structured grids are a common choice in climate modeling due to their ease of implementation. However, they come with many limitations, which have been a driving factor in the adoption of unstructured grids.
Limited Flexibility in Handling Complex Geometries#
Structured grids, with their regular and grid-aligned cells, struggle to capture complex geometries. This limitation forces modelers to use a higher number of grid points to approximate complex features, which can lead to increased computational costs and reduced simulation efficiency. Additionally, the inability to precisely represent intricate boundaries can result in inaccuracies in modeling climate processes that are influenced by these geographical features, potentially affecting the reliability of the simulation outcomes.
Difficulty in Local Grid Refinement#
Climate phenomena like tropical cyclones, atmospheric fronts, and localized convection require high-resolution grids to be accurately modeled. Structured grids make it challenging to refine the grid locally in regions where such detailed resolution is needed without uniformly increasing the grid resolution across the entire model domain. This uniform refinement results in inefficient use of computational resources, as large areas of the model may have unnecessarily high resolution where it is not required.
Pole Point Singularities#
Pole point singularities are a significant challenge in climate models that utilize structured latitude-longitude grids. In such grid systems, the lines of longitude converge as they approach the Earth’s poles, causing the grid cells to become increasingly small near these regions. This convergence leads to several issues:
Numerical Instability: The drastically reduced cell size near the poles can cause numerical methods to become unstable, requiring smaller time steps to maintain accuracy and stability in simulations.
Accuracy Problems: The distortion of grid cell shapes near the poles can lead to inaccuracies in representing physical processes, such as atmospheric circulation and ocean currents, which are critical for realistic climate simulations.
uxds = ux.UxDataset.from_structured(ds_ersstv5)
Below is a plot of Sea Surface Temperature on a structured grid.
uxds["sst"][0].plot(
projection=ccrs.Orthographic(central_latitude=60),
periodic_elements="split",
coastline=True,
grid=True,
title="SST Near North Pole",
)
If we expose the grid structure, we can observe the singularity at the poles.
uxds.uxgrid.plot(
projection=ccrs.Orthographic(central_latitude=60),
periodic_elements="split",
width=600,
linewidth=1,
coastline=True,
title="Grid Near NorthPole",
)