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.,
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-Ioannou-2003, Srinivasan-Young-2012, and Constantinou-etal-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
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:
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.
The equation is time-stepped forward in Fourier space:
The state variable sol
is the Fourier transform of vorticity, [ζh
](@ref GeophysicalFlows.BarotropicQGQL.Vars).
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
and the nonlinear terms are computed via
GeophysicalFlows.BarotropicQGQL.calcN!
which in turn calls [calcN_advection!
](@ref GeophysicalFlows.BarotropicQGQL.calcN_advection!)
and [addforcing!
](@ref GeophysicalFlows.BarotropicQGQL.addforcing!).
All required parameters are included inside [Params
](@ref GeophysicalFlows.BarotropicQGQL.Params)
and all module variables are included inside [Vars
](@ref GeophysicalFlows.BarotropicQGQL.Vars).
For the decaying case (no forcing, F = 0
), variables are constructed with [Vars
](@ref GeophysicalFlows.BarotropicQGQL.Vars).
For the forced case (F \ne 0
) variables are constructed with either [ForcedVars
](@ref GeophysicalFlows.BarotropicQGQL.ForcedVars)
or [StochasticForcedVars
](@ref GeophysicalFlows.BarotropicQGQL.StochasticForcedVars).
GeophysicalFlows.BarotropicQGQL.updatevars!
GeophysicalFlows.BarotropicQGQL.set_zeta!
The kinetic energy of the fluid is obtained via:
GeophysicalFlows.BarotropicQGQL.energy
while the enstrophy via:
GeophysicalFlows.BarotropicQGQL.enstrophy
Other diagnostic include: [dissipation
](@ref GeophysicalFlows.BarotropicQGQL.dissipation),
[drag
](@ref GeophysicalFlows.BarotropicQGQL.drag), and [work
](@ref GeophysicalFlows.BarotropicQGQL.work).
- [
examples/barotropicqgql_betaforced.jl
](@ref barotropicqgql_betaforced_example): 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 radiusk_f
in wavenumber space. This example demonstrates that the anisotropic inverse energy cascade is not required for zonation.