Functions

This page lists all Modules, Functions and Structures available in this package.

QGDipoles.jl Module

QGDipolesModule

Package 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.

source

Vortex Structures

LQG

QGDipoles.LQGParamsType
LQGParams

Stores the parameters for an LQG dipolar vortex solution

Fields:

  • U: Vortex speed
  • : Vortex radius
  • R: Rossby radius
  • β: background PV gradient
  • ActiveLayers: 1 => layer contains vortex region
  • H: thickness of each layer
  • x₀: Vortex position
  • α: Direction of vortex propagation
  • M: number of coefficients in Zernike expansion
  • tol: maximum error in solution evaluation
  • K₀: initial guess for eigenvalue
  • a₀: initial guess for coefficients
  • UseAnalytic: use analytic solution (1-layer only)
  • CalcVelocity: flag to determine if velocity is calculated
  • CalcVorticity: flag to determine if vorticity is calculated
  • CalcEnergy: flag to determine if energy is calculated
  • CalcEnstrophy: flag to determine if enstrophy is calculated
source
QGDipoles.LQGVortexType
LQGVortex

Stores fields and diagnostics for an LQG dipolar vortex solution

Fields:

  • params: Vortex params
  • ψ: streamfunction
  • q: potential vorticity anomaly
  • K: eigenvalue
  • a: coefficient matrix
  • u: x velocity
  • v: y velocity
  • ζ: vertical vorticity
  • KE: kinetic energy
  • PE: potential energy
  • EN: enstrophy
