BarotropicQGQL
Basic Equations
This module solves the quasi-linear quasi-geostrophic barotropic vorticity equation on a beta plane of variable fluid depth $H - h(x, y)$. Quasi-linear refers to the dynamics that neglects the eddy–eddy interactions in the eddy evolution equation after an eddy–mean flow decomposition, e.g.,
\[\phi(x, y, t) = \overline{\phi}(y, t) + \phi'(x, y, t) ,\]
where overline above denotes a zonal mean, $\overline{\phi}(y, t) = \int \phi(x, y, t) \, 𝖽x / L_x$, and prime denotes deviations from the zonal mean. This approximation is used in many process-model studies of zonation, e.g., Farrell and Ioannou (2003), Srinivasan and Young (2012), and Constantinou et al. (2014).
As in the SingleLayerQG module, 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{f_0 h}{H}}_{\text{topographic PV}} .\]
The dynamical variable is the component of the vorticity of the flow normal to the plane of motion, $\zeta \equiv \partial_x v - \partial_y u = \nabla^2 \psi$. Also, we denote the topographic PV with $\eta \equiv f_0 h / H$. After we apply the eddy-mean flow decomposition above, the QGPV dynamics are:
\[\begin{aligned} \partial_t \overline{\zeta} + \mathsf{J}(\overline{\psi}, \overline{\zeta} + \overline{\eta}) + \overline{\mathsf{J}(\psi', \zeta' + \eta')} & = \underbrace{- \left[\mu + \nu(-1)^{n_\nu} \nabla^{2n_\nu} \right] \overline{\zeta} }_{\textrm{dissipation}} , \\ \partial_t \zeta' + \mathsf{J}(\psi', \overline{\zeta} + \overline{\eta}) + \mathsf{J}(\overline{\psi}, \zeta' + \eta') + & \underbrace{\mathsf{J}(\psi', \zeta' + \eta') - \overline{\mathsf{J}(\psi', \zeta' + \eta')}}_{\textrm{EENL}} + \beta \partial_x \psi' = \\ & = \underbrace{-\left[\mu + \nu(-1)^{n_\nu} \nabla^{2n_\nu} \right] \zeta'}_{\textrm{dissipation}} + F . \end{aligned}\]
where $\mathsf{J}(a, b) = (\partial_x a)(\partial_y b) - (\partial_y a)(\partial_x b)$. On the right hand side, $F(x, y, t)$ is forcing (which is assumed to have zero zonal mean, $\overline{F} = 0$), $\mu$ is linear drag, and $\nu$ is hyperviscosity. Plain old viscosity corresponds to $n_{\nu} = 1$.
Quasi-linear dynamics neglects the term eddy-eddy nonlinearity (EENL) term above.
Implementation
The equation is time-stepped forward in Fourier space:
\[\partial_t \widehat{\zeta} = - \widehat{\mathsf{J}(\psi, \zeta + \eta)}^{\textrm{QL}} + \beta \frac{i k_x}{|𝐤|^2} \widehat{\zeta} - \left ( \mu + \nu |𝐤|^{2n_\nu} \right ) \widehat{\zeta} + \widehat{F} .\]
The state variable sol
is the Fourier transform of vorticity, ζh
.
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 superscript QL on the Jacobian term above denotes that triad interactions that correspond to the EENL term are removed.
The linear operator is constructed in Equation
GeophysicalFlows.BarotropicQGQL.Equation
— FunctionEquation(params, grid)
Return the equation for two-dimensional barotropic QG QL problem with parameters params
and on grid
. Linear operator $L$ includes bottom drag $μ$, (hyper)-viscosity of order $n_ν$ with coefficient $ν$ and the $β$ term:
\[L = - μ - ν |𝐤|^{2 n_ν} + i β k_x / |𝐤|² .\]
Nonlinear term is computed via calcN!
function.
and the nonlinear terms are computed via
GeophysicalFlows.BarotropicQGQL.calcN!
— FunctioncalcN!(N, sol, t, clock, vars, params, grid)
Calculate the nonlinear term, that is the advection term and the forcing,
\[N = - \widehat{𝖩(ψ, ζ + η)}^{\mathrm{QL}} + 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
GeophysicalFlows.BarotropicQGQL.updatevars!
— Functionupdatevars!(sol, vars, params, grid)
updatevars!(prob)
Update the vars
of a problem prob
that has grid
and params
with the solution in sol
.
GeophysicalFlows.BarotropicQGQL.set_zeta!
— Functionset_zeta!(prob, zeta)
set_zeta!(sol, vars, grid, zeta)
Set the solution sol
as the transform of zeta
and update variables vars
on the grid
.
Diagnostics
The kinetic energy of the fluid is obtained via:
GeophysicalFlows.BarotropicQGQL.energy
— Functionenergy(sol, grid)
energy(prob)
Return the domain-averaged kinetic energy of sol
.
while the enstrophy via:
GeophysicalFlows.BarotropicQGQL.enstrophy
— Functionenstrophy(sol, grid, vars)
enstrophy(prob)
Return the domain-averaged enstrophy of sol
.
Other diagnostic include: dissipation
, drag
, and work
.
Examples
examples/barotropicqgql_betaforced.jl
: Simulate forced-dissipative quasi-linear quasi-geostrophic flow on a beta plane demonstrating zonation. The forcing is temporally delta-correlated and its spatial structure is isotropic with power in a narrow annulus of total radius $k_f$ in wavenumber space. This example demonstrates that the anisotropic inverse energy cascade is not required for zonation.