Public Documentation
Legendre.Legendre — ModuleCollections 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
Legendre.AbstractLegendreNorm — Typeabstract type AbstractLegendreNorm endAbstract trait supertype for normalization conditions of the Associated Legendre polynomials.
Legendre.LegendreNormCoeff — Typestruct LegendreNormCoeff{N<:AbstractLegendreNorm,T<:Real} <: AbstractLegendreNormPrecomputed 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]Legendre.LegendreSphereCoeff — TypeLegendreSphereCoeff{T}Table type of precomputed recursion relation coefficients for the spherical harmonic normalization. Alias for LegendreNormCoeff{LegendreSphereNorm,T}.
Legendre.LegendreSphereNorm — Typestruct LegendreSphereNorm <: AbstractLegendreNorm endTrait type denoting the spherical-harmonic normalization of the associated Legendre polynomials.
Legendre.LegendreUnitCoeff — TypeLegendreUnitCoeff{T}Precomputed recursion relation coefficients for the standard unit normalization. Alias for LegendreNormCoeff{LegendreUnitNorm,T}.
Legendre.LegendreUnitNorm — Typestruct LegendreUnitNorm <: AbstractLegendreNorm endTrait type denoting the unit normalization of the associated Legendre polynomials.
Legendre.Nlm — MethodN = Nlm([T=Float64], l, m)Computes the normalization constant
which defines the Spherical Harmonic normalized functions $λ_ℓ^m(x)$ in terms of the standard unit normalized $P_ℓ^m(x)$
using numbers of type T.
Legendre.Pl! — MethodPl!(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).
Legendre.Pl — Methodp = 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).
Legendre.Plm! — MethodPlm!(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).
Legendre.Plm — Methodp = Plm(l::Integer, m::Integer, x::Number)Computes the associated Legendre polynomials using unit normalization; equivalent to p = legendre(LegendreUnitNorm(), l, m, x).
Legendre.legendre! — Methodlegendre!(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 atxfor degreeland orderm. - If
ndims(Λ) == ndims(x) + 1, thenlis interpreted aslmax, andΛfilled with polynomial values for all degrees0 ≤ l ≤ lmaxof orderm. - If
ndims(Λ) == ndims(x) + 2, thenlis interpreted aslmaxandmasmmax, andΛis filled with polynomial values for all degrees0 ≤ l ≤ lmaxand orders0 ≤ m ≤ min(mmax, l).
Legendre.legendre — Methodp = 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 outputPhas the same shape asxand is filled with the polynomial values of orderland degreem. - If
l isa UnitRange && m isa Integer, thenlis interpreted aslmax, and the outputPhas one more dimension thanxwith the trailing dimension spanning the degrees0 ≤ l ≤ lmax. - If
l isa UnitRange && m isa UnitRange, thenlis interpreted aslmaxandmasmmax, and the outputPhas two more dimensions thanxwith the trailing dimensions spanning the degrees0 ≤ l ≤ lmaxand orders0 ≤ m ≤ mmax, respectively.
Note that in second and third case, the UnitRanges must satisify first(l) == 0 and first(m) == 0.
Legendre.legendre — Methodp = 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).
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).
Legendre.λlm — Methodλ = λlm(l::Integer, m::Integer, x::Number)Computes the associated Legendre polynomials using spherical-harmonic normalization; equivalent to λ = legendre(LegendreSphereNorm(), l, m, x).