Functions
This page lists all Modules, Functions and Structures available in this package.
QGDipoles.jl Module
QGDipoles
— ModulePackage for creating steady modon solutions to the layered quasi-geostrophic equations and the surface quasi geostrophic equations.
See examples/
for example scripts and here for documentation.
Vortex Structures
LQG
QGDipoles.LQGParams
— TypeLQGParams
Stores the parameters for an LQG dipolar vortex solution
Fields:
U
: Vortex speedℓ
: Vortex radiusR
: Rossby radiusβ
: background PV gradientActiveLayers
: 1 => layer contains vortex regionH
: thickness of each layerx₀
: Vortex positionα
: Direction of vortex propagationM
: number of coefficients in Zernike expansiontol
: maximum error in solution evaluationK₀
: initial guess for eigenvaluea₀
: initial guess for coefficientsUseAnalytic
: use analytic solution (1-layer only)CalcVelocity
: flag to determine if velocity is calculatedCalcVorticity
: flag to determine if vorticity is calculatedCalcEnergy
: flag to determine if energy is calculatedCalcEnstrophy
: flag to determine if enstrophy is calculated
QGDipoles.LQGVortex
— TypeLQGVortex
Stores fields and diagnostics for an LQG dipolar vortex solution
Fields:
params
: Vortex paramsψ
: streamfunctionq
: potential vorticity anomalyK
: eigenvaluea
: coefficient matrixu
: x velocityv
: y velocityζ
: vertical vorticityKE
: kinetic energyPE
: potential energyEN
: enstrophy
QGDipoles.DefLQGParams
— FunctionDefLQGParams(; U=1, ℓ=1, R=Inf, β=0, ActiveLayers=1, H=1, x₀=[0,0], α=0, M=8, tol=1e-6, K₀=nothing, a₀=nothing, UseAnalytic=false, CalcVelocity=false, CalcVorticity=false, CalcEnergy=false, CalcEnstrophy=false)
Defines an LQGParams
structure using the given inputs
Keyword arguments:
U
: Vortex speedℓ
: Vortex radiusR
: Rossby radiusβ
: background PV gradientActiveLayers
: 1 => layer contains vortex regionH
: thickness of each layerx₀
: Vortex positionα
: Direction of vortex propagationM
: number of coefficients in Zernike expansiontol
: maximum error in solution evaluationK₀
: initial guess for eigenvaluea₀
: initial guess for coefficientsUseAnalytic
: use analytic solution (1-layer only)CalcVelocity
: flag to determine if velocity is calculatedCalcVorticity
: flag to determine if vorticity is calculatedCalcEnergy
: flag to determine if energy is calculatedCalcEnstrophy
: flag to determine if enstrophy is calculated
QGDipoles.DefLQGVortex
— FunctionDefLQGVortex(grid; U=1, ℓ=1, R=Inf, β=0, ActiveLayers=1, H=1, x₀=[0,0], α=0, M=8, tol=1e-6, K₀=nothing, a₀=nothing, UseAnalytic=false, CalcVelocity=false, CalcVorticity=false, CalcEnergy=false, CalcEnstrophy=false)
Defines an LQGVortex
solution structure using the given inputs
Arguments:
grid
: grid structure
Keyword arguments:
U
: Vortex speedℓ
: Vortex radiusR
: Rossby radiusβ
: background PV gradientActiveLayers
: 1 => layer contains vortex regionH
: thickness of each layerx₀
: Vortex positionα
: Direction of vortex propagationM
: number of coefficients in Zernike expansiontol
: maximum error in solution evaluationK₀
: initial guess for eigenvaluea₀
: initial guess for coefficientsUseAnalytic
: use analytic solution (1-layer only)CalcVelocity
: flag to determine if velocity is calculatedCalcVorticity
: flag to determine if vorticity is calculatedCalcEnergy
: flag to determine if energy is calculatedCalcEnstrophy
: flag to determine if enstrophy is calculated
DefLQGVortex(grid, params)
Defines an LQGVortex
solution structure using the given inputs
Arguments:
grid
: grid structureparams
: vortex parameters,LQGParams
structure
SQG
QGDipoles.SQGParams
— TypeSQGParams
Stores the parameters for an SQG dipolar vortex solution
Fields:
U
: Vortex speedℓ
: Vortex radiusR
: Rossby radiusβ
: background PV gradientx₀
: Vortex positionα
: Direction of vortex propagationM
: number of coefficients in Zernike expansiontol
: maximum error in solution evaluationK₀
: initial guess for eigenvaluea₀
: initial guess for coefficientsCalcVelocity
: flag to determine if velocity is calculatedCalcVorticity
: flag to determine if vorticity is calculatedCalcEnergy
: flag to determine if energy is calculated
QGDipoles.SQGVortex
— TypeSQGVortex
Stores fields and diagnostics for an SQG dipolar vortex solution
Fields:
params
: Vortex paramsψ
: surface streamfunctionb
: surface buoyancyK
: eigenvaluea
: coefficient matrixu
: x velocityv
: y velocityζ
: vertical vorticityE
: domain integrated energySPE
: surface potential energy
QGDipoles.DefSQGParams
— FunctionDefSQGParams(; U=1, ℓ=1, R=[Inf,Inf], β=0, x₀=[0,0], α=0, M=12, tol=1e-6, K₀=nothing, a₀=nothing, CalcVelocity=false, CalcVorticity=false, CalcEnergy=false)
Defines an SQGParams
structure using the given inputs
Keyword arguments:
U
: Vortex speedℓ
: Vortex radiusR
: Rossby radiusβ
: background PV gradientx₀
: Vortex positionα
: Direction of vortex propagationM
: number of coefficients in Zernike expansiontol
: maximum error in solution evaluationK₀
: initial guess for eigenvaluea₀
: initial guess for coefficientsCalcVelocity
: flag to determine if velocity is calculatedCalcVorticity
: flag to determine if vorticity is calculatedCalcEnergy
: flag to determine if energy is calculated
QGDipoles.DefSQGVortex
— FunctionDefSQGVortex(grid; U=1, ℓ=1, R=[Inf,Inf], β=0, x₀=[0,0], α=0, M=12, tol=1e-6, K₀=nothing, a₀=nothing, CalcVelocity=false, CalcVorticity=false, CalcEnergy=false)
Defines an SQGVortex
solution structure using the given inputs
Arguments:
grid
: grid structure
Keyword arguments:
U
: Vortex speedℓ
: Vortex radiusR
: Rossby radiusβ
: background PV gradientx₀
: Vortex positionα
: Direction of vortex propagationM
: number of coefficients in Zernike expansiontol
: maximum error in solution evaluationK₀
: initial guess for eigenvaluea₀
: initial guess for coefficientsCalcVelocity
: flag to determine if velocity is calculatedCalcVorticity
: flag to determine if vorticity is calculatedCalcEnergy
: flag to determine if energy is calculated
DefSQGVortex(grid, params)
Defines an SQGVortex
solution structure using the given inputs
Arguments:
grid
: grid structureparams
: vortex parameters,LQGParams
structure
Shared
QGDipoles.UpdateParams
— FunctionUpdateParams(params::LQGParams; kwargs...)
Creates an LQGParams
structure by replacing parameters in params
with the given keywords
Arguments:
params
:LQGParams
parameter structure
Keyword arguments:
kwargs...
: keyword arguments forDefLQGParams
UpdateParams(params::SQGParams; kwargs...)
Creates an SQGParams
structure by replacing parameters in params
with the given keywords
Arguments:
params
:SQGParams
parameter structure
Keyword arguments:
kwargs...
: keyword arguments forDefSQGParams
QGDipoles.UpdateVortex
— FunctionUpdateVortex(grid, vortex::LQGVortex; kwargs...)
Creates an LQGVortex
structure by replacing parameters in vortex.params
with the given keywords
Arguments:
grid
: grid structurevortex
:LQGVortex
structure
Keyword arguments:
kwargs...
: keyword arguments forDefLQGParams
UpdateVortex(grid, vortex::SQGVortex; kwargs...)
Creates an SQGVortex
structure by replacing parameters in vortex.params
with the given keywords
Arguments:
grid
: grid structurevortex
:SQGVortex
structure
Keyword arguments:
kwargs...
: keyword arguments forDefSQGParams
High-Level Functions
LQG
QGDipoles.CreateModonLQG
— FunctionCreateModonLQG(grid; U=1, ℓ=1, R=Inf, β=0, ActiveLayers=1, x₀=[0, 0], α=0, M=8, tol=1e-6, K₀=nothing, a₀=nothing)
High level wrapper function for calculating $ψ$ and $q$ for the Layered QG model using given parameters
Arguments:
grid
: grid structure containingx
,y
, andKrsq
Keyword arguments:
- (
U
,ℓ
): vortex speed and radius, Numbers (default: (1
,1
)) - (
R
,β
): Rossby radii and (y) PV gradients in each layer, Numbers or Vectors, (default: (Inf
,0
)) ActiveLayers
: vector of 1s or 0s where 1 denotes an active layer, Number or Vector, (default:[1,..,1]
)x₀
: position of vortex center, vector (default:[0, 0]
)α
: initial angle of vortex, Number (default: 0)M
: number of coefficient to solve for, Integer (default:8
)tol
: error tolerance passed toQuadGK
andNLSolve
functions, Number (default:1e-6
)K₀
,a₀
: initial guesses for $K$ and $a$, Arrays or nothings (default:nothing
)
Note: provide values of K₀ and a₀ for active layers ONLY.
QGDipoles.CreateLCD
— FunctionCreateLCD(grid; U=1, ℓ=1, x₀=[0, 0], α=0)
High level wrapper function for calculating $ψ$ and $q$ for the Lamb-Chaplygin dipole using given parameters
Arguments:
grid
: grid structure containingx
,y
, andKrsq
Keyword arguments:
- (
U
,ℓ
): vortex speed and radius, Numbers (default: (1
,1
)) x₀
: position of vortex center, vector (default:[0, 0]
)α
: initial angle of vortex, Number (default:0
)
Note: This function uses the analytic solution for the LCD to calculate $ψ$ and $q$.
QGDipoles.CreateLRD
— FunctionCreateLRD(grid; U=1, ℓ=1, R=Inf, β=0, x₀=[0, 0], α=0)
High level wrapper function for calculating $ψ$ and $q$ for the Larichev-Reznik dipole using given parameters
Arguments:
grid
: grid structure containingx
,y
, andKrsq
Keyword arguments:
- (
U
,ℓ
): vortex speed and radius, Numbers (default: (1
,1
)) - (
R
,β
): Rossby radii and (y) PV gradient, Numbers, (default: (Inf
,0
)) x₀
: position of vortex center, vector (default:[0, 0]
)α
: initial angle of vortex, Number (default:0
)
Note: This function uses the analytic solution for the LRD to calculate $ψ$ and $q$.
SQG
QGDipoles.CreateModonSQG
— FunctionCreateModonSQG(grid; U=1, ℓ=1, R=[Inf, Inf], β=0, x₀=[0, 0], α=0, M=12, tol=1e-6, K₀=nothing, a₀=nothing)
High level wrapper function for calculating $ψ$ and $b$ for the SQG model using given parameters
Arguments:
grid
: grid structure containingx
,y
, andKrsq
Keyword arguments:
- (
U
,ℓ
): vortex speed and radius, Numbers (default: (1
,1
)) R
: vector of $[R, R']$, Vector (default:[Inf, Inf]
)β
: beta-plane (y) PV gradient, Number (default:0
)x₀
: position of vortex center, vector (default:[0, 0]
)α
: initial angle of vortex, Number (default:0
)M
: number of coefficient to solve for, Integer (default:12
)tol
: error tolerance passed toQuadGK
andNLSolve
functions, Number (default:1e-6
)K₀
,a₀
: initial guesses for $K$ and $a$, Arrays or nothings (default:nothing
)
Note: Here R is the baroclinic Rossby radius, R = NH/f, and R' = R₀²/R where R₀ is the barotropic Rossby radius, R₀ = √(gH)/f. For infinite depth, R' = g/(fN).
QGDipoles.Eval_ψ_SQG
— FunctionEval_ψ_SQG(grid, ψ; z=[0], U=1, R=[Inf, Inf], β=0)
Evaluates $ψ$ at specified depths, $z ∈ [-R, 0]$, for the SQG problem
Arguments:
grid
: grid structure containingx
,y
, andKrsq
ψ
: surface streamfunction, calculated usingCalc_ψb
orCreateModonSQG
Keyword arguments:
z
: vector of depths (default:[0]
)U
: vortex speed, Number (default:1
)R
: vector of $[R, R']$, Vector (default:[Inf, Inf]
)β
: beta-plane (y) PV gradient, Number (default:0
)
Note: Here $R$ is the baroclinic Rossby radius, $R = NH/f$, and $R' = R₀²/R$ where $R₀$ is the barotropic Rossby radius, $R₀ = √(gH)/f$. For infinite depth, $R' = g/(fN)$.
QGDipoles.Eval_q_SQG
— FunctionEval_q_SQG(grid, ψ; z=[0], U=1, R=[Inf, Inf], β=0)
Evaluates $q$ at specified depths, $z ∈ [-R, 0]$, for the SQG problem
Arguments:
grid
: grid structure containingx
,y
, andKrsq
ψ
: surface streamfunction, calculated usingCalc_ψb
orCreateModonSQG
Keyword arguments:
z
: vector of depths (default:[0]
)U
: vortex speed, Number (default:1
)R
: vector of $[R, R']$, Vector (default:[Inf, Inf]
)β
: beta-plane (y) PV gradient, Number (default:0
)
Note: Here $R$ is the baroclinic Rossby radius, $R = NH/f$, and $R' = R₀²/R$ where $R₀$ is the barotropic Rossby radius, $R₀ = √(gH)/f$. For infinite depth, $R' = g/(fN)$.
QGDipoles.Eval_b_SQG
— FunctionEval_b_SQG(grid, ψ; z=[0], U=1, R=[Inf, Inf], β=0)
Evaluates $b$ at specified depths, $z ∈ [-R, 0]$, for the SQG problem
Arguments:
grid
: grid structure containingx
,y
, andKrsq
ψ
: surface streamfunction, calculated usingCalc_ψb
orCreateModonSQG
Keyword arguments:
z
: vector of depths (default:[0]
)U
: vortex speed, Number (default:1
)R
: vector of $[R, R']$, Vector (default:[Inf, Inf]
)β
: beta-plane (y) PV gradient, Number (default:0
)
Note: Here $R$ is the baroclinic Rossby radius, $R = NH/f$, and $R' = R₀²/R$ where $R₀$ is the barotropic Rossby radius, $R₀ = √(gH)/f$. For infinite depth, $R' = g/(fN)$.
QGDipoles.Eval_w_SQG
— FunctionEval_w_SQG(grid, ψ; z=[0], U=1, R=[Inf, Inf], β=0)
Evaluates N²w at specified depths, $z ∈ [-R, 0]$, for the SQG problem using $N²w = -J[ψ + Uy, b]$
Arguments:
grid
: grid structure containingx
,y
, andKrsq
ψ
: surface streamfunction, calculated usingCalc_ψb
orCreateModonSQG
Keyword arguments:
z
: vector of depths (default:[0]
)U
: vortex speed, Number (default:1
)R
: vector of $[R, R']$, Vector (default:[Inf, Inf]
)β
: beta-plane (y) PV gradient, Number (default:0
)
Note: Here $R$ is the baroclinic Rossby radius, $R = NH/f$, and $R' = R₀²/R$ where $R₀$ is the barotropic Rossby radius, $R₀ = √(gH)/f$. For infinite depth, $R' = g/(fN)$.
Note: this function is not accurate at the surface as $∇b$ is discontinuous there. Instead use $w = -U∂η/∂x$ where $η = fψ/g$ is the surface elevation, or $w = 0$ if $R' = ∞$.
Utilities
QGDipoles.CreateGrid
— FunctionCreateGrid(Nx, Ny, Lx, Ly; cuda=false)
Define the numerical grid as a GridStruct
Arguments:
Nx
,Ny
: number of gridpoints in x and y directions, IntegersLx
,Ly
: x and y domains, either vectors of endpoints or lengths, Vectors or Numbers
Keyword arguments:
cuda
:true
; use CUDA CuArray for fields (default:false
)
CreateGrid(; Nx=512, Ny=512, Lx=[-5,5], Ly=[-5,5], cuda=false)
Define the numerical grid as a GridStruct
using a keyword-based method
Keyword arguments:
Nx
,Ny
: number of gridpoints in x and y directions, IntegersLx
,Ly
: x and y domains, either vectors of endpoints or lengths, Vectors or Numberscuda
:true
; useCUDA
CuArray for fields (default:false
)
QGDipoles.Calc_uv
— FunctionCalc_uv(grid, ψ)
Calculate the velocity fields from $ψ$ using $(u, v) = (-∂ψ/∂y, ∂ψ/∂x)$
Arguments:
grid
: grid structure containingkr
andl
ψ
: streamfunction, Array
QGDipoles.Calc_∇
— FunctionCalc_∇(grid, f)
Calculate the gradient $∇f$ for a given field $f$
Arguments:
grid
: grid structure containingkr
andl
f
: function, Array
QGDipoles.Calc_ζ
— FunctionCalc_ζ(grid, ψ)
Calculate the vertical vorticity using $ζ = ∂v/∂x - ∂u/∂y = ∇²ψ$
Arguments:
grid
: grid structure containingKrsq
ψ
: streamfunction, Array
Diagnostics
LQG
QGDipoles.EnergyLQG
— FunctionEnergyLQG(grid, ψ; R=Inf, H=[1])
Calculates the kinetic and potential energy for the LQG system
Arguments:
grid
: grid structure containingKrsq
ψ
: streamfunction in each layer, Array or CuArray
Keyword arguments:
R
: Rossby radius in each layer, Number or Vector (default:Inf
)H
: Thickness of each layer, Number or Vector (default:[1]
)
QGDipoles.EnstrophyLQG
— FunctionEnstrophyLQG(grid, q; H=[1])
Calculates the enstrophy for the LQG system
Arguments:
grid
: grid structure containingKrsq
q
: potential vorticity anomaly in each layer, Array or CuArray
Keyword arguments:
H
: Thickness of each layer, Number or Vector (default:[1]
)
SQG
QGDipoles.EnergySQG
— FunctionEnergySQG(grid, ψ, b; R′)
Calculates the energies for the SQG system; the total domain integrated energy and the surface potential energy
Arguments:
grid
: grid structure containingKrsq
ψ
: surface streamfunction, Array or CuArrayb
: surface buoyancy, , Array or CuArray
Keyword arguments:
R′
: reduced barotropic Rossby radius, Number (default:Inf
)
Note: the surface potential energy is sometimes referred to as the generalised enstrophy or the buoyancy variance.
Low-Level Functions
LQG
QGDipoles.BuildLinSysLQG
— FunctionBuildLinSysLQG(M, λ, μ; tol=1e-6)
Builds the terms in the inhomogeneous eigenvalue problem; $A$, $B$, $c$, $d$ for the LQG problem
Arguments:
M
: number of coefficient to solve for, Integerλ
: ratio of vortex radius to Rossby radius in each layer, Number or Vectorμ
: nondimensional (y) vorticity gradient in each layer, Number or Vector
Keyword arguments:
tol
: error tolerance forQuadGK
viaJJ_int
, Number (default:1e-6
)
QGDipoles.ApplyPassiveLayers
— FunctionApplyPassiveLayers(A, B, c, d, ActiveLayers)
Removes rows and columns corresponding to passive layers from the system
Arguments:
A
,B
,c
,d
: inhomogeneous eigenvalue problem terms, ArraysActiveLayers
: vector of 1s or 0s where 1 denotes an active layer, Number or Vector
QGDipoles.IncludePassiveLayers
— FunctionIncludePassiveLayers(K, a, ActiveLayers)
Includes columns corresponding to passive layers in the eigenvalue and coefficient arrays
Arguments:
K
,a
: eigenvalue and coefficient arrays describing system solution, ArraysActiveLayers
: vector of 1s or 0s where 1 denotes an active layer, Number or Vector
QGDipoles.Calc_ψq
— FunctionCalc_ψq(grid, a; U, ℓ, R, β, x₀=[0, 0], α=0)
Calculate $ψ$ and $q$ in a layered QG model using coefficients and vortex parameters
Arguments:
grid
: grid structure containing x, y, and Krsqa
: M x N array of coefficients, Array
Keyword arguments:
- (
U
,ℓ
): vortex speed and radius, Numbers (default:(1, 1)
) - (
R
,β
): Rossby radii and (y) PV gradients in each layer, Numbers or Vectors (default:(Inf, 0)
) x₀
: position of vortex center, vector (default:[0, 0]
)α
: initial angle of vortex, Number (default:0
)
SQG
QGDipoles.BuildLinSysSQG
— FunctionBuildLinSysSQG(M, λ, μ; tol=1e-6)
Builds the terms in the inhomogeneous eigenvalue problem; $A$, $B$, $c$, $d$ for the SQG problem
Arguments:
M
: number of coefficient to solve for, Integerλ
: ratio of vortex radius to Rossby radius in each layer, Number or Vectorμ
: nondimensional (y) vorticity gradient in each layer, Number or Vector
Keyword arguments:
tol
: error tolerance forQuadGK
viaJJ_int
, Number (default:1e-6
)
QGDipoles.Calc_ψb
— FunctionCalc_ψb(grid, a; U, ℓ, R, β, x₀=[0, 0], α=0)
Calculate SQG fields $ψ$ and $b$ using coefficients and vortex parameters
Arguments:
grid
: grid structure containingx
,y
, andKrsq
a
: M x 1 array of coefficients, Array
Keyword arguments:
- (
U
,ℓ
): vortex speed and radius, Numbers (default:1
) R
: vector of $[R, R']$, Vector (default:[Inf, Inf]
)β
: beta-plane (y) PV gradient, Number (default:1
)x₀
: position of vortex center, vector (default:[0, 0]
)α
: initial angle of vortex, Number (default:0
)
Note: Here R is the baroclinic Rossby radius, R = NH/f, and R' = R₀²/R where R₀ is the barotropic Rossby radius, R₀ = √(gH)/f. For infinite depth, R' = g/(fN).
Shared
QGDipoles.SolveInhomEVP
— FunctionSolveInhomEVP(A, B, c, d; K₀=nothing, a₀=nothing, tol=1e-6, method=:eigensolve, m=2, warn=true)
Solves the inhomogeneous eigenvalue problem using nonlinear root finding
Arguments:
A
,B
,c
,d
: inhomogeneous eigenvalue problem terms, Arrays
Keyword arguments:
K₀
,a₀
: initial guesses for $K$ and $a$, Arrays or nothings (default:nothing
)tol
: error tolerance fornlsolve
, Number (default:1e-6
)method
::eigensolve
or:nlsolve
, for $N > 1$:nlsolve
is used automatically (default::eigensolve
)m
: exponent of $K$ in eignevalue problem (default:2
)warn
: iftrue
displays warning if solution includes unextracted passive layers (default:true
)
Extras
Monopoles
QGDipoles.CreateRankine
— FunctionCreateRankine(grid; ℓ=1, Γ=2π, x₀=[0, 0])
Calculates the Rankine vortex for a 1 layer system. This vortex appears as a point vortex in the far field but consists of solid body rotation within the region $r < ℓ$.
Arguments:
grid
: grid structure containingx
,y
, andKrsq
Keyword arguments:
ℓ
: vortex speed and radius, Numbers (default:1
)Γ
: vortex circulation (default:2π
)x₀
: position of vortex center, vector (default:[0, 0]
)
Note: This function outputs $(u, v)$ directly since the solution has discontinuous velocity at the vortex boundary, $r = ℓ$, so derivatives evaluated with Fourier transforms exhibit Gibbs phenomenon.
QGDipoles.Create1LMonopole
— FunctionCreate1LMonopole(grid; ℓ=1, Γ=2π, R=Inf, x₀=[0, 0])
Calculates a monopolar vortex satisfying a Long's model assumption $q = F(ψ)$ where $q = [∇²-1/R²]ψ$. We take $F(z) = -(K²+1/R²)(z-z₀)$ for $r < ℓ$ and $F(z) = 0$ for $r > ℓ$ and $z₀ = ψ(r=ℓ)$. These solutions exist only on an f-plane $(β = 0)$.
Arguments:
grid
: grid structure containingx
,y
, andKrsq
Keyword arguments:
ℓ
: vortex speed and radius, Numbers (default:1
)Γ
: vortex circulation (default:2π
)R
: Rossby radius (default:Inf
)x₀
: position of vortex center, vector (default:[0, 0]
)
Note: This vortex has a continuous vorticity distribution so calculating (u, v)
from ψ
with Fourier transforms will work. This function outputs (u, v)
from the analytical expressions for consistency with CreateRankine
.
QGDipoles.CreateLQGMonopole
— FunctionCreateLQGMonopole(grid; ℓ=1, E=1, R=Inf, x₀=[0, 0])
Calculates a monopolar vortex in the LQG model using a numerical approach. We assume that $qⱼ + βⱼ = Fⱼ(ψⱼ + Uy)$ and write $Fⱼ(z) = -Kⱼ² z + Eⱼ$. Expanding the expression gives $qⱼ + βⱼ = -Kⱼ²(ψⱼ + Uy) + Eⱼ$ which by linearity can be split into a dipole equation $qⱼ + βⱼ = -Kⱼ²(ψⱼ + Uy)$ and a monopole equation $qⱼ = Eⱼ$. Outside the vortex, we take $qⱼ = 0$.
Arguments:
grid
: grid structure containingx
,y
, andKrsq
Keyword arguments:
ℓ
: vortex speed and radius, Numbers (default:1
)E
: vector of $Eⱼ$ values, Number or Vector (default:[1, ... , 1]
)R
: Rossby radius (default:Inf
)x₀
: position of vortex center, vector (default:[0, 0]
)
QGDipoles.InvertVorticity1LQG
— FunctionInvertVorticity1LQG(grid, q; R=Inf)
This function inverts the potential vorticity relation $q = [∇²-1/R²]ψ$ for 1-layer QG
Arguments:
grid
: grid structure containingx
,y
, andKrsq
q
: potential vorticity field, Array
Keyword arguments:
R
: Rossby radius (default:Inf
)
Note: This function is designed to be used for creating periodic streamfunctions using the vorticity fields generated by Create1LMonopole
. It does not support multi-layer QG and is only valid on an f-plane (β = 0).
QGDipoles.InvertVorticityLQG
— FunctionInvertVorticityLQG(grid, q; R=Inf)
This function inverts the potential vorticity relation $q = ΔN ψ$ for the LQG model
Arguments:
grid
: grid structure containingx
,y
, andKrsq
q
: potential vorticity field, Array
Keyword arguments:
R
: Rossby radius, Number or Vector (default:Inf
)
Internal
QGDipoles.A_func
— FunctionA_func(ξ, λ, μ)
Evaluates the matrix function $A(ξ, λ, μ) = K(ξ) [K(ξ) + D(μ)]⁻¹ ξ⁻¹$
Arguments:
ξ
: point in $[0, ∞)$, Numberλ
: ratio of vortex radius to Rossby radius in each layer, Number or Vectorμ
: nondimensional (y) vorticity gradient in each layer, Number or Vector
QGDipoles.B_func
— FunctionB_func(ξ, λ, μ)
Evaluates the matrix function $B(ξ, λ, μ) = [K(ξ) + D(μ)]⁻¹ ξ⁻¹$
Arguments:
ξ
: point in $[0, ∞)$, Numberλ
: ratio of vortex radius to Rossby radius in each layer, Number or Vectorμ
: nondimensional (y) vorticity gradient in each layer, Number or Vector
QGDipoles.JJ_int
— FunctionJJ_int(F, j, k, tol=1e-6)
Evaluates the integral $I = ∫_0^∞ F(ξ) J_{2j+2}(ξ) J_{2k+2}(ξ) \mathrm{d}ξ$
Arguments:
F
: function to integrate, typicallyA_func
orB_func
, Functionj
: first Bessel function index, Integerk
: second Bessel function index, Integertol
: error tolerance forQuadGK
, Number (default:1e-6
)
Note: This integral is performed by deforming the contour of integration into the complex plane where the Bessel function decays exponentially in the imaginary direction.
QGDipoles.InhomEVP_F!
— FunctionInhomEVP_F!(F, J, x, A, B, c, d, e)
Calculates the function $F$ and it's derivatives, $J$, at a given point $x$
Arguments:
F
,J
: values of $F$ and it's derivatives, updated by functionx
: evaluation point, ArrayA
,B
,c
: inhomogeneous eigenvalue problem terms, Arrayse
: basis spanning the space perpendicular to the $d[n]$, Array
QGDipoles.OrthogSpace
— FunctionOrthogSpace(v)
Extends the input to an orthonormal basis over $R^n$ using the Gram-Schmidt method
Arguments:
v
: array with vectors as columns, Array
QGDipoles.ZernikeR
— FunctionZernikeR(n, x)
Define the Zernike radial function using the jacobi
function from SpecialFunctions
Arguments:
n
: order, Integerx
: evaluation point, Number or Array
Note: this function is defined on $[-1, 1]$ and is set to $0$ for $|x| > 1$
QGDipoles.GridStruct
— TypeGridStruct
Stores the grid variables in physical and Fourier space
Fields:
x
,y
: x and y points in physical space, Rangeskr
,l
: x and y points in Fourier space, ArraysKrsq
:kr²+l²
in Fourier space, Array
QGDipoles.ΔNCalc
— FunctionΔNCalc(K², R, β, U=1)
Defines the $Δ_N(β)$ matrix used to invert for $ψ$ and $q$
Arguments:
K²
: value of $k²+l²$ in Fourier space, Array- (
R
,β
): Rossby radii and (y) PV gradients in each layer, Numbers or Vectors U
: vortex speed, Number (default:1
)
QGDipoles.CartesianGrid
— FunctionCartesianGrid(grid)
Formats the $(x, y)$ ranges from grid
as two-dimensional Arrays
Arguments:
grid
: grid structure containingkr
andl
QGDipoles.PolarGrid
— FunctionPolarGrid(x, y, x₀)
Calculates the polar coordinates from (x
, y
) as two-dimensional Array centred on x₀
Arguments:
x
,y
: 2D Arrays for $x$ and $y$, created usingCartesianGrid
x₀
: Vector
QGDipoles.AreaInteg2
— FunctionAreaInteg2(f, grid)
Calculates the integral $I = ∫_A f^2 \mathrm{d}A$ where $A$ is the 2D domain described by grid
.
Arguments:
f
: input Array in real or Fourier spacegrid
: grid structure
Note: f can be entered in real space or Fourier space, we use the rfft function to calculate the Fourier transform so array sizes can distinguish the two.
Base & RecipesBase
Base.summary
— FunctionBase.summary function for custom type GridStruct
Base.summary
Summary method for custom type LQGParams
Base.summary
Summary method for custom type LQGVortex
Base.summary
Summary method for custom type SQGParams
Base.summary
Summary method for custom type SQGVortex
Base.show
— FunctionBase.show function for custom type GridStruct
Base.show
Show method for custom type LQGParams
Base.show
Show method for custom type LQGVortex
Base.show
Show method for custom type SQGParams
Base.show
Show method for custom type SQGVortex
RecipesBase.apply_recipe
— Functionf(grid, F::AbstractArray; layer=1)
Recipe for plotting field F
on grid
Arguments:
grid
: grid objectF
: field, may have mltiple layers
Keyword arguments:
layer
: layer to plot (default:1
)
f(grid, F::CuArray; layer=1)
Recipe for plotting field F
on grid
Arguments:
grid
: grid objectF
: field, may have mltiple layers
Keyword arguments:
layer
: layer to plot (default:1
)