Library
Documentation for the public user interface.
Oceanostics.jl
Oceanostics.TKEBudgetTerms
Oceanostics.TKEBudgetTerms.IsotropicKineticEnergyDissipationRate — MethodIsotropicKineticEnergyDissipationRate(
model;
U,
V,
W,
location
)
Calculate the Viscous Dissipation Rate, defined as
ε = 2 ν SᵢⱼSᵢⱼ,where Sᵢⱼ is the strain rate tensor, for a fluid with an isotropic turbulence closure (i.e., a turbulence closure where ν (eddy or not) is the same for all directions.
Oceanostics.TKEBudgetTerms.KineticEnergy — MethodKineticEnergy(model, u, v, w; location, kwargs...)
Calculate the kinetic energy of model manually specifying u, v and w.
Oceanostics.TKEBudgetTerms.KineticEnergy — MethodKineticEnergy(model; kwargs...)
Calculate the kinetic energy of model.
Oceanostics.TKEBudgetTerms.KineticEnergyDiffusiveTerm — MethodKineticEnergyDiffusiveTerm(model; location)
Return a KernelFunctionOperation that computes the diffusive term of the KE prognostic equation:
DIFF = uᵢ∂ⱼτᵢⱼwhere uᵢ are the velocity components and τᵢⱼ is the diffusive flux of i momentum in the j-th direction.
Oceanostics.TKEBudgetTerms.KineticEnergyDissipationRate — MethodKineticEnergyDissipationRate(model; U, V, W, location)
Calculate the pseudo viscous Dissipation Rate, defined as
ε = ν (∂uᵢ/∂xⱼ) (∂uᵢ/∂xⱼ)for a fluid with an isotropic turbulence closure (i.e., a turbulence closure where ν (eddy or not) is the same for all directions.
Oceanostics.TKEBudgetTerms.KineticEnergyForcingTerm — MethodKineticEnergyForcingTerm(model; location)
Return a KernelFunctionOperation that computes the forcing term of the KE prognostic equation:
FORC = uᵢFᵤᵢwhere uᵢ are the velocity components and Fᵤᵢ is the forcing term(s) in the uᵢ prognostic equation (i.e. the forcing for uᵢ).
Oceanostics.TKEBudgetTerms.KineticEnergyTendency — MethodKineticEnergyTendency(model; location)
Return a KernelFunctionOperation that computes the tendency of the KE except for the nonhydrostatic pressure:
Oceanostics.TKEBudgetTerms.TurbulentKineticEnergy — MethodTurbulentKineticEnergy(model, u, v, w; U, V, W, location)
Calculate the turbulent kinetic energy of model manually specifying u, v, w and optionally background velocities U, V and W.
Oceanostics.TKEBudgetTerms.TurbulentKineticEnergy — MethodTurbulentKineticEnergy(model; kwargs...)
Calculate the turbulent kinetic energy of model.
Oceanostics.TKEBudgetTerms.XPressureRedistribution — MethodXPressureRedistribution(model, u′, p′)
Calculate the pressure redistribution term in the x direction. Here u′ and p′ are the fluctuations around a mean.
Oceanostics.TKEBudgetTerms.XShearProductionRate — MethodXShearProductionRate(model, u, v, w, U, V, W; location)
Calculate the shear production rate in the model's x direction, considering velocities u, v, w and background (or average) velocities U, V and W.
Oceanostics.TKEBudgetTerms.XShearProductionRate — MethodXShearProductionRate(model; U, V, W, kwargs...)
Calculate the shear production rate in the model's x direction. At least one of the mean velocities U, V and W must be specified otherwise the output will be zero.
Oceanostics.TKEBudgetTerms.YPressureRedistribution — MethodYPressureRedistribution(model, v′, p′)
Calculate the pressure redistribution term in the y direction. Here v′ and p′ are the fluctuations around a mean.
Oceanostics.TKEBudgetTerms.YShearProductionRate — MethodYShearProductionRate(model, u, v, w, U, V, W; location)
Calculate the shear production rate in the model's y direction, considering velocities u, v, w and background (or average) velocities U, V and W.
Oceanostics.TKEBudgetTerms.YShearProductionRate — MethodYShearProductionRate(model; U, V, W, kwargs...)
Calculate the shear production rate in the model's y direction. At least one of the mean velocities U, V and W must be specified otherwise the output will be zero.
Oceanostics.TKEBudgetTerms.ZPressureRedistribution — MethodZPressureRedistribution(model, w′, p′)
Calculate the pressure redistribution term in the z direction. Here w′ and p′ are the fluctuations around a mean.
Oceanostics.TKEBudgetTerms.ZShearProductionRate — MethodZShearProductionRate(model, u, v, w, U, V, W; location)
Calculate the shear production rate in the model's z direction, considering velocities u, v, w and background (or average) velocities U, V and W.
Oceanostics.TKEBudgetTerms.ZShearProductionRate — MethodZShearProductionRate(model; U, V, W, kwargs...)
Calculate the shear production rate in the model's z direction. At least one of the mean velocities U, V and W must be specified otherwise the output will be zero.
Oceanostics.TracerVarianceBudgetTerms
Oceanostics.TracerVarianceBudgetTerms.TracerVarianceDiffusiveTerm — MethodTracerVarianceDiffusiveTerm(model, tracer_name; location)
Return a KernelFunctionOperation that computes the diffusive term of the tracer variance prognostic equation using Oceananigans' diffusive tracer flux divergence kernel:
DIFF = 2 c ∂ⱼFⱼwhere c is the tracer, and Fⱼ is the tracer's diffusive flux in the j-th direction.
grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1))
model = NonhydrostaticModel(grid=grid, tracers=:b, closure=SmagorinskyLilly())
DIFF = TracerVarianceDiffusiveTerm(model, :b)Oceanostics.TracerVarianceBudgetTerms.TracerVarianceDissipationRate — MethodTracerVarianceDissipationRate(
model,
tracer_name;
tracer,
location
)
Return a KernelFunctionOperation that computes the isotropic variance dissipation rate for tracer_name in model.tracers. The isotropic variance dissipation rate is defined as
χ = 2 ∂ⱼc ⋅ Fⱼwhere Fⱼ is the diffusive flux of c in the j-th direction and ∂ⱼ is the gradient operator. χ is implemented in its conservative formulation based on the equation above.
Note that often χ is written as χ = 2κ (∇c ⋅ ∇c), which is the special case for Fickian diffusion (κ is the tracer diffusivity).
Here tracer_name is needed even when passing tracer in order to get the appropriate tracer_index. When passing tracer, this function should be used as
grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1))
model = NonhydrostaticModel(grid=grid, tracers=:b, closure=SmagorinskyLilly())
b̄ = Field(Average(model.tracers.b, dims=(1,2)))
b′ = model.tracers.b - b̄
χb = TracerVarianceDissipationRate(model, :b, tracer=b′)Oceanostics.TracerVarianceBudgetTerms.TracerVarianceTendency — MethodTracerVarianceTendency(model, tracer_name; location)
Return a KernelFunctionOperation that computes the tracer variance tendency:
TEND = 2 c ∂ₜcwhere c is the tracer and ∂ₜc is the tracer tendency (computed using Oceananigans' tracer tendency kernel).
grid = RectilinearGrid(size=(4, 4, 4), extent=(1, 1, 1))
model = NonhydrostaticModel(grid=grid, tracers=:b, closure=SmagorinskyLilly())
DIFF = TracerVarianceTendency(model, :b)Oceanostics.FlowDiagnostics
Oceanostics.FlowDiagnostics.DirectionalErtelPotentialVorticity — MethodDirectionalErtelPotentialVorticity(
model,
direction;
location
)
Calculate the contribution from a given direction to the Ertel Potential Vorticity basde on a model and a direction. The Ertel Potential Vorticity is defined as
EPV = ωₜₒₜ ⋅ ∇bwhere ωₜₒₜ is the total (relative + planetary) vorticity vector, b is the buoyancy and ∇ is the gradient operator.
Oceanostics.FlowDiagnostics.ErtelPotentialVorticity — MethodErtelPotentialVorticity(model; location, add_background)
Calculate the Ertel Potential Vorticty for model, where the characteristics of the Coriolis rotation are taken from model.coriolis. The Ertel Potential Vorticity is defined as
EPV = ωₜₒₜ ⋅ ∇bwhere ωₜₒₜ is the total (relative + planetary) vorticity vector, b is the buoyancy and ∇ is the gradient operator.
Oceanostics.FlowDiagnostics.QVelocityGradientTensorInvariant — MethodQVelocityGradientTensorInvariant(model; location)
Calculate the value of the Q velocity gradient tensor invariant. This is usually just called Q and it is generally used for identifying and visualizing vortices in fluid flow.
The definition and nomenclature comes from the equation for the eigenvalues λ of the velocity gradient tensor ∂ⱼuᵢ:
λ³ + P λ² + Q λ + T = 0from where Q is defined as
Q = ½ ( ΩᵢⱼΩᵢⱼ - SᵢⱼSᵢⱼ)and where Sᵢⱼ= ½(∂ⱼuᵢ + ∂ᵢuⱼ) and Ωᵢⱼ= ½(∂ⱼuᵢ - ∂ᵢuⱼ). More info about it can be found in doi:10.1063/1.5124245.
Oceanostics.FlowDiagnostics.RichardsonNumber — MethodRichardsonNumber(model; location, add_background)
Calculate the Richardson Number as
Ri = (∂b/∂z) / (|∂u⃗ₕ/∂z|²)where z is the true vertical direction (ie anti-parallel to gravity).
Oceanostics.FlowDiagnostics.RossbyNumber — MethodRossbyNumber(
model;
location,
add_background,
dWdy_bg,
dVdz_bg,
dUdz_bg,
dWdx_bg,
dUdy_bg,
dVdx_bg
)
Calculate the Rossby number using the vorticity in the rotation axis direction according to model.coriolis. Rossby number is defined as
Ro = ωᶻ / fwhere ωᶻ is the vorticity in the Coriolis axis of rotation and f is the Coriolis rotation frequency.
Oceanostics.FlowDiagnostics.StrainRateTensorModulus — MethodStrainRateTensorModulus(model; location)
Calculate the modulus (absolute value) of the strain rate tensor S, which is defined as the symmetric part of the velocity gradient tensor:
Sᵢⱼ = ½(∂ⱼuᵢ + ∂ᵢuⱼ)Its modulus is then defined (using Einstein summation notation) as
|| Sᵢⱼ || = √( Sᵢⱼ Sᵢⱼ)Oceanostics.FlowDiagnostics.ThermalWindPotentialVorticity — MethodThermalWindPotentialVorticity(model; f, location)
Calculate the Potential Vorticty assuming thermal wind balance for model, where the characteristics of the Coriolis rotation are taken from model.coriolis. The Potential Vorticity in this case is defined as
TWPV = (f + ωᶻ) ∂b/∂z - f ((∂U/∂z)² + (∂V/∂z)²)where f is the Coriolis frequency, ωᶻ is the relative vorticity in the z direction, b is the buoyancy, and ∂U/∂z and ∂V/∂z comprise the thermal wind shear.
Oceanostics.FlowDiagnostics.VorticityTensorModulus — MethodVorticityTensorModulus(model; location)
Calculate the modulus (absolute value) of the vorticity tensor Ω, which is defined as the antisymmetric part of the velocity gradient tensor:
Ωᵢⱼ = ½(∂ⱼuᵢ - ∂ᵢuⱼ)Its modulus is then defined (using Einstein summation notation) as
|| Ωᵢⱼ || = √( Ωᵢⱼ Ωᵢⱼ)