Functions

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

QGDipoles.jl

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

JJ_integ.jl

QGDipoles.A_funcFunction

Function: 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

Function: 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

Function: 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, 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

lin_sys.jl

QGDipoles.BuildLinSysFunction

Function: 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 Vector
  • tol: error tolerance for QuadGK via JJ_int, Number (default: 1e-6)
  • sqg: false; creates layered QG system, true; creates SQG system (default: false)
source
QGDipoles.ApplyPassiveLayersFunction

Function: 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

Function: 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.SolveInhomEVPFunction

Function: 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, Arrays
  • K₀, a₀: initial guesses for K and a, Arrays or nothings (default: nothing)
  • tol: error tolerance for nlsolve, Number (default: 1e-6)
  • method: 0 - eigensolve for N = 1 and nlsolve for N > 1, 1 - nlsolve (default: 0)
  • m: exponent of K in eignevalue problem (default: 2)
  • sqg: false, uses m value specified; true, sets m=1 (default: false)
  • warn: if true 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.

source
QGDipoles.InhomEVP_F!Function

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

Function: 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

create_modon.jl

QGDipoles.ZernikeRFunction

Function: 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

Structure: GridStruct

Stores the grid variables in physical and Fourier space

Arguments:

  • 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.CreateGridFunction

Function: 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
  • cuda: true; use CUDA CuArray for fields (default: false)
source
QGDipoles.Calc_ψqFunction

Function: 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 Krsq
  • x₀: position of vortex center, vector (default: [0, 0])
  • α: initial angle of vortex, Number (default: 0)
source
QGDipoles.Calc_ψbFunction

Function: 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, Number
  • grid: grid structure containing x, y, and Krsq
  • 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
QGDipoles.Calc_uvFunction

Function: Calc_uv(ψ, grid)

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

Arguments:

  • ψ: streamfunction, Array
  • grid: grid structure containing kr and l
source
QGDipoles.ΔNCalcFunction

Function: Δ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.CreateModonLQGFunction

Function: 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 Krsq
  • M: 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 to QuadGK and NLSolve functions, Number (default: 1e-6)

Note: provide values of K₀ and a₀ for active layers ONLY.

source
QGDipoles.CreateModonSQGFunction

Function: 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 Krsq
  • M: 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 to QuadGK and NLSolve 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).

source
QGDipoles.CreateLCDFunction

Function: 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.

source
QGDipoles.CreateLRDFunction

Function: 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.

source
QGDipoles.Eval_ψ_SQGFunction

Function: 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 using Calc_ψb or CreateModonSQG
  • 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

Function: 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 using Calc_ψb or CreateModonSQG
  • 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

Function: 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 using Calc_ψb or CreateModonSQG
  • 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

Function: 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 using Calc_ψb or CreateModonSQG
  • 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
QGDipoles.Calc_∇Function

Function: Calc_∇(f, grid)

Calculate the gradient ∇f for a given field f

Arguments:

  • f: function, Array
  • grid: grid structure containing kr and l
source
QGDipoles.CartesianGridFunction

Function: CartesianGrid(grid)

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

Arguments:

  • grid: grid structure containing kr and l
source
QGDipoles.PolarGridFunction

Function: 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

energetics.jl

QGDipoles.EnergyLQGFunction

Function: 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 CuArray
  • R: Rossby radius in each layer, Number or Vector
  • H: Thickness of each layer, Number or Vector
source
QGDipoles.EnstrophyLQGFunction

Function: 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
  • H: Thickness of each layer, Number or Vector
source
QGDipoles.EnergySQGFunction

Function: 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
  • 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
QGDipoles.AreaInteg2Function

Function: 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 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

vortex_types.jl

QGDipoles.LQGParamsType

Structure: LQGParams

Stores the parameters for an LQG dipolar vortex solution

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
  • CalcEnergy: flag to determine if energy is calculated
  • CalcEnstrophy: flag to determine if enstrophy is calculated
source
QGDipoles.SQGParamsType

Structure: SQGParams

Stores the parameters for an SQG dipolar vortex solution

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
  • CalcEnergy: flag to determine if energy is calculated
source
QGDipoles.LQGVortexType

Structure: LQGVortex

Stores fields and diagnostics for an LQG dipolar vortex solution

Arguments:

  • params: Vortex params
  • ψ: streamfunction
  • q: potential vorticity anomaly
  • K: eigenvalue
  • a: coefficient matrix
  • u: x velocity
  • v: y velocity
  • KE: kinetic energy
  • PE: potential energy
  • EN: enstrophy
source
QGDipoles.SQGVortexType

Structure: SQGVortex

Stores fields and diagnostics for an SQG dipolar vortex solution

Arguments:

  • params: Vortex params
  • ψ: surface streamfunction
  • b: surface buoyancy
  • K: eigenvalue
  • a: coefficient matrix
  • u: x velocity
  • v: y velocity
  • E: domain integrated energy
  • SPE: surface potential energy
source
QGDipoles.DefLQGParamsFunction

Function: DefLQGParams

Defines an LQGParams structure using the given inputs

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
  • CalcEnergy: flag to determine if energy is calculated
  • CalcEnstrophy: flag to determine if enstrophy is calculated
source
QGDipoles.DefSQGParamsFunction

Function: DefSQGParams

Defines an SQGParams structure using the given inputs

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
  • CalcEnergy: flag to determine if energy is calculated
source
QGDipoles.DefLQGVortexFunction

Function: DefLQGVortex

Defines an LQGVortex solution structure using the given inputs

Arguments:

  • grid: grid structure
  • 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
  • CalcEnergy: flag to determine if energy is calculated
  • CalcEnstrophy: flag to determine if enstrophy is calculated
source

Function: DefLQGVortex

Defines an LQGVortex solution structure using the given inputs

Arguments:

  • grid: grid structure
  • params: vortex parameters, LQGParams structure
source
QGDipoles.DefSQGVortexFunction

Function: DefSQGVortex

Defines an SQGVortex solution structure using the given inputs

Arguments:

  • grid: grid structure
  • 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
  • CalcEnergy: flag to determine if energy is calculated
source

Function: DefSQGVortex

Defines an SQGVortex solution structure using the given inputs

Arguments:

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

monopoles.jl

QGDipoles.CreateRankineFunction

Function: 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: )
  • 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.CreateMonopoleFunction

Function: 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: )
  • 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.InvertVorticity1LQGFunction

Function: 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
  • R: 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).

source

Base

Base.summaryFunction

Base.summary function for custom type GridStruct

source

Base.summary function for custom type LQGParams

source

Base.summary function for custom type SQGParams

source

Base.summary function for custom type LQGVortex

source

Base.summary function for custom type SQGVortex

source
Base.showFunction

Base.show function for custom type GridStruct

source

Base.show function for custom type LQGParams

source

Base.show function for custom type SQGParams

source

Base.show function for custom type LQGVortex

source

Base.show function for custom type SQGVortex

source