source
QGDipoles.DefLQGParamsFunction
DefLQGParams(; 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 radius
  • R: Rossby radius
  • β: background PV gradient
  • ActiveLayers: 1 => layer contains vortex region
  • H: thickness of each layer
  • x₀: Vortex position
  • α: Direction of vortex propagation
  • M: number of coefficients in Zernike expansion
  • tol: maximum error in solution evaluation
  • K₀: initial guess for eigenvalue
  • a₀: initial guess for coefficients
  • UseAnalytic: use analytic solution (1-layer only)
  • CalcVelocity: flag to determine if velocity is calculated
  • CalcVorticity: flag to determine if vorticity is calculated
  • CalcEnergy: flag to determine if energy is calculated
  • CalcEnstrophy: flag to determine if enstrophy is calculated
source
QGDipoles.DefLQGVortexFunction
DefLQGVortex(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 radius
  • R: Rossby radius
  • β: background PV gradient
  • ActiveLayers: 1 => layer contains vortex region
  • H: thickness of each layer
  • x₀: Vortex position
  • α: Direction of vortex propagation
  • M: number of coefficients in Zernike expansion
  • tol: maximum error in solution evaluation
  • K₀: initial guess for eigenvalue
  • a₀: initial guess for coefficients
  • UseAnalytic: use analytic solution (1-layer only)
  • CalcVelocity: flag to determine if velocity is calculated
  • CalcVorticity: flag to determine if vorticity is calculated
  • CalcEnergy: flag to determine if energy is calculated
  • CalcEnstrophy: flag to determine if enstrophy is calculated
source
DefLQGVortex(grid, params)

Defines an LQGVortex solution structure using the given inputs

Arguments:

  • grid: grid structure
  • params: vortex parameters, LQGParams structure
source

SQG

QGDipoles.SQGParamsType
SQGParams

Stores the parameters for an SQG dipolar vortex solution

Fields:

  • U: Vortex speed
  • : Vortex radius
  • R: Rossby radius
  • β: background PV gradient
  • x₀: Vortex position
  • α: Direction of vortex propagation
  • M: number of coefficients in Zernike expansion
  • tol: maximum error in solution evaluation
  • K₀: initial guess for eigenvalue
  • a₀: initial guess for coefficients
  • CalcVelocity: flag to determine if velocity is calculated
  • CalcVorticity: flag to determine if vorticity is calculated
  • CalcEnergy: flag to determine if energy is calculated
source
QGDipoles.SQGVortexType
SQGVortex

Stores fields and diagnostics for an SQG dipolar vortex solution

Fields:

  • params: Vortex params
  • ψ: surface streamfunction
  • b: surface buoyancy
  • K: eigenvalue
  • a: coefficient matrix
  • u: x velocity
  • v: y velocity
  • ζ: vertical vorticity
  • E: domain integrated energy
  • SPE: surface potential energy
source
QGDipoles.DefSQGParamsFunction
DefSQGParams(; 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 radius
  • R: Rossby radius
  • β: background PV gradient
  • x₀: Vortex position
  • α: Direction of vortex propagation
  • M: number of coefficients in Zernike expansion
  • tol: maximum error in solution evaluation
  • K₀: initial guess for eigenvalue
  • a₀: initial guess for coefficients
  • CalcVelocity: flag to determine if velocity is calculated
  • CalcVorticity: flag to determine if vorticity is calculated
  • CalcEnergy: flag to determine if energy is calculated
source
QGDipoles.DefSQGVortexFunction
DefSQGVortex(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 radius
  • R: Rossby radius
  • β: background PV gradient
  • x₀: Vortex position
  • α: Direction of vortex propagation
  • M: number of coefficients in Zernike expansion
  • tol: maximum error in solution evaluation
  • K₀: initial guess for eigenvalue
  • a₀: initial guess for coefficients
  • CalcVelocity: flag to determine if velocity is calculated
  • CalcVorticity: flag to determine if vorticity is calculated
  • CalcEnergy: flag to determine if energy is calculated
source
DefSQGVortex(grid, params)

Defines an SQGVortex solution structure using the given inputs

Arguments:

  • grid: grid structure
  • params: vortex parameters, LQGParams structure
source

Shared

QGDipoles.UpdateParamsFunction
UpdateParams(params::LQGParams; kwargs...)

Creates an LQGParams structure by replacing parameters in params with the given keywords

Arguments:

Keyword arguments:

source
UpdateParams(params::SQGParams; kwargs...)

Creates an SQGParams structure by replacing parameters in params with the given keywords

Arguments:

Keyword arguments:

source
QGDipoles.UpdateVortexFunction
UpdateVortex(grid, vortex::LQGVortex; kwargs...)

Creates an LQGVortex structure by replacing parameters in vortex.params with the given keywords

Arguments:

  • grid: grid structure
  • vortex: LQGVortex structure

Keyword arguments:

source
UpdateVortex(grid, vortex::SQGVortex; kwargs...)

Creates an SQGVortex structure by replacing parameters in vortex.params with the given keywords

Arguments:

  • grid: grid structure
  • vortex: SQGVortex structure

Keyword arguments:

source

High-Level Functions

LQG

QGDipoles.CreateModonLQGFunction
CreateModonLQG(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 containing x, y, and Krsq

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 to QuadGK and NLSolve 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.

source
QGDipoles.CreateLCDFunction
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

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$.

source
QGDipoles.CreateLRDFunction
CreateLRD(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 containing x, y, and Krsq

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$.

source

SQG

QGDipoles.CreateModonSQGFunction
CreateModonSQG(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 containing x, y, and Krsq

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 to QuadGK and NLSolve 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).

source
QGDipoles.Eval_ψ_SQGFunction
Eval_ψ_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 containing x, y, and Krsq
  • ψ: surface streamfunction, calculated using Calc_ψb or CreateModonSQG

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)$.

source
QGDipoles.Eval_q_SQGFunction
Eval_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 containing x, y, and Krsq
  • ψ: surface streamfunction, calculated using Calc_ψb or CreateModonSQG

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)$.

source
QGDipoles.Eval_b_SQGFunction
Eval_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 containing x, y, and Krsq
  • ψ: surface streamfunction, calculated using Calc_ψb or CreateModonSQG

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)$.

source
QGDipoles.Eval_w_SQGFunction
Eval_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 containing x, y, and Krsq
  • ψ: surface streamfunction, calculated using Calc_ψb or CreateModonSQG

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' = ∞$.

source

Utilities

QGDipoles.CreateGridFunction
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, Integers
  • Lx, 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)
source
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, Integers
  • Lx, Ly: x and y domains, either vectors of endpoints or lengths, Vectors or Numbers
  • cuda: true; use CUDA CuArray for fields (default: false)
source
QGDipoles.Calc_uvFunction
Calc_uv(grid, ψ)

Calculate the velocity fields from $ψ$ using $(u, v) = (-∂ψ/∂y, ∂ψ/∂x)$

Arguments:

  • grid: grid structure containing kr and l
  • ψ: streamfunction, Array
source
QGDipoles.Calc_∇Function
Calc_∇(grid, f)

Calculate the gradient $∇f$ for a given field $f$

Arguments:

  • grid: grid structure containing kr and l
  • f: function, Array
source
QGDipoles.Calc_ζFunction
Calc_ζ(grid, ψ)

Calculate the vertical vorticity using $ζ = ∂v/∂x - ∂u/∂y = ∇²ψ$

Arguments:

  • grid: grid structure containing Krsq
  • ψ: streamfunction, Array
source

Diagnostics

LQG

QGDipoles.EnergyLQGFunction
EnergyLQG(grid, ψ; R=Inf, H=[1])

Calculates the kinetic and potential energy for the LQG system

Arguments:

  • grid: grid structure containing Krsq
  • ψ: 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])
source
QGDipoles.EnstrophyLQGFunction
EnstrophyLQG(grid, q; H=[1])

Calculates the enstrophy for the LQG system

Arguments:

  • grid: grid structure containing Krsq
  • q: potential vorticity anomaly in each layer, Array or CuArray

Keyword arguments:

  • H: Thickness of each layer, Number or Vector (default: [1])
