Plotting with Matplotlib#
Although UXarray’s primary plotting API leverages the HoloViz ecosystem, users can still create visualizations using Matplotlib by converting UXarray objects into compatible Matplotlib collections, such as LineCollection and PolyCollection.
This user guide will cover:
Converting a
Grid
to aLineCollection
Converting a
UxDataArray
to aPolyCollection
Using Geographic Projections & Elements
Handling periodic elements along the antimeridian
import uxarray as ux
import matplotlib.pyplot as plt
import cartopy.crs as ccrs
import cartopy.feature as cfeature
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)
Visualize Grid Topology with LineCollection
#
The Grid.to_linecollection()
method can be used to convert a Grid
instance into a matplotlib.collections.LineCollection
instance. It represents a collection of lines that represent the edges of an unstructured grid.
lc = uxds.uxgrid.to_linecollection(colors="black", linewidths=0.5)
lc
<matplotlib.collections.LineCollection at 0x7333b9919be0>
Once we have converted our Grid
to a LineCollection
, we can directly use Matplotlib.
fig, ax = plt.subplots(
1,
1,
figsize=(10, 10),
constrained_layout=True,
subplot_kw={"projection": ccrs.PlateCarree()},
)
ax.add_feature(cfeature.LAND)
ax.add_feature(cfeature.COASTLINE)
ax.add_collection(lc)
ax.set_global()
ax.set_title("LineCollection Plot")
plt.show()
/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_land.zip
warnings.warn(f'Downloading: {url}', DownloadWarning)

We can also specify a projection directly when constructing a LineCollection
, which provides better performance compared to re-projecting the data with Matplotlib during figure creation.
projection = ccrs.Robinson()
lc_direct_projection = uxds.uxgrid.to_linecollection(
override=True, colors="black", linewidths=0.5, projection=projection
)
fig, ax = plt.subplots(
1,
1,
figsize=(10, 10),
constrained_layout=True,
subplot_kw={"projection": projection},
)
ax.add_feature(cfeature.LAND)
ax.add_feature(cfeature.COASTLINE)
ax.add_collection(lc_direct_projection)
ax.set_global()
ax.set_title("LineCollection Plot (Explicit Projection)")
plt.show()

Visualize Data with PolyCollection
#
The Grid.to_polycollection()
method can be used to convert a UxDataArray
containing a face-centered data variable into a matplotlib.collections.PolyCollection
instance. It represents a collection of polygons that represent the faces of an unstructured grid, shaded using the values of the face-centered data variable.
pc = uxds["psi"].to_polycollection()
pc
<matplotlib.collections.PolyCollection at 0x7333553d39d0>
# disables grid lines
pc.set_antialiased(False)
pc.set_cmap("plasma")
fig, ax = plt.subplots(
1,
1,
figsize=(10, 5),
facecolor="w",
constrained_layout=True,
subplot_kw=dict(projection=ccrs.PlateCarree()),
)
ax.add_feature(cfeature.COASTLINE)
ax.add_feature(cfeature.BORDERS)
ax.add_collection(pc)
ax.set_global()
plt.title("PolyCollection Plot with Projection & Features")
Text(0.5, 1.0, 'PolyCollection Plot with Projection & Features')
/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_cultural/ne_110m_admin_0_boundary_lines_land.zip
warnings.warn(f'Downloading: {url}', DownloadWarning)

projection = ccrs.Orthographic(central_longitude=-90, central_latitude=41)
pc = uxds["psi"].to_polycollection(projection=projection, override=True)
pc.set_antialiased(False)
pc.set_cmap("plasma")
fig, ax = plt.subplots(
1,
1,
figsize=(10, 5),
facecolor="w",
constrained_layout=True,
subplot_kw=dict(projection=projection),
)
ax.add_feature(cfeature.COASTLINE)
ax.add_feature(cfeature.BORDERS)
ax.add_collection(pc)
ax.set_global()
plt.title("PolyCollection Plot with Projection & Features")
Text(0.5, 1.0, 'PolyCollection Plot with Projection & Features')

