API Reference

Functions

AssociatedLegendrePolynomials.legendreFunction
p = legendre(norm::AbstractLegendreNorm, l::DimOrInd, m::DimOrInd, x)

Computes the associated Legendre polynomials $N_ℓ^m P_ℓ^m(x)$ of degree(s) l and order(s) m at x for the normalization scheme norm.

  • If l and m are integers, returns a single value.
  • If l is a range 0:lmax and m an integer, returns the vector of values for order m and all degrees 0 ≤ l ≤ lmax.
  • If l is a range 0:lmax and m is a range 0:mmax, returns the matrix of values for all degrees 0 ≤ l ≤ lmax and orders 0 ≤ m ≤ mmax.

Note that in both the second and third cases, the ranges must have a first index of 0.

AssociatedLegendrePolynomials.legendre!Function
legendre!(norm::AbstractLegendreNorm, Λ, l::Integer, m::Integer, x)

Fills the array Λ with the Legendre polynomial values $N_ℓ^m P_ℓ^m(x)$ up to/of degree(s) l and order(s) m for the normalization scheme norm. Λ must be an array with between 0 and 2 more dimensions than x, with the leading dimensions having the same shape as x.

  • If ndims(Λ) == ndims(x), then Λ is filled with the polynomial values at x for degree l and order m.
  • If ndims(Λ) == ndims(x) + 1, then l is interpreted as lmax, and Λ filled with polynomial values for all degrees 0 ≤ l ≤ lmax of order m.
  • If ndims(Λ) == ndims(x) + 2, then l is interpreted as lmax and m as mmax, and Λ is filled with polynomial values for all degrees 0 ≤ l ≤ lmax and orders 0 ≤ m ≤ min(mmax, l).
AssociatedLegendrePolynomials.NlmFunction
N = Nlm([T=Float64], l, m)

Computes the normalization constant

\[ N_\ell^m \equiv \sqrt{\frac{2\ell+1}{4\pi} \frac{(\ell-m)!}{(\ell+m)!}}\]

which defines the Spherical Harmonic normalized functions $\lambda_\ell^m(x)$ in terms of the standard unit normalized $P_\ell^m(x)$

\[ \lambda_\ell^m(x) \equiv N_\ell^m P_\ell^m(x)\]

using numbers of type T.

See also Plm and λlm.

Normalizations

AssociatedLegendrePolynomials.LegendreUnitNormType
struct LegendreUnitNorm <: AbstractLegendreNorm end

Trait type denoting the unnormalized associated Legendre functions $P_\ell^m(x)$ which solve the colatitude $\theta$ part of Laplace's equation in spherical coordinates where $x=\cos(\theta)$. The degree $\ell=0$, order $m=0$ constant function is normalized to be unity:

\[ P_0^0(x) = 1\]

AssociatedLegendrePolynomials.LegendreFourPiNormType
struct LegendreFourPiNorm <: AbstractLegendreNorm end

Trait type denoting the $4\pi$ (geodesy) normalization of the associated Legendre functions $F_\ell^m(x)$. The orthonormal normalization is defined such that

\[ \int_{-1}^{1} \left[ F_\ell^m(x) \right]^2 \,dx = 2\]

The normalization factor with respect to the standard (unnormalized) Legendre functions $P_\ell^m(x)$ (LegendreUnitNorm) is given by

\[ F_\ell^m(x) \equiv \sqrt{2\pi(2\ell+1) \frac{(\ell-m)!}{(\ell+m)!}} P_\ell^m(x)\]

AssociatedLegendrePolynomials.LegendreOrthoNormType
struct LegendreOrthoNorm <: AbstractLegendreNorm end

Trait type denoting the orthonormal (full) normalization of the associated Legendre functions $O_\ell^m(x)$. The orthonormal normalization is defined such that

\[ \int_{-1}^{1} \left[ O_\ell^m(x) \right]^2 \,dx = 1\]

The normalization factor with respect to the standard (unnormalized) Legendre functions $P_\ell^m(x)$ (LegendreUnitNorm) is given by

\[ O_\ell^m(x) \equiv \sqrt{\frac{2\ell+1}{2} \frac{(\ell-m)!}{(\ell+m)!}} P_\ell^m(x)\]

AssociatedLegendrePolynomials.LegendreSphereNormType
struct LegendreSphereNorm <: AbstractLegendreNorm end

Trait type denoting the spherical-harmonic normalization of the associated Legendre functions $\lambda_\ell^m(x)$. The spherical-harmonic normalization is defined such that

\[ \int_{-1}^{1} \left[ \lambda_\ell^m(x) \right]^2 \,dx = \frac{1}{2\pi}\]

The normalization factor with respect to the standard (unnormalized) Legendre functions $P_\ell^m(x)$ (LegendreUnitNorm) is given by

