Definition

Legendre polynomials

Legendre polynomials satisfy the differential equation

\[\frac{d}{dx}\left(\left(1-x^{2}\right)\frac{d}{dx}\right)P_{\ell}\left(x\right)=-\ell\left(\ell+1\right)P_{\ell}\left(x\right)\]

By default this package evaluates the standard polynomials that satisfy $P_\ell(1) = 1$ and $P_0(x) = 1$. These are normalized as

\[\int_{-1}^1 P_m(x) P_n(x) dx = \frac{2}{2n+1} \delta_{mn}.\]

Associated Legendre polynomials

Associated Legendre polynomials satisfy the differential equation

\[\left[\frac{d}{dx}\left(\left(1-x^{2}\right)\frac{d}{dx}\right)-\frac{m^{2}}{1-x^{2}}\right]P_{\ell}^{m}\left(x\right)=-\ell\left(\ell+1\right)P_{\ell}^{m}\left(x\right),\]

and are related to the Legendre polynomials through

\[P_{\ell}^{m}\left(x\right)=\left(-1\right)^{m}\left(1-x^{2}\right)^{m/2}\frac{d^{m}}{dx^{m}}\left(P_\ell\left(x\right)\right).\]

For $m=0$, we obtain

\[P_{\ell}^{0}\left(x\right)=P_\ell\left(x\right).\]

The Condon-Shortley phase is included by default, and may be toggled by using the keyword argument csphase.

Quick Start

Evaluate the Legendre polynomial for one l at an argumentx as Pl(x, l):

julia> p = Pl(0.5, 3)
-0.4375

julia> p ≈ -7/16 # analytical value
true

Evaluate the associated Legendre Polynomial one l,m pair as Plm(x, l, m):

julia> p = Plm(0.5, 3, 2)
5.624999999999997

julia> p ≈ 45/8 # analytical value
true

Evaluate the nth derivative of Pl for one l as dnPl(x, l, n):

julia> dnPl(0.5, 3, 2)
7.5

Evaluate the Legendre polynomials for l in lmin:lmax as collectPl(x; lmin, lmax). By default lmin is chosen to be 0, and may be omitted.

julia> collectPl(0.5, lmax = 3)
4-element OffsetArray(::Vector{Float64}, 0:3) with eltype Float64 with indices 0:3:
  1.0
  0.5
 -0.125
 -0.4375

Evaluate the associated Legendre Polynomials for order m and degree l in lmin:lmax as collectPlm(x; m, lmin, lmax). By default lmin is chosen to be abs(m), and may be omitted.

julia> collectPlm(0.5, lmax = 5, m = 3)
3-element OffsetArray(::Vector{Float64}, 3:5) with eltype Float64 with indices 3:5:
  -9.742785792574933
 -34.09975027401223
 -42.62468784251535

Evaluate all the nth derivatives of Pl as collectdnPl(x; lmax, n):

julia> collectdnPl(0.5, lmax = 5, n = 3)
6-element OffsetArray(::Vector{Float64}, 0:5) with eltype Float64 with indices 0:5:
  0.0
  0.0
  0.0
 15.0
 52.5
 65.625

Recursion relations

Compute Legendre polynomials, Associated Legendre polynomials, and their derivatives using a 3-term recursion relation (Bonnet’s recursion formula).

\[P_\ell(x) = \left((2\ell-1) x P_{\ell-1}(x) - (\ell-1)P_{\ell - 2}(x)\right)/\ell\]

Analogous to Legendre polynomials, one may evaluate associated Legendre polynomials using a 3-term recursion relation. This is evaluated by iterating over the normalized associated Legendre functions, and multiplying the norm at the final stage. Such an iteration avoids floating-point overflow.

The relation used in evaluating normalized associated Legendre polynomials $\bar{P}_{\ell}^{m}\left(x\right)$ is

\[\bar{P}_{\ell}^{m}\left(x\right)=\alpha_{\ell m}\left(x\bar{P}_{\ell-1}^{m}\left(x\right)-\frac{1}{\alpha_{\ell-1m}}\bar{P}_{\ell-2}^{m}\left(x\right)\right),\]

