Cross-Sections#

This section demonstrates how to extract cross-sections from an unstructured grid using UXarray, which allows the analysis and visualization across slices of grids. Cross-sections can be performed directly on a ux.Grid object or on a ux.UxDataArray

import cartopy.crs as ccrs
import geoviews as gv
import geoviews.feature as gf

import uxarray as ux

projection = ccrs.Robinson()
base_path = "../../test/meshfiles/ugrid/outCSne30/"
grid_path = base_path + "outCSne30.ug"
data_path = base_path + "outCSne30_vortex.nc"

uxds = ux.open_dataset(grid_path, data_path)
uxds["psi"].plot(
    cmap="inferno",
    periodic_elements="split",
    projection=projection,
    title="Global Plot",
)

Constant Latitude#

Cross-sections along constant latitude lines can be obtained by using the .cross_section.constant_latitude(lat) method. The sliced grid will be made up of the faces that contain at least one edge that intersects with a line of constant latitude.

For example, we can obtain a cross-section at 0 degrees latitude by doing the following:

lat = 0

uxda_constant_lat = uxds["psi"].cross_section.constant_latitude(lat)

Since the result is a new UxDataArray, we can directly plot the result to see the cross-section.

(
    uxda_constant_lat.plot(
        rasterize=False,
        backend="bokeh",
        cmap="inferno",
        projection=projection,
        global_extent=True,
        coastline=True,
        title=f"Cross Section at {lat} degrees latitude",
    )
    * gf.grid(projection=projection)
)
/home/docs/checkouts/readthedocs.org/user_builds/uxarray/conda/latest/lib/python3.13/site-packages/cartopy/io/__init__.py:241: DownloadWarning: Downloading: https://naturalearth.s3.amazonaws.com/110m_physical/ne_110m_graticules_30.zip
  warnings.warn(f'Downloading: {url}', DownloadWarning)
/home/docs/checkouts/readthedocs.org/user_builds/uxarray/conda/latest/lib/python3.13/site-packages/cartopy/io/__init__.py:241: DownloadWarning: Downloading: https://naturalearth.s3.amazonaws.com/110m_physical/ne_110m_coastline.zip
  warnings.warn(f'Downloading: {url}', DownloadWarning)

You can also perform operations on the cross-section, such as taking the mean.

print(f"Global Mean: {uxds['psi'].data.mean()}")
print(f"Mean at {lat} degrees lat: {uxda_constant_lat.data.mean()}")
Global Mean: 1.0000000015332655
Mean at 0 degrees lat: 1.0000000052455802

Constant Longitude#

Cross-sections along constant longitude lines can be obtained using the .cross_section.constant_longitude(lon) method. The sliced grid will be made up of the faces that contain at least one edge that intersects with a line of constant longitude.

lon = 90

uxda_constant_lon = uxds["psi"].cross_section.constant_longitude(lon)
(
    uxda_constant_lon.plot(
        rasterize=False,
        backend="bokeh",
        cmap="inferno",
        projection=projection,
        global_extent=True,
        coastline=True,
        title=f"Cross Section at {lon} degrees longitude",
        periodic_elements="split",
    )
    * gf.grid(projection=projection)
)

Constant Latitude Interval#

Cross-sections between two lines of latitudes can be obtained using the .cross_section.constant_lats_interval(lats) method. The sliced grid will contain faces that are strictly between the latitude interval.

lats = [-20, 20]

uxda_constant_lat_interval = uxds["psi"].cross_section.constant_latitude_interval(lats)
(
    uxda_constant_lat_interval.plot(
        rasterize=False,
        backend="bokeh",
        cmap="inferno",
        projection=projection,
        global_extent=True,
        coastline=True,
        title=f"Cross Section between {lats[0]} and {lats[1]} degrees latitude",
        periodic_elements="split",
    )
    * gf.grid(projection=projection)
)

Constant Longitude Interval#

Cross-sections between two lines of longitude can be obtained using the .cross_section.constant_lons_interval(lons) method. The sliced grid will contain faces that are strictly between the longitude interval.

lons = [-25, 25]

uxda_constant_lon_interval = uxds["psi"].cross_section.constant_longitude_interval(lats)
(
    uxda_constant_lon_interval.plot(
        rasterize=False,
        backend="bokeh",
        cmap="inferno",
        projection=projection,
        global_extent=True,
        coastline=True,
        title=f"Cross Section between {lons[0]} and {lons[1]} degrees longitude",
        periodic_elements="split",
    )
    * gf.grid(projection=projection)
)

Arbitrary Great Circle Arc (GCA)#

Warning

Arbitrary great circle arc cross sections are not yet implemented.