Poly
Polynomial functions and conveniences that are needed for the spectral or spatial discretization.
Polynomials
Limace.Poly.jacobi
— Functionjacobi(n,a,b,x)
Jacobi polynomial
\[J_n^{(a,b)}(x)\]
Limace.Poly.ylm
— Functionylm(ℓ::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.
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.adamgaunt
— Functionadamgaunt(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}$.
Limace.Poly.elsasser
— Functionelsasser(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}$.
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.inners
— Functioninners(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.
Limace.Poly.innert
— Functioninnert(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.
These functions are used for the discrete integration in radius. The derivatives are computed using automatic differentiation (based on ForwardDiff.jl).
Limace.Quadrature.rquad
— Functionrquad(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)
The radial integrals of the projections are all integrated using
Limace.Quadrature.∫dr
— Function∫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 ∫.
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.p
— Functionp(l) -> Any
$p(l) = l(l+1)$ following notation of Ivers and Phillips (2008)
Limace.Poly._∂ll
— Function_∂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).
Limace.Poly.D
— FunctionD(f, l, r) -> Any
\[D_l(f) = \frac{\partial^2 f}{\partial r^2} +\frac{2}{r}\frac{\partial f}{\partial r} - \frac{l(l+1)}{r^2}f\]
Below equation (31) in Ivers and Phillips (2008).
Limace.Poly.dylmdθ
— Functiondylmdθ(l, m, θ, ϕ) -> Any
Derivative of spherical harmonic $Y_l^m$ in $\theta$.
Limace.Poly.dylmdϕ
— Functiondylmdϕ(l, m, θ, ϕ) -> Any
Derivative of spherical harmonic $Y_l^m$ in $\phi$.