Momentum equation
The UMomentumEquation, VMomentumEquation, and WMomentumEquation modules provide diagnostics for every term on the right-hand side of the momentum equations. The momentum equation for the $i$-th velocity component $u_i$ is
\[\partial_t u_i = \underbrace{-\partial_j(u_j u_i)}_{\text{advection}} + \underbrace{\hat{g}_i \, b}_{\text{buoyancy}} + \underbrace{-\epsilon_{ijk} f_j u_k}_{\text{Coriolis}} + \underbrace{-\partial_i p}_{\text{pressure}} + \underbrace{-\partial_j \tau_{ij}}_{\text{viscous}} + \underbrace{(\nabla \times \mathbf{u}^S) \times \mathbf{u}}_{\text{Stokes shear}} + \underbrace{\partial_t u_i^S}_{\text{Stokes tendency}} + \underbrace{F_{u_i}}_{\text{forcing}}\]
where $\hat{g}_i$ is the $i$-th component of the gravitational unit vector, $b$ is the buoyancy, $f_j$ is the Coriolis frequency, $p$ is the pressure, $\tau_{ij}$ is the viscous/subgrid stress tensor, $\mathbf{u}^S$ is the Stokes drift, and $F_{u_i}$ is the forcing. This decomposition lets the user compute each contribution independently, build diagnostics like budget closure, or analyse the energetics of individual processes.
Each module wraps the corresponding Oceananigans velocity-tendency kernel and provides diagnostics at the natural grid location for that velocity component:
| Module | Location | Tendency wrapped |
|---|---|---|
UMomentumEquation | (Face, Center, Center) | u_velocity_tendency (NH) / hydrostatic_free_surface_u_velocity_tendency (HFS) |
VMomentumEquation | (Center, Face, Center) | v_velocity_tendency (NH) / hydrostatic_free_surface_v_velocity_tendency (HFS) |
WMomentumEquation | (Center, Center, Face) | w_velocity_tendency (NH only) |
Constructors accept either a full Oceananigans model object (the convenience form) or the individual arguments expected by the underlying kernel function. The model-only form is enough for most analysis workflows.
Example
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid; buoyancy=BuoyancyTracer(), tracers=:b, coriolis=FPlane(f=1e-4));
julia> ADV = UMomentumEquation.Advection(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: div_𝐯u (generic function with 10 methods)
└── arguments: ("Centered", "NamedTuple", "Field")
julia> BUOY = UMomentumEquation.BuoyancyAcceleration(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: x_dot_g_bᶠᶜᶜ (generic function with 10 methods)
└── arguments: ("BuoyancyForce", "NamedTuple")
julia> COR = VMomentumEquation.CoriolisAcceleration(model)
KernelFunctionOperation at (Center, Face, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: y_f_cross_U (generic function with 10 methods)
└── arguments: ("FPlane", "NamedTuple")Budget closure
The combined RHS — Tendency(model) — wraps Oceananigans' internal *_velocity_tendency kernel directly. The sum of the individual term diagnostics, taken with the correct signs from Oceananigans' source, equals Tendency to machine precision. For NonhydrostaticModel the u/v formula is
\[\text{Tendency} = -\text{Advection} + \text{BuoyancyAcceleration} -\text{CoriolisAcceleration} -\text{PressureGradient} -\text{TotalViscousDissipation} +\text{StokesShear} + \text{StokesTendency} + \text{Forcing}\]
The w-momentum equation is similar but with one subtlety: Oceananigans' w_velocity_tendency has no explicit -\partial_z p line. The vertical hydrostatic balance is treated as a property of the pressure-projection step, and whether the buoyancy term appears in the tendency depends on whether the model uses hydrostatic/nonhydrostatic pressure splitting:
hydrostatic_pressure_anomaly | Buoyancy term in Tendency | Budget formula |
|---|---|---|
Field (default, splitting on) | absorbed by $pHY'$ | -ADV - COR - TVISC + SS + ST + FORC |
nothing (splitting off) | included as $+b$ | -ADV + BUOY - COR - TVISC + SS + ST + FORC |
HydrostaticFreeSurfaceModel caveats
WMomentumEquation.Tendency,WMomentumEquation.Forcing(model, Val(:w)), and the U/V/WStokesShearandStokesTendencythrow anArgumentErroronHydrostaticFreeSurfaceModel:wis diagnosed from continuity rather than evolved by a prognostic equation, and HFS has nostokes_driftfield.Advection(::HydrostaticFreeSurfaceModel)wrapsdiv_𝐯u/div_𝐯v(the flux-form kernel). HFS's defaultVectorInvariantmomentum advection usesU_dot_∇uinstead, so to use theAdvectiondiagnostic with HFS you must build the model with a flux-form scheme, e.g.momentum_advection = Centered().- For HFS,
Tendencywrapshydrostatic_free_surface_*_velocity_tendency, which differs from the NH tendency: it has no Stokes terms, includes a barotropic free-surface pressure gradient, and a grid-slope contribution. Budget closure for HFS is not yet supported by this module — it would require additional diagnostics for those terms.
Forcing dispatch quirk
Forcing aliases KernelFunctionOperation directly (no narrowing on the kernel function), so its constructor methods are shared across the three modules. To avoid clobbering UMomentumEquation.Forcing(model) when the V and W modules are also loaded, the V/W single-argument convenience form takes an explicit Val tag:
UMomentumEquation.Forcing(model) # u-momentum forcing
VMomentumEquation.Forcing(model, Val(:v)) # v-momentum forcing
WMomentumEquation.Forcing(model, Val(:w)) # w-momentum forcingThe full explicit form Forcing(model, forcing, clock, model_fields, Val(:u/:v/:w)) is available in all three modules.
U-momentum
Oceanostics.UMomentumEquation.Advection — Type
Advection(model, u, v, w, advection_scheme; location)
Calculate the advection of u-momentum as
ADV = ∂ⱼ (uⱼ u)using Oceananigans' kernel div_𝐯u.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> ADV = UMomentumEquation.Advection(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: div_𝐯u (generic function with 10 methods)
└── arguments: ("Centered", "NamedTuple", "Field")Oceanostics.UMomentumEquation.BuoyancyAcceleration — Type
BuoyancyAcceleration(model, buoyancy, tracers; location)
Calculate the buoyancy acceleration in the x-direction as
BUOY = ĝₓ bwhere ĝₓ is the x-component of the gravitational unit vector and b is the buoyancy.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid; buoyancy=BuoyancyTracer(), tracers=:b);
julia> BUOY = UMomentumEquation.BuoyancyAcceleration(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: x_dot_g_bᶠᶜᶜ (generic function with 10 methods)
└── arguments: ("BuoyancyForce", "NamedTuple")Oceanostics.UMomentumEquation.CoriolisAcceleration — Type
CoriolisAcceleration(model, coriolis, velocities; location)
Calculate the Coriolis acceleration in the x-direction as
COR = f × uwhere f is the Coriolis parameter vector and u is the velocity vector.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid; coriolis=FPlane(f=1e-4));
julia> COR = UMomentumEquation.CoriolisAcceleration(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: x_f_cross_U (generic function with 10 methods)
└── arguments: ("FPlane", "NamedTuple")Oceanostics.UMomentumEquation.PressureGradient — Type
PressureGradient(model, hydrostatic_pressure; location)
Calculate the hydrostatic pressure gradient force in the x-direction as
PRES = ∂p/∂xwhere p is the hydrostatic pressure anomaly.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> PRES = UMomentumEquation.PressureGradient(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: hydrostatic_pressure_gradient_x (generic function with 2 methods)
└── arguments: ("Nothing",)Oceanostics.UMomentumEquation.ViscousDissipation — Type
ViscousDissipation(
model,
closure,
diffusivities,
clock,
model_fields,
buoyancy;
location
)
Calculate the viscous dissipation term (excluding immersed boundaries) as
VISC = ∂ⱼ τ₁ⱼ,where τ₁ⱼ is the viscous stress tensor for the x-momentum equation.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> VISC = UMomentumEquation.ViscousDissipation(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: ∂ⱼ_τ₁ⱼ (generic function with 10 methods)
└── arguments: ("Nothing", "Nothing", "Clock", "NamedTuple", "Nothing")Oceanostics.UMomentumEquation.ImmersedViscousDissipation — Type
ImmersedViscousDissipation(
model,
velocities,
u_immersed_bc,
closure,
diffusivities,
clock,
model_fields;
location
)
Calculate the viscous dissipation term due to immersed boundaries as
VISC = ∂ⱼ τ₁ⱼ,where τ₁ⱼ is the immersed boundary viscous stress tensor for the x-momentum equation.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> VISC = UMomentumEquation.ImmersedViscousDissipation(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: immersed_∂ⱼ_τ₁ⱼ (generic function with 2 methods)
└── arguments: ("NamedTuple", "Nothing", "Nothing", "Nothing", "Clock", "NamedTuple")Oceanostics.UMomentumEquation.TotalViscousDissipation — Type
TotalViscousDissipation(
model,
velocities,
u_immersed_bc,
closure,
diffusivities,
clock,
model_fields,
buoyancy;
location
)
Calculate the total viscous dissipation term as
VISC = ∂ⱼ τ₁ⱼ + ∂ⱼ τ₁ⱼ_immersed,where τ₁ⱼ is the interior viscous stress tensor and τ₁ⱼ_immersed is the immersed boundary viscous stress tensor for the x-momentum equation.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> VISC = UMomentumEquation.TotalViscousDissipation(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: total_∂ⱼ_τ₁ⱼ (generic function with 1 method)
└── arguments: ("NamedTuple", "Nothing", "Nothing", "Nothing", "Clock", "NamedTuple", "Nothing")Oceanostics.UMomentumEquation.StokesShear — Type
StokesShear(model, stokes_drift, velocities, time; location)
Calculate the Stokes shear term as
STOKES_SHEAR = (∇ × uˢ) × uwhere uˢ is the Stokes drift velocity and u is the velocity vector.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> STOKES = UMomentumEquation.StokesShear(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: x_curl_Uˢ_cross_U (generic function with 10 methods)
└── arguments: ("Nothing", "NamedTuple", "Float64")Oceanostics.UMomentumEquation.StokesTendency — Type
StokesTendency(model, stokes_drift, time; location)
Calculate the Stokes tendency term as
STOKES_TEND = ∂uˢ/∂twhere uˢ is the Stokes drift velocity.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> STOKES = UMomentumEquation.StokesTendency(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: ∂t_uˢ (generic function with 10 methods)
└── arguments: ("Nothing", "Float64")Oceanostics.UMomentumEquation.Forcing — Type
Forcing(
model,
forcing_func,
clock,
model_fields,
;
location
)
Calculate the forcing term Fᵘ on the x-momentum equation for model.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> FORC = UMomentumEquation.Forcing(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: Returns (generic function with 1 method)
└── arguments: ("Clock", "NamedTuple")Oceanostics.UMomentumEquation.Tendency — Type
Tendency(
model,
advection_scheme,
coriolis,
closure,
u_immersed_bc,
velocities,
free_surface,
tracers,
buoyancy,
closure_fields,
hydrostatic_pressure_anomaly,
auxiliary_fields,
vertical_coordinate,
clock,
forcing_func;
location
)
Calculate the total tendency of the u-momentum equation as computed by Oceananigans.
For NonhydrostaticModel, this includes:
- Advection: -∇⋅(𝐯u)
- Background advection terms
- Buoyancy: ĝₓ b
- Coriolis: -f × u
- Pressure gradient: -∂p/∂x
- Viscous dissipation: -∇⋅τ₁
- Immersed viscous dissipation
- Stokes shear: (∇ × uˢ) × u
- Stokes tendency: ∂uˢ/∂t
- Forcing: Fᵘ
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> TEND = UMomentumEquation.Tendency(model)
KernelFunctionOperation at (Face, Center, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: u_velocity_tendency (generic function with 1 method)
└── arguments: ("Centered", "Nothing", "Nothing", "Nothing", "Nothing", "Nothing", "Oceananigans.Models.NonhydrostaticModels.BackgroundFields", "NamedTuple", "NamedTuple", "NamedTuple", "Nothing", "Nothing", "Clock", "Returns")V-momentum
Oceanostics.VMomentumEquation.Advection — Type
Advection(model, u, v, w, advection_scheme; location)
Calculate the advection of v-momentum as
ADV = ∂ⱼ (uⱼ v)using Oceananigans' kernel div_𝐯v.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> ADV = VMomentumEquation.Advection(model)
KernelFunctionOperation at (Center, Face, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: div_𝐯v (generic function with 10 methods)
└── arguments: ("Centered", "NamedTuple", "Field")Oceanostics.VMomentumEquation.BuoyancyAcceleration — Type
BuoyancyAcceleration(model, buoyancy, tracers; location)
Calculate the buoyancy acceleration in the y-direction as
BUOY = ĝᵧ bwhere ĝᵧ is the y-component of the gravitational unit vector and b is the buoyancy.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid; buoyancy=BuoyancyTracer(), tracers=:b);
julia> BUOY = VMomentumEquation.BuoyancyAcceleration(model)
KernelFunctionOperation at (Center, Face, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: y_dot_g_bᶜᶠᶜ (generic function with 10 methods)
└── arguments: ("BuoyancyForce", "NamedTuple")Oceanostics.VMomentumEquation.CoriolisAcceleration — Type
CoriolisAcceleration(model, coriolis, velocities; location)
Calculate the Coriolis acceleration in the y-direction as
COR = (f × u)ᵧwhere f is the Coriolis parameter vector and u is the velocity vector.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid; coriolis=FPlane(f=1e-4));
julia> COR = VMomentumEquation.CoriolisAcceleration(model)
KernelFunctionOperation at (Center, Face, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: y_f_cross_U (generic function with 10 methods)
└── arguments: ("FPlane", "NamedTuple")Oceanostics.VMomentumEquation.PressureGradient — Type
PressureGradient(model, hydrostatic_pressure; location)
Calculate the hydrostatic pressure gradient force in the y-direction as
PRES = ∂p/∂ywhere p is the hydrostatic pressure anomaly.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> PRES = VMomentumEquation.PressureGradient(model)
KernelFunctionOperation at (Center, Face, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: hydrostatic_pressure_gradient_y (generic function with 2 methods)
└── arguments: ("Nothing",)Oceanostics.VMomentumEquation.ViscousDissipation — Type
ViscousDissipation(
model,
closure,
diffusivities,
clock,
model_fields,
buoyancy;
location
)
Calculate the viscous dissipation term (excluding immersed boundaries) as
VISC = ∂ⱼ τ₂ⱼ,where τ₂ⱼ is the viscous stress tensor for the y-momentum equation.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> VISC = VMomentumEquation.ViscousDissipation(model)
KernelFunctionOperation at (Center, Face, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: ∂ⱼ_τ₂ⱼ (generic function with 10 methods)
└── arguments: ("Nothing", "Nothing", "Clock", "NamedTuple", "Nothing")Oceanostics.VMomentumEquation.ImmersedViscousDissipation — Type
ImmersedViscousDissipation(
model,
velocities,
v_immersed_bc,
closure,
diffusivities,
clock,
model_fields;
location
)
Calculate the viscous dissipation term due to immersed boundaries as
VISC = ∂ⱼ τ₂ⱼ,where τ₂ⱼ is the immersed boundary viscous stress tensor for the y-momentum equation.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> VISC = VMomentumEquation.ImmersedViscousDissipation(model)
KernelFunctionOperation at (Center, Face, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: immersed_∂ⱼ_τ₂ⱼ (generic function with 2 methods)
└── arguments: ("NamedTuple", "Nothing", "Nothing", "Nothing", "Clock", "NamedTuple")Oceanostics.VMomentumEquation.TotalViscousDissipation — Type
TotalViscousDissipation(
model,
velocities,
v_immersed_bc,
closure,
diffusivities,
clock,
model_fields,
buoyancy;
location
)
Calculate the total viscous dissipation term as
VISC = ∂ⱼ τ₂ⱼ + ∂ⱼ τ₂ⱼ_immersed,where τ₂ⱼ is the interior viscous stress tensor and τ₂ⱼ_immersed is the immersed boundary viscous stress tensor for the y-momentum equation.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> VISC = VMomentumEquation.TotalViscousDissipation(model)
KernelFunctionOperation at (Center, Face, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: total_∂ⱼ_τ₂ⱼ (generic function with 1 method)
└── arguments: ("NamedTuple", "Nothing", "Nothing", "Nothing", "Clock", "NamedTuple", "Nothing")Oceanostics.VMomentumEquation.StokesShear — Type
StokesShear(model, stokes_drift, velocities, time; location)
Calculate the Stokes shear term as
STOKES_SHEAR = ((∇ × uˢ) × u)ᵧwhere uˢ is the Stokes drift velocity and u is the velocity vector.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> STOKES = VMomentumEquation.StokesShear(model)
KernelFunctionOperation at (Center, Face, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: y_curl_Uˢ_cross_U (generic function with 10 methods)
└── arguments: ("Nothing", "NamedTuple", "Float64")Oceanostics.VMomentumEquation.StokesTendency — Type
StokesTendency(model, stokes_drift, time; location)
Calculate the Stokes tendency term as
STOKES_TEND = ∂vˢ/∂twhere vˢ is the y-component of the Stokes drift velocity.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> STOKES = VMomentumEquation.StokesTendency(model)
KernelFunctionOperation at (Center, Face, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: ∂t_vˢ (generic function with 10 methods)
└── arguments: ("Nothing", "Float64")Oceanostics.VMomentumEquation.Forcing — Type
Forcing(
model,
forcing_func,
clock,
model_fields,
;
location
)
Calculate the forcing term Fᵛ on the y-momentum equation for model.
Forcing is a type alias for the generic KernelFunctionOperation (no narrowing on the kernel function), so a constructor on Forcing(model) would clash across the U/V/W modules. To disambiguate, the V-momentum convenience constructor takes an explicit Val(:v) tag.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> FORC = VMomentumEquation.Forcing(model, Val(:v))
KernelFunctionOperation at (Center, Face, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: Returns (generic function with 1 method)
└── arguments: ("Clock", "NamedTuple")Oceanostics.VMomentumEquation.Tendency — Type
Tendency(
model,
advection_scheme,
coriolis,
closure,
v_immersed_bc,
velocities,
free_surface,
tracers,
buoyancy,
closure_fields,
hydrostatic_pressure_anomaly,
auxiliary_fields,
vertical_coordinate,
clock,
forcing_func;
location
)
Calculate the total tendency of the v-momentum equation as computed by Oceananigans.
For NonhydrostaticModel, this includes:
- Advection: -∇⋅(𝐯v)
- Background advection terms
- Buoyancy: ĝᵧ b
- Coriolis: -(f × u)ᵧ
- Pressure gradient: -∂p/∂y
- Viscous dissipation: -∇⋅τ₂
- Immersed viscous dissipation
- Stokes shear: ((∇ × uˢ) × u)ᵧ
- Stokes tendency: ∂vˢ/∂t
- Forcing: Fᵛ
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> TEND = VMomentumEquation.Tendency(model)
KernelFunctionOperation at (Center, Face, Center)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: v_velocity_tendency (generic function with 1 method)
└── arguments: ("Centered", "Nothing", "Nothing", "Nothing", "Nothing", "Nothing", "Oceananigans.Models.NonhydrostaticModels.BackgroundFields", "NamedTuple", "NamedTuple", "NamedTuple", "Nothing", "Nothing", "Clock", "Returns")W-momentum
Oceanostics.WMomentumEquation.Advection — Type
Advection(model, u, v, w, advection_scheme; location)
Calculate the advection of w-momentum as
ADV = ∂ⱼ (uⱼ w)using Oceananigans' kernel div_𝐯w.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> ADV = WMomentumEquation.Advection(model)
KernelFunctionOperation at (Center, Center, Face)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: div_𝐯w (generic function with 10 methods)
└── arguments: ("Centered", "NamedTuple", "Field")Oceanostics.WMomentumEquation.BuoyancyAcceleration — Type
BuoyancyAcceleration(model, buoyancy, tracers; location)
Calculate the buoyancy acceleration in the z-direction as
BUOY = ĝ_z bwhere ĝ_z is the z-component of the gravitational unit vector and b is the buoyancy.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid; buoyancy=BuoyancyTracer(), tracers=:b);
julia> BUOY = WMomentumEquation.BuoyancyAcceleration(model)
KernelFunctionOperation at (Center, Center, Face)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: z_dot_g_bᶜᶜᶠ (generic function with 10 methods)
└── arguments: ("BuoyancyForce", "NamedTuple")Oceanostics.WMomentumEquation.CoriolisAcceleration — Type
CoriolisAcceleration(model, coriolis, velocities; location)
Calculate the Coriolis acceleration in the z-direction as
COR = (f × u)_zwhere f is the Coriolis parameter vector and u is the velocity vector.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid; coriolis=FPlane(f=1e-4));
julia> COR = WMomentumEquation.CoriolisAcceleration(model)
KernelFunctionOperation at (Center, Center, Face)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: z_f_cross_U (generic function with 10 methods)
└── arguments: ("FPlane", "NamedTuple")Oceanostics.WMomentumEquation.ViscousDissipation — Type
ViscousDissipation(
model,
closure,
diffusivities,
clock,
model_fields,
buoyancy;
location
)
Calculate the viscous dissipation term (excluding immersed boundaries) as
VISC = ∂ⱼ τ₃ⱼ,where τ₃ⱼ is the viscous stress tensor for the z-momentum equation.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> VISC = WMomentumEquation.ViscousDissipation(model)
KernelFunctionOperation at (Center, Center, Face)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: ∂ⱼ_τ₃ⱼ (generic function with 10 methods)
└── arguments: ("Nothing", "Nothing", "Clock", "NamedTuple", "Nothing")Oceanostics.WMomentumEquation.ImmersedViscousDissipation — Type
ImmersedViscousDissipation(
model,
velocities,
w_immersed_bc,
closure,
diffusivities,
clock,
model_fields;
location
)
Calculate the viscous dissipation term due to immersed boundaries as
VISC = ∂ⱼ τ₃ⱼ,where τ₃ⱼ is the immersed boundary viscous stress tensor for the z-momentum equation.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> VISC = WMomentumEquation.ImmersedViscousDissipation(model)
KernelFunctionOperation at (Center, Center, Face)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: immersed_∂ⱼ_τ₃ⱼ (generic function with 2 methods)
└── arguments: ("NamedTuple", "Nothing", "Nothing", "Nothing", "Clock", "NamedTuple")Oceanostics.WMomentumEquation.TotalViscousDissipation — Type
TotalViscousDissipation(
model,
velocities,
w_immersed_bc,
closure,
diffusivities,
clock,
model_fields,
buoyancy;
location
)
Calculate the total viscous dissipation term as
VISC = ∂ⱼ τ₃ⱼ + ∂ⱼ τ₃ⱼ_immersed,where τ₃ⱼ is the interior viscous stress tensor and τ₃ⱼ_immersed is the immersed boundary viscous stress tensor for the z-momentum equation.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> VISC = WMomentumEquation.TotalViscousDissipation(model)
KernelFunctionOperation at (Center, Center, Face)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: total_∂ⱼ_τ₃ⱼ (generic function with 1 method)
└── arguments: ("NamedTuple", "Nothing", "Nothing", "Nothing", "Clock", "NamedTuple", "Nothing")Oceanostics.WMomentumEquation.StokesShear — Type
StokesShear(model, stokes_drift, velocities, time; location)
Calculate the Stokes shear term as
STOKES_SHEAR = ((∇ × uˢ) × u)_zwhere uˢ is the Stokes drift velocity and u is the velocity vector.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> STOKES = WMomentumEquation.StokesShear(model)
KernelFunctionOperation at (Center, Center, Face)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: z_curl_Uˢ_cross_U (generic function with 10 methods)
└── arguments: ("Nothing", "NamedTuple", "Float64")Oceanostics.WMomentumEquation.StokesTendency — Type
StokesTendency(model, stokes_drift, time; location)
Calculate the Stokes tendency term as
STOKES_TEND = ∂wˢ/∂twhere wˢ is the z-component of the Stokes drift velocity.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> STOKES = WMomentumEquation.StokesTendency(model)
KernelFunctionOperation at (Center, Center, Face)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: ∂t_wˢ (generic function with 10 methods)
└── arguments: ("Nothing", "Float64")Oceanostics.WMomentumEquation.Forcing — Type
Forcing(
model,
forcing_func,
clock,
model_fields,
;
location
)
Calculate the forcing term Fʷ on the z-momentum equation for model.
Forcing is a type alias for the generic KernelFunctionOperation (no narrowing on the kernel function), so a constructor on Forcing(model) would clash across the U/V/W modules. To disambiguate, the W-momentum convenience constructor takes an explicit Val(:w) tag.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> FORC = WMomentumEquation.Forcing(model, Val(:w))
KernelFunctionOperation at (Center, Center, Face)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: Returns (generic function with 1 method)
└── arguments: ("Clock", "NamedTuple")Oceanostics.WMomentumEquation.Tendency — Type
Tendency(
model,
advection_scheme,
coriolis,
stokes_drift,
closure,
w_immersed_bc,
buoyancy,
background_fields,
velocities,
tracers,
auxiliary_fields,
diffusivities,
hydrostatic_pressure,
clock,
forcing_func;
location
)
Calculate the total tendency of the w-momentum equation as computed by Oceananigans.
For NonhydrostaticModel, this includes:
- Advection: -∇⋅(𝐯w)
- Background advection terms
- Buoyancy: ĝ_z b
- Coriolis: -(f × u)_z
- Pressure gradient: -∂p/∂z
- Viscous dissipation: -∇⋅τ₃
- Immersed viscous dissipation
- Stokes shear: ((∇ × uˢ) × u)_z
- Stokes tendency: ∂wˢ/∂t
- Forcing: Fʷ
HydrostaticFreeSurfaceModel does not have a prognostic w-momentum equation (w is diagnosed from continuity), so Tendency is only defined for NonhydrostaticModel.
julia> using Oceananigans, Oceanostics
julia> grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1));
julia> model = NonhydrostaticModel(grid);
julia> TEND = WMomentumEquation.Tendency(model)
KernelFunctionOperation at (Center, Center, Face)
├── grid: 4×4×4 RectilinearGrid{Float64, Periodic, Periodic, Bounded} on CPU with 3×3×3 halo
├── kernel_function: w_velocity_tendency (generic function with 1 method)
└── arguments: ("Centered", "Nothing", "Nothing", "Nothing", "Nothing", "Nothing", "Oceananigans.Models.NonhydrostaticModels.BackgroundFields", "NamedTuple", "NamedTuple", "NamedTuple", "Nothing", "Nothing", "Clock", "Returns")