\[ \lambda_\ell^m(x) \equiv \sqrt{\frac{2\ell+1}{4\pi} \frac{(\ell-m)!}{(\ell+m)!}} P_\ell^m(x)\]

AssociatedLegendrePolynomials.LegendreNormCoeffType
struct LegendreNormCoeff{N<:AbstractLegendreNorm,T<:Real} <: AbstractLegendreNorm

Precomputed recursion relation coefficients for the normalization N and value type T.

Example

julia> LegendreNormCoeff{LegendreSphereNorm,Float64}(1)
LegendreNormCoeff{LegendreSphereNorm,Float64} for lmax = 1, mmax = 1 with coefficients:
    μ: [0.0, 1.22474]
    ν: [1.73205, 2.23607]
    α: [0.0 0.0; 1.73205 0.0]
    β: [0.0 0.0; -0.0 0.0]

Aliases

The following functors are constant aliases to the underlying normalization types which have been made callable via an appropriate call overload.

AssociatedLegendrePolynomials.PlmConstant
p = Plm(l, m, x)

Computes the associated Legendre polynomials using unit normalization; equivalent to p = legendre(LegendreUnitNorm(), l, m, x).

AssociatedLegendrePolynomials.λlmConstant
λ = λlm(l, m, x)

Computes the associated Legendre polynomials using spherical-harmonic normalization; equivalent to λ = legendre(LegendreSphereNorm(), l, m, x).

AssociatedLegendrePolynomials.Plm!Constant
Plm!(P, l, m, x)

Fills the array P with the unit-normalized associated Legendre polynomial values $P_ℓ^m(x)$; equivalent to legendre!(LegendreUnitNorm(), P, l, m, x).

AssociatedLegendrePolynomials.λlm!Constant
λlm!(Λ, l, m, x)

Fills the array Λ with the spherical-harmonic normalized associated Legendre polynomial values $λ_ℓ^m(x)$; equivalent to legendre!(LegendreSphereNorm(), P, l, m, x).

There are also aliases for pre-computed coefficients of the provided normalizations.

Normalization Interface

The following functions are unexported but considered a public API for interacting with and defining additional normalizations.

AssociatedLegendrePolynomials.initcondFunction
initcond(::N, ::Type{T}) where {N<:AbstractLegendreNorm, T}

Returns the initial condition $P_0^0(x)$ for the associated Legendre recursions based on the normalization choice N for numeric type T.

AssociatedLegendrePolynomials.coeff_μFunction
coeff_μ(norm::N, ::Type{T}, l::Integer) where {N<:AbstractLegendreNorm, T}

Returns the coefficient $\mu_\ell$ for the single-term recursion relation

\[ P_\ell^\ell(x) = -\mu_\ell \sqrt{1-x^2} P_{\ell-1}^{\ell-1}(x)\]

where $\mu_\ell$ is appropriate for the choice of normalization N.

AssociatedLegendrePolynomials.coeff_νFunction
coeff_ν(norm::N, ::Type{T}, l::Integer) where {N<:AbstractLegendreNorm, T}

Returns the coefficient $\nu_\ell$ for the single-term recursion relation

\[ P_\ell^{\ell-1}(x) = \nu_\ell x P_{\ell-1}^{\ell-1}(x)\]

where $\nu_\ell$ is appropriate for the choice of normalization N.

AssociatedLegendrePolynomials.coeff_αFunction
coeff_α(norm::N, ::Type{T}, l::Integer, m::Integer) where {N<:AbstractLegendreNorm, T}

Returns the coefficient $\alpha_\ell^m$ for the two-term recursion relation

\[ P_\ell^m(x) = \alpha_\ell^m x P_{\ell-1}^m(x) - \beta_\ell^m P_{\ell-2}^m(x)\]

where $\alpha_\ell^m$ is appropriate for the choice of normalization N.

AssociatedLegendrePolynomials.coeff_βFunction
coeff_β(norm::N, ::Type{T}, l::Integer, m::Integer) where {N<:AbstractLegendreNorm, T}

Returns the coefficient $\beta_\ell^m$ for the two-term recursion relation

\[ P_\ell^m(x) = \alpha_\ell^m x P_{\ell-1}^m(x) - \beta_\ell^m P_{\ell-2}^m(x)\]

where $\beta_\ell^m$ is appropriate for the choice of normalization N.

AssociatedLegendrePolynomials.boundscheck_hookFunction
boundscheck_hook(norm::AbstractLegendreNorm, lmax, mmax)

A bounds-checking hook executed at the beginning of each legendre! call to permit a normalization norm to validate that the given maximum $(\ell,m)$ will be within the ability to satisfy. The default case always returns nothing. A custom normalization should throw an error if lmax or mmax is out of bounds or return nothing otherwise.

For example, the precomputed coefficients of LegendreNormCoeff are limited to a given domain at time of construction and cannot be used to calculate terms to arbitrary orders/degrees.