source

SQG

QGDipoles.EnergySQGFunction
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 CuArray
  • b: 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.

source

Low-Level Functions

LQG

QGDipoles.BuildLinSysLQGFunction
BuildLinSysLQG(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 for QuadGK via JJ_int, Number (default: 1e-6)
source
QGDipoles.ApplyPassiveLayersFunction
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, Arrays
  • ActiveLayers: vector of 1s or 0s where 1 denotes an active layer, Number or Vector
source
QGDipoles.IncludePassiveLayersFunction
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, Arrays
  • ActiveLayers: vector of 1s or 0s where 1 denotes an active layer, Number or Vector
source
QGDipoles.Calc_ψqFunction
Calc_ψ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 Krsq
  • a: 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)
source

SQG

QGDipoles.BuildLinSysSQGFunction
BuildLinSysSQG(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 for QuadGK via JJ_int, Number (default: 1e-6)
source
QGDipoles.Calc_ψbFunction
Calc_ψb(grid, a; U, ℓ, R, β, x₀=[0, 0], α=0)

Calculate SQG fields $ψ$ and $b$ using coefficients and vortex parameters

Arguments:

  • grid: grid structure containing x, y, and Krsq
  • 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).

source

Shared

QGDipoles.SolveInhomEVPFunction
SolveInhomEVP(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 for nlsolve, 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: if true displays warning if solution includes unextracted passive layers (default: true)
source

Extras

Monopoles

QGDipoles.CreateRankineFunction
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

Keyword arguments:

  • : vortex speed and radius, Numbers (default: 1)
  • Γ: vortex circulation (default: )
  • 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.

source
QGDipoles.Create1LMonopoleFunction
Create1LMonopole(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

Keyword arguments:

  • : vortex speed and radius, Numbers (default: 1)
  • Γ: vortex circulation (default: )
  • 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.

source
QGDipoles.CreateLQGMonopoleFunction
CreateLQGMonopole(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 containing x, y, and Krsq

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])
source
QGDipoles.InvertVorticity1LQGFunction
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 Krsq
  • 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).

source
QGDipoles.InvertVorticityLQGFunction
InvertVorticityLQG(grid, q; R=Inf)

This function inverts the potential vorticity relation $q = ΔN ψ$ for the LQG model

Arguments:

  • grid: grid structure containing x, y, and Krsq
  • q: potential vorticity field, Array

Keyword arguments:

  • R: Rossby radius, Number or Vector (default: Inf)
source

Internal

QGDipoles.A_funcFunction
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
source
QGDipoles.B_funcFunction
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
source
QGDipoles.JJ_intFunction
JJ_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, typically A_func or B_func, Function
  • j: first Bessel function index, Integer
  • k: second Bessel function index, Integer
  • tol: 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.

source
QGDipoles.InhomEVP_F!Function
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 function
  • x: evaluation point, Array
  • A, B, c: inhomogeneous eigenvalue problem terms, Arrays
  • e: basis spanning the space perpendicular to the $d[n]$, Array
source
QGDipoles.OrthogSpaceFunction
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
source
QGDipoles.ZernikeRFunction
ZernikeR(n, x)

Define the Zernike radial function using the jacobi function from SpecialFunctions

Arguments:

  • n: order, Integer
  • x: evaluation point, Number or Array

Note: this function is defined on $[-1, 1]$ and is set to $0$ for $|x| > 1$

source
QGDipoles.GridStructType
GridStruct

Stores the grid variables in physical and Fourier space

Fields:

  • x, y: x and y points in physical space, Ranges
  • kr, l: x and y points in Fourier space, Arrays
  • Krsq: kr²+l² in Fourier space, Array
source
QGDipoles.ΔNCalcFunction
ΔNCalc(K², R, β, U=1)

Defines the $Δ_N(β)$ matrix used to invert for $ψ$ and $q$

Arguments:

  • : 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)
source
QGDipoles.CartesianGridFunction
CartesianGrid(grid)

Formats the $(x, y)$ ranges from grid as two-dimensional Arrays

Arguments:

  • grid: grid structure containing kr and l
source
QGDipoles.PolarGridFunction
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 using CartesianGrid
  • x₀: Vector
source
QGDipoles.AreaInteg2Function
AreaInteg2(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 space
  • grid: 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.

source

Base & RecipesBase

RecipesBase.apply_recipeFunction
f(grid, F::AbstractArray; layer=1)

Recipe for plotting field F on grid

Arguments:

  • grid: grid object
  • F: field, may have mltiple layers

Keyword arguments:

  • layer: layer to plot (default: 1)
source
f(grid, F::CuArray; layer=1)

Recipe for plotting field F on grid

Arguments:

  • grid: grid object
  • F: field, may have mltiple layers

Keyword arguments:

  • layer: layer to plot (default: 1)
source