SingleLayerQG
Basic Equations
This module solves the barotropic or equivalent barotropic quasi-geostrophic vorticity equation on a beta plane of variable fluid depth $H - h(x, y)$. The flow is obtained through a streamfunction $\psi$ as $(u, v) = (-\partial_y \psi, \partial_x \psi)$. All flow fields can be obtained from the quasi-geostrophic potential vorticity (QGPV). Here the QGPV is
\[ \underbrace{f_0 + \beta y}_{\text{planetary PV}} + \underbrace{\partial_x v - \partial_y u}_{\text{relative vorticity}} \underbrace{ - \frac{1}{\ell^2} \psi}_{\text{vortex stretching}} + \underbrace{\frac{f_0 h}{H}}_{\text{topographic PV}} ,\]
where $\ell$ is the Rossby radius of deformation. Purely barotropic dynamics corresponds to infinite Rossby radius of deformation ($\ell = \infty$), while a flow with a finite Rossby radius follows is said to obey equivalent-barotropic dynamics. We denote the sum of the relative vorticity and the vortex stretching contributions to the QGPV with $q \equiv \nabla^2 \psi - \psi / \ell^2$. Also, we denote the topographic PV with $\eta \equiv f_0 h / H$.
The dynamical variable is $q$. Thus, the equation solved by the module is:
\[\partial_t q + \mathsf{J}(\psi, q + \eta) + \beta \partial_x \psi = \underbrace{-\left[\mu + \nu(-1)^{n_\nu} \nabla^{2n_\nu} \right] q}_{\textrm{dissipation}} + F .\]
where $\mathsf{J}(a, b) = (\partial_x a)(\partial_y b)-(\partial_y a)(\partial_x b)$ is the two-dimensional Jacobian. On the right hand side, $F(x, y, t)$ is forcing, $\mu$ is linear drag, and $\nu$ is hyperviscosity of order $n_\nu$. Plain old viscosity corresponds to $n_\nu = 1$.
Implementation
The equation is time-stepped forward in Fourier space:
\[\partial_t \widehat{q} = - \widehat{\mathsf{J}(\psi, q + \eta)} + \beta \frac{i k_x}{|𝐤|^2 + 1/\ell^2} \widehat{q} - \left(\mu + \nu |𝐤|^{2n_\nu} \right) \widehat{q} + \widehat{F} .\]
The state variable sol
is the Fourier transform of the sum of relative vorticity and vortex stretching (when the latter is applicable), qh
.
The Jacobian is computed in the conservative form: $\mathsf{J}(f, g) = \partial_y [ (\partial_x f) g] - \partial_x[ (\partial_y f) g]$.
The linear operator is constructed in Equation
GeophysicalFlows.SingleLayerQG.Equation
— FunctionEquation(params::BarotropicQGParams, grid)
Return the equation for a barotropic QG problem with params
and grid
. Linear operator $L$ includes bottom drag $μ$, (hyper)-viscosity of order $n_ν$ with coefficient $ν$ and the $β$ term:
\[L = - μ - ν |𝐤|^{2 n_ν} + i β k_x / |𝐤|² .\]
The nonlinear term is computed via calcN!
function.
Equation(params::EquivalentBarotropicQGParams, grid)
Return the equation for an equivalent-barotropic QG problem with params
and grid
. Linear operator $L$ includes bottom drag $μ$, (hyper)-viscosity of order $n_ν$ with coefficient $ν$ and the $β$ term:
\[L = -μ - ν |𝐤|^{2 n_ν} + i β k_x / (|𝐤|² + 1/ℓ²) .\]
The nonlinear term is computed via calcN!
function.
The nonlinear terms are computed via
GeophysicalFlows.SingleLayerQG.calcN!
— FunctioncalcN!(N, sol, t, clock, vars, params, grid)
Calculate the nonlinear term, that is the advection term and the forcing,
\[N = - \widehat{𝖩(ψ, q + η)} + F̂ .\]
which in turn calls calcN_advection!
and addforcing!
.
Parameters and Variables
All required parameters are included inside Params
and all module variables are included inside Vars
.
For the decaying case (no forcing, $F = 0$), variables are constructed with Vars
. For the forced case ($F \ne 0$) variables are constructed with either ForcedVars
or StochasticForcedVars
.
Helper functions
Some helper functions included in the module are:
GeophysicalFlows.SingleLayerQG.updatevars!
— Functionupdatevars!(sol, vars, params, grid)
Update the variables in vars
with the solution in sol
.
GeophysicalFlows.SingleLayerQG.set_q!
— Functionset_q!(prob, q)
Set the solution of problem, prob.sol
as the transform of $q$ and update variables prob.vars
.
Diagnostics
The kinetic energy of the fluid is computed via:
GeophysicalFlows.SingleLayerQG.kinetic_energy
— Functionkinetic_energy(prob)
Return the problem's (prob
) domain-averaged kinetic energy of the fluid. Since $u² + v² = |{\bf ∇} ψ|²$, the domain-averaged kinetic energy is
\[\int \frac1{2} |{\bf ∇} ψ|² \frac{𝖽x 𝖽y}{L_x L_y} = \sum_{𝐤} \frac1{2} |𝐤|² |ψ̂|² .\]
while the potential energy, for an equivalent barotropic fluid, is computed via:
GeophysicalFlows.SingleLayerQG.potential_energy
— Functionpotential_energy(prob)
Return the problem's (prob
) domain-averaged potential energy of the fluid,
\[\int \frac1{2} \frac{ψ²}{ℓ²} \frac{𝖽x 𝖽y}{L_x L_y} = \sum_{𝐤} \frac1{2} \frac{|ψ̂|²}{ℓ²} .\]
The total energy is:
GeophysicalFlows.SingleLayerQG.energy
— Functionenergy(prob)
Return the problem's (prob
) domain-averaged total energy of the fluid, that is, the kinetic energy for a pure barotropic flow or the sum of kinetic and potential energies for an equivalent barotropic flow.
Other diagnostic include: energy_dissipation
, energy_drag
, energy_work
, enstrophy_dissipation
, and enstrophy_drag
, enstrophy_work
.
Examples
examples/singlelayerqg_betadecay.jl
: Simulate decaying quasi-geostrophic flow on a beta plane demonstrating zonation.examples/singlelayerqg_betaforced.jl
: Simulate forced-dissipative quasi-geostrophic flow on a beta plane demonstrating zonation. The forcing is temporally delta-correlated with isotropic spatial structure with power in a narrow annulus in wavenumber space with total wavenumber $k_f$.examples/singlelayerqg_decay_topography.jl
: Simulate two dimensional turbulence (barotropic quasi-geostrophic flow with $\beta=0$) above topography.examples/singlelayerqg_decaying_barotropic_equivalentbarotropic.jl
: Simulate two dimensional turbulence ($\beta=0$) with both infinite and finite Rossby radius of deformation and compares the evolution of the two.