Public Documentation

Legendre.LegendreModule

Collections of functions which compute the associated Legendre functions.

Based on implementation described in Limpanuparb and Milthorpe (2014) “Associated Legendre Polynomials and Spherical Harmonics Computation for Chemistry Applications” arXiv:1410.1748v1

source
Legendre.AbstractLegendreNormType
abstract type AbstractLegendreNorm end

Abstract trait supertype for normalization conditions of the Associated Legendre polynomials.

source
Legendre.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]
source
Legendre.LegendreSphereCoeffType
LegendreSphereCoeff{T}

Table type of precomputed recursion relation coefficients for the spherical harmonic normalization. Alias for LegendreNormCoeff{LegendreSphereNorm,T}.

source
Legendre.LegendreSphereNormType
struct LegendreSphereNorm <: AbstractLegendreNorm end

Trait type denoting the spherical-harmonic normalization of the associated Legendre polynomials.

source
Legendre.LegendreUnitCoeffType
LegendreUnitCoeff{T}

Precomputed recursion relation coefficients for the standard unit normalization. Alias for LegendreNormCoeff{LegendreUnitNorm,T}.

source
Legendre.LegendreUnitNormType
struct LegendreUnitNorm <: AbstractLegendreNorm end

Trait type denoting the unit normalization of the associated Legendre polynomials.

source
Legendre.NlmMethod
N = Nlm([T=Float64], l, m)

Computes the normalization constant

\[ N_ℓ^m ≡ \sqrt{\frac{2ℓ+1}{4π} \frac{(ℓ-m)!}{(ℓ+m)!}}\]

which defines the Spherical Harmonic normalized functions $λ_ℓ^m(x)$ in terms of the standard unit normalized $P_ℓ^m(x)$

\[ λ_ℓ^m(x) ≡ N_ℓ^m P_ℓ^m(x)\]

using numbers of type T.

See also Plm and λlm.

source
Legendre.Pl!Method
Pl!(P, l::Integer, x)

Fills the array P with the unit-normalized Legendre polynomial values $P_ℓ(x)$ for fixed order $m = 0$; equivalent to legendre!(LegendreUnitNorm(), P, l, 0, x).

source
Legendre.PlMethod
p = Pl(l::Integer, x::Number)

Computes the Legendre polynomials using unit normalization and for degree $m = 0$; equivalent to p = legendre(LegendreUnitNorm(), l, 0, x).

source
Legendre.Plm!Method
Plm!(P, l::Integer, m::Integer, x)

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

source
Legendre.PlmMethod
p = Plm(l::Integer, m::Integer, x::Number)

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

source
Legendre.legendre!Method
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 l and order 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).
source
Legendre.legendreMethod
p = legendre(norm::AbstractLegendreNorm, l::Integer, m::Integer, x::Number)
P = legendre.(norm::AbstractLegendreNorm, l, m, x)

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

With broadcasting syntax, the polynomials can be computed over any iterable x. Furthermore,

  • If l isa Integer && m isa Integer, then the output P has the same shape as x and is filled with the polynomial values of order l and degree m.
  • If l isa UnitRange && m isa Integer, then l is interpreted as lmax, and the output P has one more dimension than x with the trailing dimension spanning the degrees 0 ≤ l ≤ lmax.
  • If l isa UnitRange && m isa UnitRange, then l is interpreted as lmax and m as mmax, and the output P has two more dimensions than x with the trailing dimensions spanning the degrees 0 ≤ l ≤ lmax and orders 0 ≤ m ≤ mmax, respectively.

Note that in second and third case, the UnitRanges must satisify first(l) == 0 and first(m) == 0.

source
Legendre.legendreMethod
p = legendre(norm::AbstractLegendreNorm, l::Integer, x::Number)
P = legendre.(norm::AbstractLegendreNorm, l, x)

Computes the associated Legendre polynomial assuming the order $m = 0$; equivalent to legendre(norm, l, 0, x) and legendre.(norm, l, 0, x).

source
Legendre.λlm!Method
λlm!(Λ, l::Integer, m::Integer, x)

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

source
Legendre.λlmMethod
λ = λlm(l::Integer, m::Integer, x::Number)

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

source