where

\[\alpha_{\ell m}=\sqrt{\frac{\left(2\ell+1\right)\left(2\ell-1\right)}{\left(\ell-m\right)\left(\ell+m\right)}},\]

The starting value for a specific $m$ is obtained by setting $\ell = m$, in which case we use the explicit relation

\[\bar{P}_{m}^{m}\left(x\right)=\left(-1\right)^{m}\left(2m-1\right)!!\sqrt{\frac{\left(2m+1\right)}{2\left(2m\right)!}}\left(1-x^{2}\right)^{\left(m/2\right)}.\]

The polynomials, thus defined, include the Condon-Shortley phase $(-1)^m$, and may be evaluated using the standard normalization as

\[P_{\ell}^{m}\left(x\right)=\sqrt{\frac{2\left(\ell+m\right)!}{\left(2\ell+1\right)\left(\ell-m\right)!}}\bar{P}_{\ell}^{m}\left(x\right).\]

Main functions

There are six main functions:

Normalization options

By default, the Legendre and associated Legendre polynomials are not normalized. One may specify the normalization through the keyword argument norm. The normalization options accepted are

  • norm = Val(:standard): standard, unnormalized polynomials. This is the default option. Choosing this option will lead to the standard Legendre or associated Legendre polynomials that satisfy $P_\ell(0)=1$ and $P_{\ell0}(0)=1$
  • norm = Val(:normalized): fully normalized polynomials, with an L2 norm of 1. These are defined as

\[ \bar{P}_{\ell m}\left(x\right)=\sqrt{\frac{2\ell+1}{2}\frac{\left(\ell-m\right)!}{\left(\ell+m\right)!}}P_{\ell m}\left(x\right)\]

  • norm = Val(:schmidtquasi): Schmidt quasi-normalized polynomials, also known as Schmidt semi-normalized polynomials. These have an L2 norm of $\sqrt{2(2-\delta_{m0})/(2\ell+1)}$. These are defined as

\[ \bar{P}_{\ell m}\left(x\right)=\sqrt{\left(2-\delta_{m0}\right)\frac{\left(\ell-m\right)!}{\left(\ell+m\right)!}}P_{\ell m}\left(x\right)\]

  • norm = Val(:schmidt): Schmidt normalized polynomials. These have an L2 norm of $\sqrt{2(2-\delta_{m0})}$. These are defined as

\[ \bar{P}_{\ell m}\left(x\right)=\sqrt{\left(2-\delta_{m0}\right)\left(2\ell+1\right)\frac{\left(\ell-m\right)!}{\left(\ell+m\right)!}}P_{\ell m}\left(x\right)\]

Note

Irrespective of the norm specified, the 3-term recursion relations used are stable ones, so polynomials of high degrees may be computed without overflow.

Reference

LegendrePolynomials.PlMethod
Pl(x, l::Integer; [norm = Val(:standard)])

Compute the Legendre Polynomial $P_\ell(x)$ for the argument x and the degree l.

The default norm is chosen to be Val(:standard), in which case the polynomials satisfy Pl(1, l) == 1 for all l. Optionally, the norm may be set to Val(:normalized) to obtain normalized Legendre polynomials. These have an L2 norm of 1.

Examples

julia> Pl(1, 2)
1.0

julia> Pl(0.5, 4) ≈ -37/128 # analytically obtained value
true

julia> Pl(0.5, 20, norm = Val(:normalized))
-0.21895188261094017
source
LegendrePolynomials.PlmMethod
Plm(x, l::Integer, m::Integer; [kwargs...])

Compute the associated Legendre polynomial $P_\ell^m(x)$ for degree l and order $m$ lying in -l:l.

For m == 0 this function returns Legendre polynomials.

