Vector Calculus#
Computing and working with vector calculus on unstructured grids is essential for analyzing scalar and vector fields in geoscience. This notebook demonstrates UXarray’s finite-volume implementations and verifies expected identities.
We will showcase:
Gradient of face-centered scalar fields (zonal and meridional components)
Curl of vector fields (including curl of a gradient ≈ 0 and a synthetic vortex)
Divergence of vector fields (including Laplacian as ∇·∇φ, divergence of a vortex ≈ 0, and radial expansion > 0)
Scalar dot gradient (advection term v · ∇q)
Vector calculus identities verification
import holoviews as hv
import numpy as np
import uxarray as ux
hv.extension("bokeh")
Data#
This notebook uses a subset of a 30km MPAS atmosphere grid, taken centered at 45 degrees longitude and 0 degrees latitude with a radius of 2 degrees.
face_lon: Longitude at cell-centersface_lat: Latitude at cell-centersgaussian: Gaussian initialized at the center of the gridinverse_gaussian: Inverse of the gaussian above.
uxds = ux.tutorial.open_dataset("mpas-dyamond-30km-gradient")
uxds
<xarray.UxDataset> Size: 3kB
Dimensions: (n_face: 195)
Dimensions without coordinates: n_face
Data variables:
face_lat (n_face) float32 780B ...
face_lon (n_face) float32 780B ...
gaussian (n_face) float32 780B ...
inverse_gaussian (n_face) float32 780B ...1. Gradient (∇φ)#
Background#
The gradient of a scalar field φ gives both the direction of steepest increase and the rate of change in that direction. On unstructured grids, we use the Green–Gauss theorem:
Implementation#
In a finite-volume context, the gradient of a scalar field \(\phi\) is obtained by summing fluxes across each cell face and dividing by the cell’s volume.
Input |
Usage |
Output |
|---|---|---|
Scalar field \(\phi\) |
|
Vector field \(\nabla\phi\) |
Finite-volume discretization#
Discrete gradient at cell center \(C^*\)#
Usage#
Gradients can be computed using the UxDataArray.gradient() method on a face-centered data variable.
Units. By default,
gradient()andcurl()divide byuxgrid.sphere_radiusso derivatives carry physical units (e.g.[data units]/m,1/sfor velocity curl). Passscale_by_radius=Falseto keep results on the unit sphere (per radian). If the grid has nosphere_radiusattribute, the call falls back to unit-sphere output and emits aUserWarning.
grad_lat = uxds["face_lat"].gradient()
grad_lon = uxds["face_lon"].gradient()
grad_gauss = uxds["gaussian"].gradient()
grad_inv_gauss = uxds["inverse_gaussian"].gradient()
Examining one of the outputs, we find that the zonal_gradient and meridional_gradient data variables store the rate of change along longitude (east–west) and latitude (north–south), respectively.
print("Gradient components:")
print(f"Zonal gradient shape: {grad_gauss.zonal_gradient.shape}")
print(f"Meridional gradient shape: {grad_gauss.meridional_gradient.shape}")
grad_gauss
Gradient components:
Zonal gradient shape: (195,)
Meridional gradient shape: (195,)
<xarray.UxDataset> Size: 3kB
Dimensions: (n_face: 195)
Dimensions without coordinates: n_face
Data variables:
zonal_gradient (n_face) float64 2kB nan nan nan nan ... nan nan nan
meridional_gradient (n_face) float64 2kB nan nan nan nan ... nan nan nan
Attributes:
gradient: TruePlotting Gradients#
To visualize gradients, we represent them as vector fields and overlay them on the original scalar data.
def plot_gradient_vectors(uxda_grad, **kwargs):
"""
Plots gradient vectors using HoloViews
"""
uxgrid = uxda_grad.uxgrid
mag = np.hypot(uxda_grad.zonal_gradient, uxda_grad.meridional_gradient)
angle = np.arctan2(uxda_grad.meridional_gradient, uxda_grad.zonal_gradient)
return hv.VectorField(
(uxgrid.face_lon, uxgrid.face_lat, angle, mag), **kwargs
).opts(magnitude="Magnitude")
# Overlay the gradient vector field on top of the original data variable
p1 = (
uxds["face_lat"].plot(cmap="Oranges", aspect=1) * plot_gradient_vectors(grad_lat)
).opts(title="∇ Cell Latitudes")
p2 = (
uxds["face_lon"].plot(cmap="Oranges", aspect=1) * plot_gradient_vectors(grad_lon)
).opts(title="∇ Cell Longitudes")
p3 = (
uxds["gaussian"].plot(cmap="Oranges", aspect=1) * plot_gradient_vectors(grad_gauss)
).opts(title="∇ Gaussian")
p4 = (
uxds["inverse_gaussian"].plot(cmap="Oranges", aspect=1)
* plot_gradient_vectors(grad_inv_gauss)
).opts(title="∇ Inverse Gaussian")
# Compose all four plots in a 2 column layout
(p1 + p2 + p3 + p4).cols(2).opts(shared_axes=False)
2. Curl (∇ × F)#
Background#
The curl of a vector field F = (u, v) measures the local rotation or circulation. In 2D, curl produces a scalar field representing the magnitude of rotation:
Positive curl: Counter-clockwise rotation
Negative curl: Clockwise rotation
Zero curl: No local rotation (irrotational flow)
Usage#
Curl can be computed using the UxDataArray.curl() method on vector field components:
Input |
Usage |
Output |
|---|---|---|
Vector field (u, v) |
|
Scalar curl field |
Constant Fields (Mathematical Validation)#
The curl of a constant vector field should be zero everywhere (within numerical precision).
# Constant vector field: u=1, v=2 (face-centered)
u_constant = uxds["face_lat"] * 0 + 1.0
v_constant = uxds["face_lat"] * 0 + 2.0
# Compute partials via gradient
grad_u = u_constant.gradient()
grad_v = v_constant.gradient()
du_dy = grad_u["meridional_gradient"]
dv_dx = grad_v["zonal_gradient"]
curl_constant = dv_dx - du_dy
finite = np.isfinite(curl_constant.values)
vals = curl_constant.values[finite]
print(
f"Total faces: {curl_constant.size}, interior: {vals.size}, boundary NaNs: {np.isnan(curl_constant.values).sum()}"
)
if vals.size:
print(f"Finite curl range: [{vals.min():.2e}, {vals.max():.2e}]")
print(
f"Max |curl|: {np.abs(vals).max():.2e}, Mean |curl|: {np.abs(vals).mean():.2e}"
)
Total faces: 195, interior: 147, boundary NaNs: 48
Finite curl range: [0.00e+00, 0.00e+00]
Max |curl|: 0.00e+00, Mean |curl|: 0.00e+00
Gaussian Fields#
# Use Gaussian fields as vector components
u_gauss = uxds["gaussian"]
v_gauss = uxds["inverse_gaussian"]
# Compute partials via gradient
grad_u = u_gauss.gradient()
grad_v = v_gauss.gradient()
du_dy = grad_u["meridional_gradient"]
dv_dx = grad_v["zonal_gradient"]
curl_gauss = dv_dx - du_dy
finite = np.isfinite(curl_gauss.values)
vals = curl_gauss.values[finite]
print(
f"Total faces: {curl_gauss.size}, interior: {vals.size}, boundary NaNs: {np.isnan(curl_gauss.values).sum()}"
)
if vals.size:
print(f"Finite curl range: [{vals.min():.6f}, {vals.max():.6f}]")
print(f"Mean curl (finite): {vals.mean():.6f}")
Total faces: 195, interior: 147, boundary NaNs: 48
Finite curl range: [-0.000022, 0.000022]
Mean curl (finite): -0.000000
Example 1: Curl of Gradient Field#
In the continuous setting, curl(∇φ) = 0 exactly for any scalar field φ. On an unstructured mesh, however, gradient and curl are independent finite-volume stencils that do not form a discrete de Rham complex, so the identity holds only approximately. The residual is a discretization error — not a bug — and shrinks with grid refinement (O(Δx) for first-order stencils). For a Gaussian field on a coarse mesh the residual is typically O(1–10), consistent with the grid spacing.
# Extract gradient components
u_component = grad_gauss.zonal_gradient
v_component = grad_gauss.meridional_gradient
# Compute partial derivatives via gradient()
grad_u = u_component.gradient()
grad_v = v_component.gradient()
du_dy = grad_u["meridional_gradient"]
dv_dx = grad_v["zonal_gradient"]
# Curl = ∂v/∂x - ∂u/∂y
curl_of_gradient = dv_dx - du_dy
print(
f"Curl of gradient range: [{curl_of_gradient.min().values:.2e}, {curl_of_gradient.max().values:.2e}]"
)
print(f"Mean absolute curl: {abs(curl_of_gradient).mean().values:.2e}")
curl_plot = curl_of_gradient.plot(cmap="RdBu_r", aspect=1, clim=(-1e-10, 1e-10)).opts(
title="Curl of Gradient Field (Should ≈ 0)", colorbar=True
)
curl_plot
Curl of gradient range: [-1.11e-13, 1.07e-13]
Mean absolute curl: 4.35e-14
Example 2: Synthetic Vortex Field#
To better demonstrate curl, let’s create a synthetic rotating vector field (vortex):
# Get face coordinates
face_lon = uxds.uxgrid.face_lon.values
face_lat = uxds.uxgrid.face_lat.values
# Center the coordinates
lon_center = np.mean(face_lon)
lat_center = np.mean(face_lat)
x = face_lon - lon_center
y = face_lat - lat_center
# Create a vortex: u = -y, v = x (pure rotation)
u_vortex_data = -y
v_vortex_data = x
# Create UxDataArrays
u_vortex = ux.UxDataArray(
u_vortex_data, dims=["n_face"], uxgrid=uxds.uxgrid, name="u_vortex"
)
v_vortex = ux.UxDataArray(
v_vortex_data, dims=["n_face"], uxgrid=uxds.uxgrid, name="v_vortex"
)
# Compute curl via gradients
grad_u = u_vortex.gradient()
grad_v = v_vortex.gradient()
du_dy = grad_u["meridional_gradient"]
dv_dx = grad_v["zonal_gradient"]
curl_vortex = dv_dx - du_dy
print(
f"Vortex curl range: [{curl_vortex.min().values:.2f}, {curl_vortex.max().values:.2f}]"
)
print("Expected curl for pure rotation: ~2.0")
Vortex curl range: [0.00, 0.00]
Expected curl for pure rotation: ~2.0
# Plot the vortex vector field and its curl
def plot_vector_field(u, v, **kwargs):
"""Plot vector field using HoloViews"""
uxgrid = u.uxgrid
mag = np.hypot(u, v)
angle = np.arctan2(v, u)
return hv.VectorField(
(uxgrid.face_lon, uxgrid.face_lat, angle, mag), **kwargs
).opts(magnitude="Magnitude")
# Vector field plot
vortex_vectors = plot_vector_field(u_vortex, v_vortex).opts(
title="Synthetic Vortex Field", aspect=1, arrow_heads=True
)
# Curl magnitude plot
curl_magnitude = curl_vortex.plot(cmap="Reds", aspect=1).opts(
title="Curl of Vortex Field", colorbar=True
)
(vortex_vectors + curl_magnitude).opts(shared_axes=False)
Example 3: Relative Vorticity from Solid-Body Rotation#
The 2D curl on the sphere is relative vorticity \(\zeta = \partial v/\partial x - \partial u/\partial y\), a quantity meteorologists and oceanographers care about every day. With the default scale_by_radius=True, u.curl(v) returns \(\zeta\) in physical units of \(s^{-1}\).
A clean analytical check is solid-body rotation about the polar axis with angular speed \(\Omega\):
For this flow the relative vorticity on a sphere of radius \(R\) is
We construct the field on the existing MPAS subset and compare u.curl(v) against this closed form.
OMEGA = 7.292115e-5 # Earth's angular speed (rad/s)
R = uxds.uxgrid.sphere_radius
lat_rad = np.deg2rad(uxds.uxgrid.face_lat.values)
u_sbr_data = OMEGA * R * np.cos(lat_rad)
v_sbr_data = np.zeros_like(u_sbr_data)
u_sbr = ux.UxDataArray(
u_sbr_data,
dims=["n_face"],
uxgrid=uxds.uxgrid,
name="u_sbr",
attrs={"units": "m/s"},
)
v_sbr = ux.UxDataArray(
v_sbr_data,
dims=["n_face"],
uxgrid=uxds.uxgrid,
name="v_sbr",
attrs={"units": "m/s"},
)
zeta = u_sbr.curl(v_sbr)
zeta_analytic = -2.0 * OMEGA * np.sin(lat_rad)
finite = np.isfinite(zeta.values)
err = np.abs(zeta.values[finite] - zeta_analytic[finite])
print(f"curl units: {zeta.attrs['units']}")
print(f"sphere_radius (m): {R:.3e}")
print(
f"computed zeta range: [{zeta.values[finite].min():.3e}, {zeta.values[finite].max():.3e}] 1/s"
)
print(
f"analytic zeta range: [{zeta_analytic.min():.3e}, {zeta_analytic.max():.3e}] 1/s"
)
print(f"max |error|: {err.max():.3e} 1/s (≈ {err.max() / (2 * OMEGA):.1%} of 2Omega)")
curl units: (m/s)/m
sphere_radius (m): 6.371e+06
computed zeta range: [-6.673e-06, 6.388e-06] 1/s
analytic zeta range: [-4.844e-06, 5.039e-06] 1/s
max |error|: 1.112e-05 1/s (≈ 7.6% of 2Omega)
zeta.plot(cmap="RdBu_r", aspect=1).opts(
title="Relative Vorticity (1/s) — Solid-Body Rotation", colorbar=True
)
3. Divergence (∇ · F)#
Background#
The divergence of a vector field F = (u, v) measures the local expansion or contraction of the field:
Positive divergence: Expansion (source)
Negative divergence: Contraction (sink)
Zero divergence: Incompressible flow
Usage#
Divergence can be computed using the UxDataArray.divergence() method:
Input |
Usage |
Output |
|---|---|---|
Vector field (u, v) |
|
Scalar divergence field |
Constant Fields (Mathematical Validation)#
The divergence of a constant vector field should be zero everywhere (within numerical precision).
# Constant vector field: u=1, v=2 (face-centered)
u_constant = uxds["face_lat"] * 0 + 1.0
v_constant = uxds["face_lat"] * 0 + 2.0
# Compute partials via gradient
grad_u = u_constant.gradient()
grad_v = v_constant.gradient()
du_dx = grad_u["zonal_gradient"]
dv_dy = grad_v["meridional_gradient"]
div_constant = du_dx + dv_dy
finite = np.isfinite(div_constant.values)
vals = div_constant.values[finite]
print(
f"Total faces: {div_constant.size}, interior: {vals.size}, boundary NaNs: {np.isnan(div_constant.values).sum()}"
)
if vals.size:
print(f"Finite divergence range: [{vals.min():.2e}, {vals.max():.2e}]")
print(f"Max |div|: {np.abs(vals).max():.2e}, Mean |div|: {np.abs(vals).mean():.2e}")
Total faces: 195, interior: 147, boundary NaNs: 48
Finite divergence range: [0.00e+00, 0.00e+00]
Max |div|: 0.00e+00, Mean |div|: 0.00e+00
Gaussian Fields#
# Use Gaussian fields as vector components
u_gauss = uxds["gaussian"]
v_gauss = uxds["inverse_gaussian"]
# Compute partials via gradient
grad_u = u_gauss.gradient()
grad_v = v_gauss.gradient()
du_dx = grad_u["zonal_gradient"]
dv_dy = grad_v["meridional_gradient"]
div_gauss = du_dx + dv_dy
finite = np.isfinite(div_gauss.values)
vals = div_gauss.values[finite]
print(
f"Total faces: {div_gauss.size}, interior: {vals.size}, boundary NaNs: {np.isnan(div_gauss.values).sum()}"
)
if vals.size:
print(f"Finite divergence range: [{vals.min():.6f}, {vals.max():.6f}]")
print(f"Mean divergence (finite): {vals.mean():.6f}")
Total faces: 195, interior: 147, boundary NaNs: 48
Finite divergence range: [-0.000022, 0.000022]
Mean divergence (finite): 0.000000
Example 1: Divergence of Gradient Field#
# Compute divergence of the gradient (Laplacian) via partials
# ∇²φ = ∂(∇φ_x)/∂x + ∂(∇φ_y)/∂y
grad_gauss_u = grad_gauss["zonal_gradient"]
grad_gauss_v = grad_gauss["meridional_gradient"]
gxu = grad_gauss_u.gradient()["zonal_gradient"] # ∂u/∂x
gyv = grad_gauss_v.gradient()["meridional_gradient"] # ∂v/∂y
div_of_gradient = gxu + gyv
print(
f"Divergence of gradient range: [{div_of_gradient.min().values:.2e}, {div_of_gradient.max().values:.2e}]"
)
print("This is the Laplacian (∇²) of the original gaussian field")
div_plot = div_of_gradient.plot(cmap="RdBu_r", aspect=1).opts(
title="Divergence of Gradient (Laplacian)", colorbar=True
)
div_plot
Divergence of gradient range: [-1.34e-09, 2.53e-11]
This is the Laplacian (∇²) of the original gaussian field
Example 2: Divergence of Vortex Field#
Pure rotation should have zero divergence (incompressible):
# Compute divergence of the vortex via gradients: div = ∂u/∂x + ∂v/∂y
grad_u = u_vortex.gradient()
grad_v = v_vortex.gradient()
du_dx = grad_u["zonal_gradient"]
dv_dy = grad_v["meridional_gradient"]
div_vortex = du_dx + dv_dy
print(
f"Vortex divergence range: [{div_vortex.min().values:.2e}, {div_vortex.max().values:.2e}]"
)
print(f"Mean absolute divergence: {abs(div_vortex).mean().values:.2e}")
print("Pure rotation should have zero divergence")
div_vortex_plot = div_vortex.plot(cmap="RdBu_r", aspect=1, clim=(-1e-10, 1e-10)).opts(
title="Divergence of Vortex (Should ≈ 0)", colorbar=True
)
div_vortex_plot
Vortex divergence range: [-8.50e-10, 7.54e-10]
Mean absolute divergence: 2.71e-10
Pure rotation should have zero divergence
Example 3: Radial Expansion Field#
Let’s create a field that expands radially outward to demonstrate positive divergence:
# Create radial expansion field: u = x, v = y
u_radial_data = x
v_radial_data = y
u_radial = ux.UxDataArray(
u_radial_data, dims=["n_face"], uxgrid=uxds.uxgrid, name="u_radial"
)
v_radial = ux.UxDataArray(
v_radial_data, dims=["n_face"], uxgrid=uxds.uxgrid, name="v_radial"
)
# Compute curl and divergence via gradients
grad_u = u_radial.gradient()
grad_v = v_radial.gradient()
du_dy = grad_u["meridional_gradient"]
dv_dx = grad_v["zonal_gradient"]
du_dx = grad_u["zonal_gradient"]
dv_dy = grad_v["meridional_gradient"]
curl_radial = dv_dx - du_dy
div_radial = du_dx + dv_dy
print(
f"Radial field curl range: [{curl_radial.min().values:.2e}, {curl_radial.max().values:.2e}]"
)
print(
f"Radial field divergence range: [{div_radial.min().values:.2f}, {div_radial.max().values:.2f}]"
)
print("Expected: curl ≈ 0, divergence ≈ 2 (pure expansion)")
Radial field curl range: [-7.54e-10, 8.50e-10]
Radial field divergence range: [0.00, 0.00]
Expected: curl ≈ 0, divergence ≈ 2 (pure expansion)
# Plot radial expansion field and its divergence
radial_vectors = plot_vector_field(u_radial, v_radial).opts(
title="Radial Expansion Field", aspect=1, arrow_heads=True
)
radial_div = div_radial.plot(cmap="Reds", aspect=1).opts(
title="Divergence of Radial Field", colorbar=True
)
(radial_vectors + radial_div).opts(shared_axes=False)
4. Scalar Dot Gradient (v · ∇q)#
Background#
Many geophysical applications need the projection of a vector field v = (u, v) onto the gradient of a scalar field q. This quantity appears, for example, as the horizontal advection term in transport equations:
Usage#
The UxDataArray.scalardotgradient(v, q) method computes this dot product directly. The data variable it is called on (self) is treated as the zonal component u, v is the meridional component, and q is the scalar field whose gradient is taken. All three inputs must be face-centered and share the same grid and dimensions.
# Treat the Gaussian and its inverse as a vector field (u, v)
u_field = uxds["gaussian"]
v_field = uxds["inverse_gaussian"]
# Scalar field whose gradient is advected
q_field = uxds["gaussian"]
# v . grad(q) = u * dq/dx + v * dq/dy
vdotgradq = u_field.scalardotgradient(v_field, q_field)
vdotgradq
<xarray.UxDataArray 'scalar_dot_gradient' (n_face: 195)> Size: 2kB
array([ nan, nan, nan, nan,
5.63028989e-06, nan, 6.81135825e-06, 6.76469334e-06,
7.46301221e-06, 8.38914996e-06, nan, 9.28183808e-06,
9.73405536e-06, nan, 9.51194907e-06, 1.07488740e-05,
nan, 9.23667709e-06, 1.08785658e-05, 1.16340699e-05,
8.54859100e-06, 1.05414154e-05, 1.13727556e-05, nan,
9.98726160e-06, 1.06082029e-05, 1.08583973e-05, 9.18134657e-06,
6.74958334e-06, nan, 9.03085700e-06, 8.06679370e-06,
1.00024965e-05, 8.18589365e-06, nan, 9.60531999e-06,
4.50960919e-06, 4.76953606e-06, 7.94837768e-06, 5.48092557e-07,
8.38504765e-06, 9.26784188e-06, nan, 3.00196832e-06,
-3.16416704e-07, nan, 7.36219560e-06, 5.81376419e-06,
-4.10755900e-07, 6.95276139e-06, -5.06708099e-06, 6.37868147e-06,
-7.44285891e-06, -2.96034896e-07, nan, -7.89354025e-06,
-3.21274348e-06, -4.50968189e-06, 5.04660176e-06, -3.73829286e-06,
-1.28747026e-06, nan, nan, -5.09192426e-06,
6.22856714e-07, 3.80714948e-06, 5.65040700e-06, 4.40765180e-06,
8.31521845e-06, 1.22357478e-06, 2.92925673e-06, nan,
2.95949077e-06, nan, 1.42681692e-06, -2.50209042e-06,
nan, 1.36101618e-07, nan, nan,
...
9.27162913e-06, 1.04505690e-05, -5.74053716e-06, -6.72584111e-06,
-2.38791409e-06, 1.69706698e-06, -5.78712458e-06, 8.26163414e-07,
nan, nan, -7.31749309e-06, nan,
6.53408236e-06, -6.06261489e-06, -9.99628187e-06, -1.04078136e-05,
-7.22886368e-06, -7.76272850e-06, -8.92451335e-06, -9.31836072e-06,
-1.12257616e-05, -1.09655237e-05, -3.56313496e-06, -1.08229061e-05,
3.47619617e-06, -1.14263302e-05, -9.63819202e-06, 3.10740230e-06,
-9.44334196e-06, -1.11438566e-05, -1.03668144e-05, -1.00536334e-05,
-8.29278142e-06, -8.79101699e-06, -8.79634182e-06, -1.92380298e-06,
-1.81668052e-06, -1.78527112e-07, -6.30142506e-06, -3.89733318e-06,
-9.40149355e-06, -3.72690151e-06, 2.68977110e-06, -9.00885560e-06,
-6.95671190e-06, -1.05797612e-05, -7.24703875e-06, -1.17306711e-05,
-8.80584084e-06, -1.07944433e-05, -1.03301587e-05, -8.71411387e-06,
-8.83342791e-06, -9.62326758e-06, -8.66659500e-06, nan,
-5.79166917e-06, -8.03250051e-06, -8.05819981e-06, nan,
-5.67230399e-06, -6.82835905e-06, nan, -3.03963147e-06,
-5.07008262e-06, nan, -3.08346718e-06, nan,
nan, nan, nan, nan,
-8.47174992e-06, nan, -8.07671379e-06, nan,
nan, nan, nan])
Dimensions without coordinates: n_face
Attributes:
units: 1/m
long_name: scalar dot gradient
description: Dot product u * (dq/dx) + v * (dq/dy).# Verify against an explicit gradient-based computation
grad_q = q_field.gradient()
expected = u_field * grad_q["zonal_gradient"] + v_field * grad_q["meridional_gradient"]
# Boundary faces yield NaN gradients, so compare only finite values
finite = np.isfinite(vdotgradq.values) & np.isfinite(expected.values)
max_diff = np.abs(vdotgradq.values[finite] - expected.values[finite]).max()
print(f"Max difference vs. explicit u*dq/dx + v*dq/dy: {max_diff:.2e}")
print(f"Matches explicit computation: {max_diff < 1e-12}")
Max difference vs. explicit u*dq/dx + v*dq/dy: 0.00e+00
Matches explicit computation: True
# Visualize the scalar dot gradient field
vdotgradq.plot(cmap="RdBu_r", aspect=1).opts(
title="Scalar Dot Gradient v · ∇q", colorbar=True
)
5. Vector Calculus Identities#
Let’s verify some fundamental vector calculus identities using our computed fields:
Identity 1: Curl of Gradient (Discretization Residual)#
In the continuous setting: ∇ × (∇φ) = 0. UXarray’s finite-volume operators are not mimetic, so this holds only approximately. The residual below is discretization error, not a numerical bug.
# We already computed this above
max_curl_grad = np.abs(curl_of_gradient).max().values
print(f"Maximum |curl(∇φ)|: {max_curl_grad:.2e}")
print(
"Note: non-zero residual is expected discretization error from independent "
"finite-volume gradient and curl stencils (not a bug). Shrinks with grid refinement."
)
Maximum |curl(∇φ)|: 1.11e-13
Note: non-zero residual is expected discretization error from independent finite-volume gradient and curl stencils (not a bug). Shrinks with grid refinement.
Identity 2: Divergence of Curl is Zero#
For any vector field F: ∇ · (∇ × F) = 0
Note: This identity applies to 3D vector fields. In 2D, curl produces a scalar, so we can’t directly compute its divergence.
Identity 3: Properties of Special Fields#
# Summary of field properties
print("Field Properties Summary:")
print("=" * 50)
print("Gradient field:")
print(f" - Curl: {np.abs(curl_of_gradient).max().values:.2e} (≈ 0, conservative)")
print(
f" - Divergence: {div_of_gradient.min().values:.2e} to {div_of_gradient.max().values:.2e} (Laplacian)"
)
print()
print("Vortex field (pure rotation):")
print(f" - Curl: {curl_vortex.mean().values:.2f} (≈ 2, rotational)")
print(f" - Divergence: {np.abs(div_vortex).max().values:.2e} (≈ 0, incompressible)")
print()
print("Radial field (pure expansion):")
print(f" - Curl: {np.abs(curl_radial).max().values:.2e} (≈ 0, irrotational)")
print(f" - Divergence: {div_radial.mean().values:.2f} (≈ 2, expanding)")
Field Properties Summary:
==================================================
Gradient field:
- Curl: 1.11e-13 (≈ 0, conservative)
- Divergence: -1.34e-09 to 2.53e-11 (Laplacian)
Vortex field (pure rotation):
- Curl: 0.00 (≈ 2, rotational)
- Divergence: 8.50e-10 (≈ 0, incompressible)
Radial field (pure expansion):
- Curl: 8.50e-10 (≈ 0, irrotational)
- Divergence: 0.00 (≈ 2, expanding)