Poly

Polynomial functions and conveniences that are needed for the spectral or spatial discretization.

Polynomials

Limace.Poly.ylmFunction
ylm(ℓ::Int64, m::Int64, θ, φ) -> Any

Spherical harmonic $Y_l^m$ in full norm, i.e.

\[\int Y_l^mY_i^j\, \sin(\theta)\,\mathrm{d}\theta\mathrm{d}\phi = \delta_{li}\delta_{mj}\]

where $\theta$ is the colatitude and $\phi$ the azimuthal angle.

source

Integral functions

Integrals over the spherical surfaces are reduced to two specific integrals, known as the Adam-Gaunt and Elsasser integrals (James, 1973).

Limace.Poly.adamgauntFunction
adamgaunt(la,lb,lc,ma,mb,mc)

Adam-Gaunt integral (James, 1973), given by

\[A_{abc} = \oint\int Y_aY_bY_c\sin\theta\,\mathrm{d}\theta\mathrm{d}\phi,\]

where the spherical harmonics are abbreviated, so that $Y_a = Y_{l_a}^{m_a}$.

source
Limace.Poly.elsasserFunction
elsasser(la,lb,lc,ma,mb,mc)

Elsasser integral (James, 1973), given by

\[E_{abc} = \oint\int Y_c\left( \frac{\partial Y_a}{\partial \theta} \frac{\partial Y_b}{\partial \phi} - \frac{\partial Y_a}{\partial \phi}\frac{\partial Y_b}{\partial \theta} \right)\,\mathrm{d}\theta\mathrm{d}\phi,\]

where the spherical harmonics are abbreviated, so that $Y_a = Y_{l_a}^{m_a}$.

source

The radial part of the projection of poloidal and toroidal vectors can be described by two scalar functions, Limace.Poly.inners and Limace.Poly.innert.

Limace.Poly.innersFunction
inners(s, s2, l, r) -> Any

\[r\rightarrow \frac{l(l+1)}{r^2}\left( l(l+1)s(r) s_2(r) + \frac{\partial r s(r)}{\partial r}\frac{\partial r s_2(r)}{\partial r}\right)\]

Radial function to be integrated in radius when computing the inner product of two poloidal vectors.

source
Limace.Poly.innertFunction
innert(t, t2, l::Int64, r) -> Any

\[r\rightarrow l(l+1) t(r) t_2(r)\]

Radial function to be integrated in radius when computing the inner product of two toroidal vectors.

source

These functions are used for the discrete integration in radius. The derivatives are computed using automatic differentiation (based on ForwardDiff.jl).

Limace.Quadrature.rquadFunction
rquad(nr, a, b)
rquad(nr, V::Volume)

Computes nr Gauss-Legendre quadrature points and weights on a radial grid.

r0 = 0.3
r1 = 1.0

r, wr = Limace.Quadrature.rquad(50, r0, r1)

V = SphericalShell(r0,r1)

r, wr = Limace.Quadrature.rquad(50, V)
source

The radial integrals of the projections are all integrated using

Limace.Quadrature.∫drFunction
∫dr(f, r::Vector{T}, wr::Vector{T})::Complex{T} where T

Computes the integral of a function f

\[\int_{r_0}^{r_1} r^2 f(r)\,\mathrm{d}r,\]

with $r_0<r<r_1$, using quadrature points r and weights wr computed by Limace.Quadrature.rquad. This function returns a complex value by default for easier type stability.

Type \int<TAB> to write ∫.

source
using Limace.Quadrature: rquad, ∫dr
r,wr = rquad(10, 0.0, 1.0)
f = r->r^2
∫dr(f,r,wr) ≈ 1/5
true

Miscellaneous functions

The functions listed here are mostly for internal convenience. Nevertheless, they may be useful for further implementations.

Limace.Poly._∂llFunction
_∂ll(f, l, l1, r) -> Any

\[\partial_l^{l_1} = \begin{cases} \frac{\partial f}{\partial r} + \frac{l+1}{r}f \quad \mathrm{if}\, l_1 = l-1\\ \frac{\partial f}{\partial r} - \frac{l}{r}f \quad \mathrm{if}\, l_1 = l+1 \end{cases}\]

Equation (25) in Ivers and Phillips (2008).

source