Functions
This page lists all Modules, Functions and Structures available in this package.
QGDipoles.jl
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.
JJ_integ.jl
QGDipoles.A_func
— FunctionFunction: A_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
— FunctionFunction: B_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
— FunctionFunction: JJ_int(F, j, k, tol=1e-6)
Evaluates the integral $I = \int_0^\infty F(\xi) J_{2j+2}(\xi) J_{2k+2}(\xi) \mathrm{d}\xi$
Arguments:
F
: function to integrate, typicallyA_func
orB_func
, Functionj
: first Bessel function index, Integerk
: second Bessel function index, Integertol
: error tolerance for QuadGK, 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.
lin_sys.jl
QGDipoles.BuildLinSys
— FunctionFunction: BuildLinSys(M, λ, μ; tol=1e-6, sqg=false)
Builds the terms in the inhomogeneous eigenvalue problem; A, B, c, d
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 Vectortol
: error tolerance for QuadGK viaJJ_int
, Number (default:1e-6
)sqg
:false
; creates layered QG system,true
; creates SQG system (default:false
)
QGDipoles.ApplyPassiveLayers
— FunctionFunction: ApplyPassiveLayers(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
— FunctionFunction: IncludePassiveLayers(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.SolveInhomEVP
— FunctionFunction: SolveInhomEVP(A, B, c, d; K₀=nothing, a₀=nothing, tol=1e-6, method=0, m=2, sqg=false, warn=true)
Solves the inhomogeneous eigenvalue problem using nonlinear root finding
Arguments:
A
,B
,c
,d
: inhomogeneous eigenvalue problem terms, ArraysK₀
,a₀
: initial guesses for K and a, Arrays or nothings (default:nothing
)tol
: error tolerance fornlsolve
, Number (default:1e-6
)method
:0
- eigensolve for N = 1 andnlsolve
for N > 1,1
-nlsolve
(default:0
)m
: exponent of K in eignevalue problem (default:2
)sqg
:false
, usesm
value specified;true
, setsm=1
(default:false
)warn
: iftrue
displays warning if solution includes unextracted passive layers (default:true
)
Note: setting sqg=true
overwrites the value of m
and is equivalent to setting m=1
. The option to set both is included for consistency with BuildLinSys
and more generality with the value of m
.
QGDipoles.InhomEVP_F!
— FunctionFunction: InhomEVP_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
— FunctionFunction: OrthogSpace(v)
Extends the input to an orthonormal basis over R^n using the Gram-Schmidt method
Arguments:
v
: array with vectors as columns, Array
create_modon.jl
QGDipoles.ZernikeR
— FunctionFunction: ZernikeR(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
— TypeStructure: GridStruct
Stores the grid variables in physical and Fourier space
Arguments:
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.CreateGrid
— FunctionFunction: CreateGrid(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 Numberscuda
:true
; use CUDA CuArray for fields (default:false
)
QGDipoles.Calc_ψq
— FunctionFunction: Calc_ψq(a, U, ℓ, R, β, grid, x₀=[0, 0], α=0)
Calculate ψ and q in a layered QG model using coefficients and vortex parameters
Arguments:
a
: M x N array of coefficients, Array- (
U
,ℓ
): vortex speed and radius, Numbers - (
R
,β
): Rossby radii and (y) PV gradients in each layer, Numbers or Vectors grid
: grid structure containing x, y, and Krsqx₀
: position of vortex center, vector (default:[0, 0]
)α
: initial angle of vortex, Number (default:0
)
QGDipoles.Calc_ψb
— FunctionFunction: Calc_ψb(a, U, ℓ, R, β, grid, x₀=[0, 0], α=0)
Calculate SQG fields ψ and b using coefficients and vortex parameters
Arguments:
a
: M x 1 array of coefficients, Array- (
U
,ℓ
): vortex speed and radius, Numbers R
: vector of [R, R'], Vectorβ
: beta-plane (y) PV gradient, Numbergrid
: grid structure containing x, y, and Krsqx₀
: 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).
QGDipoles.Calc_uv
— FunctionFunction: Calc_uv(ψ, grid)
Calculate the velocity fields from ψ using (u, v) = (-∂ψ/∂y, ∂ψ/∂x)
Arguments:
ψ
: streamfunction, Arraygrid
: grid structure containing kr and l
QGDipoles.ΔNCalc
— FunctionFunction: Δ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.CreateModonLQG
— FunctionFunction: CreateModonLQG(grid, M, U=1, ℓ=1, R=1, β=0, ActiveLayers=1, x₀=[0, 0], α=0; K₀=nothing, a₀=nothing, tol=1e-6)
High level wrapper function for calculating ψ and q for the Layered QG model using given parameters
Arguments:
grid
: grid structure containing x, y, and KrsqM
: number of coefficient to solve for, Integer (default:8
)- (
U
,ℓ
): vortex speed and radius, Numbers (default: (1
,1
)) - (
R
,β
): Rossby radii and (y) PV gradients in each layer, Numbers or Vectors, (default: (1
,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)K₀
,a₀
: initial guesses for K and a, Arrays or nothings (default:nothing
)tol
: error tolerance passed toQuadGK
andNLSolve
functions, Number (default:1e-6
)
Note: provide values of K₀ and a₀ for active layers ONLY.
QGDipoles.CreateModonSQG
— FunctionFunction: CreateModonSQG(grid, M, U=1, ℓ=1, R=[Inf, Inf], β=0, x₀=[0, 0], α=0; K₀=nothing, a₀=nothing, tol=1e-6)
High level wrapper function for calculating ψ and b for the SQG model using given parameters
Arguments:
grid
: grid structure containing x, y, and KrsqM
: number of coefficient to solve for, Integer (default:12
)- (
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
)K₀
,a₀
: initial guesses for K and a, Arrays or nothings (default:nothing
)tol
: error tolerance passed toQuadGK
andNLSolve
functions, Number (default:1e-6
)
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.CreateLCD
— FunctionFunction: CreateLCD(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 containing x, y, and Krsq- (
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
— FunctionFunction: CreateLRD(grid, U=1, ℓ=1, R=1, β=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 containing x, y, and Krsq- (
U
,ℓ
): vortex speed and radius, Numbers (default: (1
,1
)) - (
R
,β
): Rossby radii and (y) PV gradient, Numbers, (default: (1
,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.
QGDipoles.Eval_ψ_SQG
— FunctionFunction: Eval_ψ_SQG(grid, ψ, z=[0], U=1, R=[Inf, Inf], β=0)
Evaluates ψ at specified depths, z in [-R, 0], for the SQG problem
Arguments:
grid
: grid structure containing x, y, and Krsqψ
: surface streamfunction, calculated usingCalc_ψb
orCreateModonSQG
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
— FunctionFunction: Eval_q_SQG(grid, ψ, z=[0], U=1, R=[Inf, Inf], β=0)
Evaluates q at specified depths, z in [-R, 0], for the SQG problem
Arguments:
grid
: grid structure containing x, y, and Krsqψ
: surface streamfunction, calculated usingCalc_ψb
orCreateModonSQG
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
— FunctionFunction: Eval_b_SQG(grid, ψ, z=[0], U=1, R=[Inf, Inf], β=0)
Evaluates b at specified depths, z in [-R, 0], for the SQG problem
Arguments:
grid
: grid structure containing x, y, and Krsqψ
: surface streamfunction, calculated usingCalc_ψb
orCreateModonSQG
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
— FunctionFunction: Eval_w_SQG(grid, ψ, z=[0], U=1, R=[Inf, Inf], β=0)
Evaluates N²w at specified depths, z in [-R, 0], for the SQG problem using N²w = -J[ψ + Uy, b]
Arguments:
grid
: grid structure containing x, y, and Krsqψ
: surface streamfunction, calculated usingCalc_ψb
orCreateModonSQG
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' = ∞.
QGDipoles.Calc_∇
— FunctionFunction: Calc_∇(f, grid)
Calculate the gradient ∇f for a given field f
Arguments:
f
: function, Arraygrid
: grid structure containing kr and l
QGDipoles.CartesianGrid
— FunctionFunction: CartesianGrid(grid)
Formats the (x, y) ranges from grid
as two-dimensional Arrays
Arguments:
grid
: grid structure containing kr and l
QGDipoles.PolarGrid
— FunctionFunction: PolarGrid(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
energetics.jl
QGDipoles.EnergyLQG
— FunctionFunction: EnergyLQG(grid, ψ, R, H=[1])
Calculates the kinetic and potential energy for the LQG system
Arguments:
grid
: grid structure containing Krsqψ
: streamfunction in each layer, Array or CuArrayR
: Rossby radius in each layer, Number or VectorH
: Thickness of each layer, Number or Vector
QGDipoles.EnstrophyLQG
— FunctionFunction: EnstrophyLQG(grid, q, H=[1])
Calculates the enstrophy for the LQG system
Arguments:
grid
: grid structure containing Krsqq
: potential vorticity anomaly in each layer, Array or CuArrayH
: Thickness of each layer, Number or Vector
QGDipoles.EnergySQG
— FunctionFunction: EnergySQG(grid, ψ, b, R′)
Calculates the energies for the SQG system; the total domain integrated energy and the surface potential energy
Arguments:
grid
: grid structure containing Krsqψ
: surface streamfunction, Array or CuArrayb
: surface buoyancy, , Array or CuArrayR′
: reduced barotropic Rossby radius, Number (default:Inf
)
Note: the surface potential energy is sometimes referred to as the generalised enstrophy or the buoyancy variance.
QGDipoles.AreaInteg2
— FunctionFunction: AreaInteg2(f, grid)
Calculates the integral $I = \int_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.
vortex_types.jl
QGDipoles.LQGParams
— TypeStructure: LQGParams
Stores the parameters for an LQG dipolar vortex solution
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 calculatedCalcEnergy
: flag to determine if energy is calculatedCalcEnstrophy
: flag to determine if enstrophy is calculated
QGDipoles.SQGParams
— TypeStructure: SQGParams
Stores the parameters for an SQG dipolar vortex solution
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 calculatedCalcEnergy
: flag to determine if energy is calculated
QGDipoles.LQGVortex
— TypeStructure: LQGVortex
Stores fields and diagnostics for an LQG dipolar vortex solution
Arguments:
params
: Vortex paramsψ
: streamfunctionq
: potential vorticity anomalyK
: eigenvaluea
: coefficient matrixu
: x velocityv
: y velocityKE
: kinetic energyPE
: potential energyEN
: enstrophy
QGDipoles.SQGVortex
— TypeStructure: SQGVortex
Stores fields and diagnostics for an SQG dipolar vortex solution
Arguments:
params
: Vortex paramsψ
: surface streamfunctionb
: surface buoyancyK
: eigenvaluea
: coefficient matrixu
: x velocityv
: y velocityE
: domain integrated energySPE
: surface potential energy
QGDipoles.DefLQGParams
— FunctionFunction: DefLQGParams
Defines an LQGParams structure using the given inputs
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 calculatedCalcEnergy
: flag to determine if energy is calculatedCalcEnstrophy
: flag to determine if enstrophy is calculated
QGDipoles.DefSQGParams
— FunctionFunction: DefSQGParams
Defines an SQGParams structure using the given inputs
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 calculatedCalcEnergy
: flag to determine if energy is calculated
QGDipoles.DefLQGVortex
— FunctionFunction: DefLQGVortex
Defines an LQGVortex solution structure using the given inputs
Arguments:
grid
: grid structureU
: 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 calculatedCalcEnergy
: flag to determine if energy is calculatedCalcEnstrophy
: flag to determine if enstrophy is calculated
Function: DefLQGVortex
Defines an LQGVortex solution structure using the given inputs
Arguments:
grid
: grid structureparams
: vortex parameters, LQGParams structure
QGDipoles.DefSQGVortex
— FunctionFunction: DefSQGVortex
Defines an SQGVortex solution structure using the given inputs
Arguments:
grid
: grid structureU
: 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 calculatedCalcEnergy
: flag to determine if energy is calculated
Function: DefSQGVortex
Defines an SQGVortex solution structure using the given inputs
Arguments:
grid
: grid structureparams
: vortex parameters, LQGParams structure
monopoles.jl
QGDipoles.CreateRankine
— FunctionFunction: CreateRankine(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 containing x, y, and Krsqℓ
: 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.CreateMonopole
— FunctionFunction: CreateMonopole(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 containing x, y, and Krsqℓ
: 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.InvertVorticity1LQG
— FunctionFunction: InvertVorticity1LQG(grid, q, R=Inf)
This function inverts the potential vorticity relation q = [∇²-1/R²]ψ for 1-layer QG
Arguments:
grid
: grid structure containing x, y, and Krsqq
: potential vorticity field, ArrayR
: Rossby radius (default:Inf
)
Note: This function is designed to be used for creating periodic streamfunctions using the vorticity fields generated by CreateMonopole
. It does not support multi-layer QG and is only valid on an f-plane (β = 0).
Base
Base.summary
— FunctionBase.summary function for custom type GridStruct
Base.summary function for custom type LQGParams
Base.summary function for custom type SQGParams
Base.summary function for custom type LQGVortex
Base.summary function for custom type SQGVortex
Base.show
— FunctionBase.show function for custom type GridStruct
Base.show function for custom type LQGParams
Base.show function for custom type SQGParams
Base.show function for custom type LQGVortex
Base.show function for custom type SQGVortex