Cross-Sections#
UXarray allows for cross-sections to be performed along arbitrary great-circle arcs (GCAs) and lines of constant latitude or longitude.
This enables analysis of workflows such as Vertical Cross-Sections when vertical dimensions (i.e., levels) are present, temporal visualizations (e.g., Hovmöller diagram), and even visualizations of arbitrary dimensions along GCAs or constant lat/lon lines.
import matplotlib.pyplot as plt
import numpy as np
import uxarray as ux
Overview#
Cross-sections can be performed by using the .cross_section
method on a ux.UxDataArray
. The original data variable is sampled along steps
number of evenly spaced points, with the result stored as a xr.DataArray
, since there is no longer any grid information associated with the sampled result.
Note
The cross_section.bounding_box
, cross_section.bounding_circle
and cross_section.nearest_neighbor
methods will be deprecated in a future release, as they are now housed under the subset
accessor.
Let’s open up our data and define a small helper function for labeling the sampled coordinates.
def set_lon_lat_xticks(ax, cross_section, n_ticks=6):
"""Utility function to draw stacked lat/lon points along the sampled cross-section"""
da = cross_section
N = da.sizes["steps"]
tick_pos = np.linspace(0, N - 1, n_ticks, dtype=int)
lons = da["lon"].values[tick_pos]
lats = da["lat"].values[tick_pos]
tick_labels = []
for lon, lat in zip(lons, lats):
lon_dir = "E" if lon >= 0 else "W"
lat_dir = "N" if lat >= 0 else "S"
tick_labels.append(f"{abs(lon):.2f}°{lon_dir}\n{abs(lat):.2f}°{lat_dir}")
ax.set_xticks(tick_pos)
ax.set_xticklabels(tick_labels)
ax.set_xlabel("Longitude\nLatitude")
plt.tight_layout()
return fig, ax
grid_path = "../../test/meshfiles/scrip/ne30pg2/grid.nc"
data_path = "../../test/meshfiles/scrip/ne30pg2/data.nc"
uxds = ux.open_dataset(grid_path, data_path)
uxds["RELHUM"]
<xarray.UxDataArray 'RELHUM' (lev: 72, n_face: 21600)> Size: 6MB [1555200 values with dtype=float32] Coordinates: * lev (lev) float64 576B 0.1238 0.1828 0.2699 ... 986.2 993.8 998.5 time object 8B ... Dimensions without coordinates: n_face Attributes: mdims: 1 units: percent long_name: Relative humidity standard_name: relative_humidity cell_methods: time: mean
Arbitrary Great Circle Arc (GCA)#
A cross‑section can be performed between two arbitrary (lon,lat) points, which will form a geodesic arc.
start_point = (-45, -45)
end_point = (45, 45)
cross_section_gca = uxds["RELHUM"].cross_section(
start=start_point, end=end_point, steps=100
)
fig, ax = plt.subplots()
cross_section_gca.plot(ax=ax)
set_lon_lat_xticks(ax, cross_section_gca)
ax.set_title(f"Cross Section between {start_point} and {end_point}")
ax.invert_yaxis()

Constant Latitude#
A constant‐latitude cross‐section samples data along a horizontal line at a fixed latitude.
lat = 45
cross_section_const_lat = uxds["RELHUM"].cross_section(lat=lat, steps=100)
fig, ax = plt.subplots()
cross_section_const_lat.plot(ax=ax)
set_lon_lat_xticks(ax, cross_section_const_lat)
ax.set_title(f"Cross Section at {lat}° latitude.")
ax.invert_yaxis()

Constant Longitude#
A constant‐longitude cross‐section samples data along a vertical line at a fixed longitude
lon = 0
cross_section_const_lon = uxds["RELHUM"].cross_section(lon=lon, steps=100)
fig, ax = plt.subplots()
cross_section_const_lon.plot(ax=ax)
set_lon_lat_xticks(ax, cross_section_const_lon)
ax.set_title(f"Cross Section at {lon}° longitude.")
ax.invert_yaxis()