Handling Periodic Elements#
Global Data#
If your grid contains elements that cross the antimeridian, plotting them without any corrections will lead to artifacts, as can be observed in the first plot below.
UXarray provides two ways of handling these elements:
Exclusion: Periodic Elements will be excluded from the plot, with no other corrections being done, indicated by setting
periodic_elements='exclude'
, this is the default.Splitting: Each periodic element is split into two across the antimeridian, indicated by setting
periodic_elements='split'
Ignore: Periodic Elements will be included in the plot, without any processing done to them, indicated by setting
periodic_elements='ignore'
Warning
Setting periodic_elements='split'
will lead to roughly a 20 times perfromance hit compared to the other method, so it is suggested to only use this option for small grids.
methods = ["ignore", "exclude", "split"]
poly_collections = [
uxds["psi"].to_polycollection(periodic_elements=method) for method in methods
]
fig, axes = plt.subplots(
nrows=3, figsize=(20, 10), subplot_kw={"projection": ccrs.PlateCarree()}
)
for ax, pc, method in zip(axes, poly_collections, methods):
pc.set_linewidth(0)
pc.set_cmap("plasma")
ax.set_xlim((-180, 180))
pc.set_antialiased(False)
ax.set_ylim((-90, 90))
ax.add_collection(pc)
ax.set_title(f"periodic_elements='{method}'")

projection = ccrs.Orthographic(central_longitude=-180, central_latitude=-41)
# collection with split polygons, will be much slower
pc_split = uxds["psi"].to_polycollection(periodic_elements="split")
# collection with excluded periodic polygons with explicit projection
pc_exclude = uxds["psi"].to_polycollection(
periodic_elements="exclude", projection=projection
)
pc_split.set_antialiased(False)
pc_split.set_cmap("plasma")
pc_exclude.set_antialiased(False)
pc_exclude.set_cmap("plasma")
fig, axes = plt.subplots(
1,
2,
figsize=(10, 5),
constrained_layout=True,
subplot_kw=dict(projection=projection),
)
ax1, ax2 = axes
ax1.add_feature(cfeature.COASTLINE)
ax1.add_feature(cfeature.BORDERS)
ax1.add_collection(pc_split)
ax1.set_global()
ax1.set_title("Split Polygons (Projected with Matplotlib)")
ax2.add_feature(cfeature.COASTLINE)
ax2.add_feature(cfeature.BORDERS)
ax2.add_collection(pc_exclude)
ax2.set_global()
ax2.set_title("Excluded Polygons (Explicit Projection)")
Text(0.5, 1.0, 'Excluded Polygons (Explicit Projection)')

Regional Data#
If you grid doesn’t contain any periodic elements, it is always suggested to keep periodic_elements='ignore'
for the best performance, as there is no difference in the resulting plots.
methods = ["ignore", "exclude", "split"]
poly_collections = [
uxds["psi"]
.subset.bounding_circle((0, 0), 20)
.to_polycollection(periodic_elements=method)
for method in methods
]
fig, axes = plt.subplots(
nrows=3, figsize=(10, 10), subplot_kw={"projection": ccrs.PlateCarree()}
)
for ax, pc, method in zip(axes, poly_collections, methods):
pc.set_linewidth(0)
pc.set_cmap("plasma")
pc.set_antialiased(False)
ax.set_xlim((-20, 20))
ax.set_ylim((-20, 20))
ax.add_collection(pc)
ax.set_title(f"periodic_elements='{method}'")
/home/docs/checkouts/readthedocs.org/user_builds/uxarray/conda/latest/lib/python3.13/site-packages/uxarray/grid/coordinates.py:255: UserWarning: This cannot be guaranteed to work correctly on concave polygons
warnings.warn("This cannot be guaranteed to work correctly on concave polygons")
