Diagnostic Quantities
Pointwise and domain-integrated diagnostic quantities on Arakawa C-grids.
Kinetic Energy
finitevolx.kinetic_energy(u, v, mask=None)
Kinetic energy at T-points (cell centers) on an Arakawa C-grid.
Eq: ke[j, i] = 0.5 * (u²_on_T[j, i] + v²_on_T[j, i])
where u² and v² are averaged from face-points to T-points: u²_on_T[j, i] = 0.5 * (u[j, i+1/2]² + u[j, i-1/2]²) = 0.5 * (u[j, i]² + u[j, i-1]²) v²_on_T[j, i] = 0.5 * (v[j+1/2, i]² + v[j-1/2, i]²) = 0.5 * (v[j, i]² + v[j-1, i]²)
Args: u (Array): x-velocity at U-points (east faces), shape [Ny, Nx]. v (Array): y-velocity at V-points (north faces), shape [Ny, Nx]. mask (Array | None): optional binary mask at T-points. If provided, KE is zeroed where mask is 0.
Returns: ke (Array): kinetic energy at T-points, shape [Ny, Nx]. Ghost ring is zero; interior is [1:-1, 1:-1].
Source code in finitevolx/_src/operators/diagnostics.py
Bernoulli Potential
finitevolx.bernoulli_potential(h, u, v, gravity=GRAVITY)
Bernoulli potential at T-points on an Arakawa C-grid.
Eq: p[j, i] = ke[j, i] + g * h[j, i]
where ke is the kinetic energy at T-points.
Args: h (Array): layer thickness at T-points, shape [Ny, Nx]. u (Array): x-velocity at U-points (east faces), shape [Ny, Nx]. v (Array): y-velocity at V-points (north faces), shape [Ny, Nx]. gravity (float): gravitational acceleration. Default = 9.81.
Returns: p (Array): Bernoulli potential at T-points, shape [Ny, Nx]. Ghost ring is zero; interior is [1:-1, 1:-1].
Example: >>> u, v, h = ... >>> p = bernoulli_potential(h=h, u=u, v=v)
Source code in finitevolx/_src/operators/diagnostics.py
Relative Vorticity
finitevolx.relative_vorticity_cgrid(u, v, dx, dy)
Relative vorticity at X-points (corners) on an Arakawa C-grid.
zeta[j+1/2, i+1/2] = (v[j+1/2, i+1] - v[j+1/2, i]) / dx - (u[j+1, i+1/2] - u[j, i+1/2]) / dy
This is a standalone functional form of the class-based
:meth:Vorticity2D.relative_vorticity. Both share the same
underlying implementation
(:func:finitevolx._src.operators.difference._curl_2d).
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
u
|
Float[Array, 'Ny Nx']
|
x-velocity at U-points. |
required |
v
|
Float[Array, 'Ny Nx']
|
y-velocity at V-points. |
required |
dx
|
float
|
Grid spacing in x. |
required |
dy
|
float
|
Grid spacing in y. |
required |
Returns:
| Type | Description |
|---|---|
Float[Array, 'Ny Nx']
|
Relative vorticity at X-points. Ghost ring is zero. |
Source code in finitevolx/_src/operators/diagnostics.py
Potential Vorticity
finitevolx.potential_vorticity(omega, f, h)
Potential vorticity (pointwise).
q = (omega + f) / h
All inputs must live on the same grid points (typically X-points
after interpolating f and h from T-points).
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
omega
|
Float[Array, 'Ny Nx']
|
Relative vorticity. |
required |
f
|
Float[Array, 'Ny Nx']
|
Coriolis parameter. |
required |
h
|
Float[Array, 'Ny Nx']
|
Layer thickness (or depth). |
required |
Returns:
| Type | Description |
|---|---|
Float[Array, 'Ny Nx']
|
Potential vorticity. Where |
Source code in finitevolx/_src/operators/diagnostics.py
QG Potential Vorticity
finitevolx.qg_potential_vorticity(psi, f0, beta, dx, dy, y, y0)
QG potential vorticity at T-points (single layer).
q = nabla^2(psi) / f0 + beta * (y - y0) / f0
This computes the per-layer (barotropic) QG PV for a single [Ny, Nx]
streamfunction. For multi-layer QG, lift this function with
:func:~finitevolx.multilayer and subtract the cross-layer stretching
term separately::
from finitevolx import multilayer, qg_potential_vorticity
# Per-layer PV (vmapped over the layer axis)
q = multilayer(lambda p: qg_potential_vorticity(p, f0, beta, dx, dy, y, y0))(
psi
) # psi: [nl, Ny, Nx] -> q: [nl, Ny, Nx]
# Cross-layer stretching (A couples layers, cannot be vmapped)
q = q - stretching_term(A, psi)
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
psi
|
Float[Array, 'Ny Nx']
|
Streamfunction at T-points. |
required |
f0
|
float
|
Reference Coriolis parameter. |
required |
beta
|
float
|
Meridional gradient of the Coriolis parameter. |
required |
dx
|
float
|
Grid spacing in x. |
required |
dy
|
float
|
Grid spacing in y. |
required |
y
|
Float[Array, 'Ny Nx']
|
Meridional coordinate at T-points. |
required |
y0
|
float
|
Reference latitude. |
required |
Returns:
| Type | Description |
|---|---|
Float[Array, 'Ny Nx']
|
QG potential vorticity at T-points.
Ghost ring is zero; interior is |
Source code in finitevolx/_src/operators/diagnostics.py
Stretching Term (Multi-Layer QG)
finitevolx.stretching_term(A, psi)
Cross-layer stretching term for multi-layer QG models.
Computes A . psi — the vertical coupling that mixes information
across layers. This is a cross-layer operation and cannot be vmapped;
it should be used alongside :func:qg_potential_vorticity (which
handles the per-layer Laplacian + beta terms via
:func:~finitevolx.multilayer).
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
A
|
Float[Array, 'nl nl']
|
Coupling (stretching) matrix. Shape |
required |
psi
|
Float[Array, 'nl Ny Nx']
|
Streamfunction at T-points for all layers. |
required |
Returns:
| Type | Description |
|---|---|
Float[Array, 'nl Ny Nx']
|
Stretching contribution, same shape as |
Source code in finitevolx/_src/operators/diagnostics.py
Strain
finitevolx.shear_strain(u, v, dx, dy)
Shear strain at X-points (corners) on an Arakawa C-grid.
Ss[j+1/2, i+1/2] = dv/dx + du/dy = (v[j+1/2, i+1] - v[j+1/2, i]) / dx + (u[j+1, i+1/2] - u[j, i+1/2]) / dy
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
u
|
Float[Array, 'Ny Nx']
|
x-velocity at U-points. |
required |
v
|
Float[Array, 'Ny Nx']
|
y-velocity at V-points. |
required |
dx
|
float
|
Grid spacing in x. |
required |
dy
|
float
|
Grid spacing in y. |
required |
Returns:
| Type | Description |
|---|---|
Float[Array, 'Ny Nx']
|
Shear strain at X-points. Ghost ring is zero. |
Source code in finitevolx/_src/operators/diagnostics.py
finitevolx.tensor_strain(u, v, dx, dy)
Tensor (normal) strain at T-points on an Arakawa C-grid.
Sn[j, i] = du/dx - dv/dy = (u[j, i+1/2] - u[j, i-1/2]) / dx - (v[j+1/2, i] - v[j-1/2, i]) / dy
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
u
|
Float[Array, 'Ny Nx']
|
x-velocity at U-points. |
required |
v
|
Float[Array, 'Ny Nx']
|
y-velocity at V-points. |
required |
dx
|
float
|
Grid spacing in x. |
required |
dy
|
float
|
Grid spacing in y. |
required |
Returns:
| Type | Description |
|---|---|
Float[Array, 'Ny Nx']
|
Tensor strain at T-points. Ghost ring is zero. |
Source code in finitevolx/_src/operators/diagnostics.py
finitevolx.strain_magnitude_squared(sn, ss)
Squared strain magnitude (pointwise).
sigma^2 = Sn^2 + Ss^2
Inputs must live on the same grid points. If sn is at T-points
and ss at X-points, interpolate one to the other before calling.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
sn
|
Float[Array, 'Ny Nx']
|
Tensor (normal) strain. |
required |
ss
|
Float[Array, 'Ny Nx']
|
Shear strain. |
required |
Returns:
| Type | Description |
|---|---|
Float[Array, 'Ny Nx']
|
Squared strain magnitude. |
Source code in finitevolx/_src/operators/diagnostics.py
Okubo-Weiss Parameter
finitevolx.okubo_weiss(sn, ss, omega)
Okubo-Weiss parameter (pointwise).
OW = Sn^2 + Ss^2 - omega^2
All inputs must live on the same grid points.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
sn
|
Float[Array, 'Ny Nx']
|
Tensor (normal) strain. |
required |
ss
|
Float[Array, 'Ny Nx']
|
Shear strain. |
required |
omega
|
Float[Array, 'Ny Nx']
|
Relative vorticity. |
required |
Returns:
| Type | Description |
|---|---|
Float[Array, 'Ny Nx']
|
Okubo-Weiss parameter. Positive = strain-dominated, negative = vorticity-dominated. |
Source code in finitevolx/_src/operators/diagnostics.py
Enstrophy
finitevolx.enstrophy(omega, mask=None)
Enstrophy (pointwise).
Z = 0.5 * omega^2
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
omega
|
Float[Array, 'Ny Nx']
|
Relative vorticity (typically at X-points). |
required |
mask
|
Float[Array, 'Ny Nx'] | None
|
Binary mask. If provided, result is zeroed where mask is 0. |
None
|
Returns:
| Type | Description |
|---|---|
Float[Array, 'Ny Nx']
|
Enstrophy at the same grid points as omega. |
Source code in finitevolx/_src/operators/diagnostics.py
finitevolx.potential_enstrophy(q, h, mask=None)
Potential enstrophy (pointwise).
PE = 0.5 * q^2 * h
Inputs must live on the same grid points.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
q
|
Float[Array, 'Ny Nx']
|
Potential vorticity. |
required |
h
|
Float[Array, 'Ny Nx']
|
Layer thickness (or depth). |
required |
mask
|
Float[Array, 'Ny Nx'] | None
|
Binary mask. If provided, result is zeroed where mask is 0. |
None
|
Returns:
| Type | Description |
|---|---|
Float[Array, 'Ny Nx']
|
Potential enstrophy. |
Source code in finitevolx/_src/operators/diagnostics.py
Available Potential Energy
finitevolx.available_potential_energy(h, H, g_prime)
Available potential energy at T-points.
APE = 0.5 * g' * (h - H)^2
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
h
|
Float[Array, 'Ny Nx']
|
Layer thickness at T-points. |
required |
H
|
Float[Array, 'Ny Nx']
|
Mean (reference) layer thickness at T-points. |
required |
g_prime
|
float
|
Reduced gravity. |
required |
Returns:
| Type | Description |
|---|---|
Float[Array, 'Ny Nx']
|
Available potential energy at T-points. |
Source code in finitevolx/_src/operators/diagnostics.py
Domain-Integrated Quantities
finitevolx.total_energy(ke, ape, dx, dy)
Domain-integrated total energy.
E = sum_{interior}(KE + APE) * dx * dy
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
ke
|
Float[Array, 'Ny Nx']
|
Kinetic energy at T-points. |
required |
ape
|
Float[Array, 'Ny Nx']
|
Available potential energy at T-points. |
required |
dx
|
float
|
Grid spacing in x. |
required |
dy
|
float
|
Grid spacing in y. |
required |
Returns:
| Type | Description |
|---|---|
Float[Array, '']
|
Scalar total energy. |
Source code in finitevolx/_src/operators/diagnostics.py
finitevolx.total_enstrophy(ens, dx, dy)
Domain-integrated enstrophy.
Z = sum_{interior}(0.5 * omega^2) * dx * dy
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
ens
|
Float[Array, 'Ny Nx']
|
Enstrophy field (e.g. from :func: |
required |
dx
|
float
|
Grid spacing in x. |
required |
dy
|
float
|
Grid spacing in y. |
required |
Returns:
| Type | Description |
|---|---|
Float[Array, '']
|
Scalar total enstrophy. |
Source code in finitevolx/_src/operators/diagnostics.py
Vertical Velocity
finitevolx.vertical_velocity(u, v, dx, dy, dz, mask=None)
Vertical velocity from the continuity equation.
Integrates the horizontal divergence from bottom (w=0) to top:
w[k+1/2] = w[k-1/2] - (du/dx + dv/dy)[k] * dz
where the horizontal divergence is computed at T-points for each z-level.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
u
|
Float[Array, 'Nz Ny Nx']
|
x-velocity at U-points, with ghost cells in all three dimensions. |
required |
v
|
Float[Array, 'Nz Ny Nx']
|
y-velocity at V-points, with ghost cells in all three dimensions. |
required |
dx
|
float
|
Grid spacing in x. |
required |
dy
|
float
|
Grid spacing in y. |
required |
dz
|
float
|
Grid spacing in z. |
required |
mask
|
Float[Array, 'Nz Ny Nx'] | None
|
Binary mask at T-points. If provided, the horizontal divergence is zeroed where mask is 0 before vertical integration. |
None
|
Returns:
| Type | Description |
|---|---|
Float[Array, 'Nzp1 Ny Nx']
|
Vertical velocity at W-points (cell interfaces in z).
Shape is |