Optional keyword arguments

  • norm: The normalization used in the function. Possible options are Val(:standard) (default), Val(:normalized), Val(:schmidtquasi) or Val(:schmidt).
    • The functions obtained with norm = Val(:normalized) have an L2 norm of $1$.
    • The functions obtained with norm = Val(:schmidtquasi) have an L2 norm of $\sqrt{\frac{2(2-\delta_{m0})}{2l+1}}$.
    • The functions obtained with norm = Val(:schmidt) have an L2 norm of $\sqrt{2(2-\delta_{m0})}$.
  • csphase::Bool: The Condon-shortley phase $(-1)^m$, which is included by default.

Examples

julia> Plm(0.5, 3, 2) ≈ 45/8 # analytically obtained value
true

julia> Plm(0.5, 4, 0) == Pl(0.5, 4)
true
source
LegendrePolynomials.collectPl!Method
collectPl!(v::AbstractVector, x; [lmin::Integer = 0], [lmax::Integer = length(v) - 1 + lmin], [kwargs...])

Compute the Legendre Polynomials $P_\ell(x)$ for the argument x and degrees l = lmin:lmax, and store them in v.

At output v[firstindex(v) + l] == Pl(x, l; kwargs...).

Optional keyword arguments

  • norm: The normalization used in the function. Possible options are Val(:standard) (default) and Val(:normalized). The functions obtained with norm = Val(:normalized) have an L2 norm of $1$.

Examples

julia> v = zeros(4);

julia> collectPl!(v, 0.5);

julia> all(zip(0:3, v)) do (l, vl); vl ≈ Pl(0.5, l); end
true

julia> collectPl!(v, 0.5, norm = Val(:normalized));

julia> all(zip(0:3, v)) do (l,vl); vl ≈ Pl(0.5, l, norm = Val(:normalized)); end
true

julia> v = zeros(0:4);

julia> collectPl!(v, 0.5, lmax = 3) # only l from 0 to 3 are populated
5-element OffsetArray(::Vector{Float64}, 0:4) with eltype Float64 with indices 0:4:
  1.0
  0.5
 -0.125
 -0.4375
  0.0
source
LegendrePolynomials.collectPlMethod
collectPl(x; lmax::Integer, [lmin::Integer = 0], [kwargs...])

Compute the Legendre Polynomial $P_\ell(x)$ for the argument x and degrees l = lmin:lmax. Return P with indices lmin:lmax, where P[l] == Pl(x, l; kwargs...)

Optional keyword arguments

  • norm: The normalization used in the function. Possible options are Val(:standard) (default) and Val(:normalized). The functions obtained with norm = Val(:normalized) have an L2 norm of $1$.

Examples

julia> v = collectPl(0.5, lmax = 4);

julia> all(l -> v[l] ≈ Pl(0.5, l), 0:4)
true

julia> v = collectPl(0.5, lmax = 4, norm = Val(:normalized));

julia> all(l -> v[l] ≈ Pl(0.5, l, norm = Val(:normalized)), 0:4)
true
source
LegendrePolynomials.collectPlm!Method
collectPlm!(v::AbstractVector, x; m::Integer, [lmin::Integer = abs(m)], [lmax::Integer = length(v) - 1 + lmin], [kwargs...])

Compute the associated Legendre polynomials $P_\ell^m(x)$ for the argument x, degrees l = lmin:lmax and order m lying in -l:l, and store the result in v.

At output, v[l + firstindex(v)] == Plm(x, l, m; kwargs...) for l = lmin:lmax.

Optional keyword arguments

  • norm: The normalization used in the function. Possible options are Val(:standard) (default), Val(:normalized), Val(:schmidtquasi) or Val(:schmidt).
    • The functions obtained with norm = Val(:normalized) have an L2 norm of $1$.
    • The functions obtained with norm = Val(:schmidtquasi) have an L2 norm of $\sqrt{\frac{2(2-\delta_{m0})}{2l+1}}$.
    • The functions obtained with norm = Val(:schmidt) have an L2 norm of $\sqrt{2(2-\delta_{m0})}$.
  • csphase::Bool: The Condon-shortley phase $(-1)^m$, which is included by default.

Examples

julia> v = zeros(2);

julia> collectPlm!(v, 0.5, lmax = 3, m = 2);

julia> all(zip(2:3, v)) do (l, vl); vl ≈ Plm(0.5, l, 2); end
true
source
LegendrePolynomials.collectPlmMethod
collectPlm(x; m::Integer, lmax::Integer, [lmin::Integer = abs(m)], [kwargs...])

Compute the associated Legendre polynomials $P_\ell^m(x)$ for the argument x, degrees l = lmin:lmax and order m lying in -l:l.

Returns v with indices lmin:lmax, with v[l] == Plm(x, l, m; kwargs...).

Optional keyword arguments

  • norm: The normalization used in the function. Possible options are Val(:standard) (default), Val(:normalized), Val(:schmidtquasi) or Val(:schmidt).
    • The functions obtained with norm = Val(:normalized) have an L2 norm of $1$.
    • The functions obtained with norm = Val(:schmidtquasi) have an L2 norm of $\sqrt{\frac{2(2-\delta_{m0})}{2l+1}}$.
    • The functions obtained with norm = Val(:schmidt) have an L2 norm of $\sqrt{2(2-\delta_{m0})}$.
  • csphase::Bool: The Condon-shortley phase $(-1)^m$, which is included by default.
  • Tnorm: Type of the normalization factor, which is dynamicall inferred by default, and is used in allocating an appropriate container. In case it is known that the norm may be expressed as a certain type without overflow (eg. Float64), this may be provided. Care must be taken to avoid overflows if Tnorm is passed as an argument.

Examples

julia> v = collectPlm(0.5, lmax = 3, m = 2);

julia> all(l -> v[l] ≈ Plm(0.5, l, 2), 2:3)
true

julia> v = collectPlm(0.5, lmax = 3, m = 2, norm = Val(:normalized));

julia> all(l -> v[l] ≈ Plm(0.5, l, 2, norm = Val(:normalized)), 2:3)
true
source
LegendrePolynomials.collectdnPl!Method
collectdnPl!(v::AbstractVector, x; lmax::Integer, n::Integer)

Compute the $n$-th derivative of a Legendre Polynomial $P_\ell(x)$ for the argument x and all degrees l = 0:lmax, and store the result in v.

The order of the derivative n must be greater than or equal to zero.

At output, v[l + firstindex(v)] == dnPl(x, l, n) for l = 0:lmax.

Examples

julia> v = zeros(4);

julia> collectdnPl!(v, 0.5, lmax = 3, n = 2)
4-element Vector{Float64}:
 0.0
 0.0
 3.0
 7.5
source
LegendrePolynomials.collectdnPlMethod
collectdnPl(x; lmax::Integer, n::Integer)

Compute the $n$-th derivative of a Legendre Polynomial $P_\ell(x)$ for the argument x and all degrees l = 0:lmax.

The order of the derivative n must be greater than or equal to zero.

Returns v with indices 0:lmax, where v[l] == dnPl(x, l, n).

Examples

julia> collectdnPl(0.5, lmax = 3, n = 2)
4-element OffsetArray(::Vector{Float64}, 0:3) with eltype Float64 with indices 0:3:
 0.0
 0.0
 3.0
 7.5
source
LegendrePolynomials.dnPlFunction
dnPl(x, l::Integer, n::Integer, [cache::AbstractVector])

Compute the $n$-th derivative $d^n P_\ell(x)/dx^n$ of the Legendre polynomial $P_\ell(x)$. Optionally a pre-allocated vector cache may be provided, which must have a minimum length of l - n + 1 and may be overwritten during the computation.

The order of the derivative n must be non-negative. For n == 0 this function just returns Legendre polynomials.

Examples

julia> dnPl(0.5, 3, 2) # second derivative of P3(x) at x = 0.5
7.5

julia> dnPl(0.5, 4, 0) == Pl(0.5, 4) # zero-th order derivative == Pl(x)
true
source