SphericalHarmonics.jl
SphericalHarmonics.NorthPole
— TypeSphericalHarmonics.NorthPole <: SphericalHarmonics.Pole
The angle $\theta = 0$ in spherical polar coordinates.
SphericalHarmonics.Pole
— TypeSphericalHarmonics.SouthPole
— TypeSphericalHarmonics.SouthPole <: SphericalHarmonics.Pole
The angle $\theta = π$ in spherical polar coordinates.
SphericalHarmonics.SphericalHarmonicsCache
— TypeSphericalHarmonicsCache
Preallocate arrays of associated Legendre polynomials and spherical harmonics. Such an object may be allocated using cache
.
SphericalHarmonics.allocate_p
— MethodSphericalHarmonics.allocate_p([T::Type = Float64], lmax::Integer)
Allocate an array large enough to store an entire set of Associated Legendre Polynomials $\bar{P}_ℓ^m(x)$ of maximum degree $ℓ$.
SphericalHarmonics.allocate_y
— FunctionSphericalHarmonics.allocate_y([T::Type = ComplexF64], lmax::Integer)
Allocate an array large enough to store an entire set of spherical harmonics $Y_{ℓ,m}(θ,ϕ)$ of maximum degree $ℓ$.
SphericalHarmonics.associatedLegendre
— MethodSphericalHarmonics.associatedLegendre(θ; l::Integer, m::Integer, [coeff = nothing], [norm = SphericalHarmonics.LMnorm()])
Evaluate the normalized associated Legendre polynomial $\bar{P}_ℓ^m(\cos \theta)$. Optionally a matrix of coefficients returned by compute_coefficients
may be provided.
The keyword argument norm
may be used to specify the how the polynomials are normalized. See computePlmcostheta
for the possible normalization options.
SphericalHarmonics.associatedLegendrex
— MethodSphericalHarmonics.associatedLegendrex(x; l::Integer, m::Integer, [coeff = nothing], [norm = SphericalHarmonics.LMnorm()])
Evaluate the normalized associated Legendre polynomial $\bar{P}_ℓ^m(x)$. Optionally a matrix of coefficients returned by compute_coefficients
may be provided.
The keyword argument norm
may be used to specify the how the polynomials are normalized. See computePlmcostheta
for the possible normalization options.
SphericalHarmonics.cache
— Methodcache([T::Type = Float64], lmax, [m_range = SphericalHarmonicModes.FullRange], [SHType = SphericalHarmonics.ComplexHarmonics()])
Allocate arrays to evaluate associated Legendre polynomials and spherical harmonics. The returned object may be passed to computePlmcostheta!
and computeYlm!
. The coefficients are cached and need not be recomputed.
Examples
julia> S = SphericalHarmonics.cache(1);
julia> computePlmcostheta!(S, pi/3, 1)
3-element normalized AssociatedLegendrePolynomials{Float64} for lmax = 1 and cosθ = 0.5:
0.3989422804014327
0.34549414947133555
-0.4231421876608172
julia> computeYlm!(S, pi/3, pi/4, 1)
4-element SHArray(::Vector{ComplexF64}, (ML(0:1, -1:1),)):
0.2820947917738782 + 0.0im
0.21157109383040865 - 0.2115710938304086im
0.24430125595146002 + 0.0im
-0.21157109383040865 - 0.2115710938304086im
Choosing a new lmax
in computePlmcostheta!
expands the cache if necessary.
julia> computePlmcostheta!(S, pi/3, 2)
6-element normalized AssociatedLegendrePolynomials{Float64} for lmax = 2 and cosθ = 0.5:
0.3989422804014327
0.34549414947133555
-0.4231421876608172
-0.11150775725954817
-0.4730873478787801
0.40970566147202964
SphericalHarmonics.computePlmcostheta!
— FunctioncomputePlmcostheta!(S::SphericalHarmonicsCache, θ, [lmax::Integer]; [norm = SphericalHarmonics.LMnorm()])
Compute an entire set of normalized Associated Legendre Polynomials $\bar{P}_ℓ^m(\cos θ)$ using the pre-computed coefficients in S
, and store the result in S
. If lmax
is not provided, the value of lmax
for which coefficients have been computed in S
is used.
The keyword argument norm
may be used to specify the how the polynomials are normalized. See computePlmcostheta
for the possible normalization options.
SphericalHarmonics.computePlmcostheta!
— MethodcomputePlmcostheta!(P::AbstractVector, θ, lmax::Integer, coeff; [norm = SphericalHarmonics.LMnorm()])
computePlmcostheta!(P::AbstractVector, θ::SphericalHarmonics.Pole, lmax::Integer; [norm = SphericalHarmonics.LMnorm()])
Compute an entire set of normalized Associated Legendre Polynomials $\bar{P}_ℓ^m(\cos θ)$ using the given coefficients, and store in the array P
. The matrix coeff
may be computed using compute_coefficients
.
The keyword argument norm
may be used to specify the how the polynomials are normalized. See computePlmcostheta
for the possible normalization options.
computePlmcostheta!(P::AbstractVector, θ, lmax::Integer, m::Integer, coeff)
Compute the Associated Legendre Polynomials $\bar{P}_ℓ^m(\cos θ)$ for for the specified $m$ and all $ℓ$ lying in $|m| ≤ ℓ ≤ ℓ_\mathrm{max}$. The array P
needs to be large enough to hold all the polynomials for $0 ≤ ℓ ≤ ℓ_\mathrm{max}$ and $0 ≤ m ≤ ℓ$, as the recursive evaluation requires the computation of extra elements. Pre-existing values in P
may be overwritten, even for azimuthal orders not equal to $m$.
SphericalHarmonics.computePlmcostheta
— MethodcomputePlmcostheta(θ; lmax::Integer, [m::Integer], [norm = SphericalHarmonics.LMnorm()])
computePlmcostheta(θ, lmax::Integer[, m::Integer]; [norm = SphericalHarmonics.LMnorm()])
Compute an entire set of normalized Associated Legendre Polynomials $\bar{P}_ℓ^m(\cos θ)$ where $0 ≤ ℓ ≤ ℓ_\mathrm{max}$ and $0 ≤ m ≤ ℓ$ for colatitude $\theta$. If m
is provided then only the polynomials corresponding to the specified m
are computed.
The polynomials are normalized as
\[\bar{P}_{\ell}^m = \sqrt{\frac{(2\ell + 1)(\ell-m)!}{2\pi (\ell+m)!}} P_{\ell m},\]
where $P_{\ell m}$ are the standard associated Legendre polynomials, and are defined in terms of Legendre polynomials $P_\ell(x)$ as
\[P_{\ell m}\left(x\right)=(-1)^m \left(1-x^{2}\right)^{m/2}\frac{d^{m}}{dx^{m}}P_{\ell}\left(x\right).\]
The normalized polynomials $\bar{P}_{\ell}^m$ satisfy
\[\int_{-1}^{1} dx\,\left| \bar{P}_{\ell}^m(x) \right|^2 = \frac{1}{\pi}\]
A different normalization may be chosen by specifying the keyword argument norm
, with optionally the Condon-Shortley phase disabled by passing the keyword argument csphase
to the constructor of the normalization specifier. The possible normalization options are:
SphericalHarmonics.LMnorm([; csphase = true])
: the default normalization described aboveSphericalHarmonics.Orthonormal([; csphase = true])
: Orthonormal polynomials that are defined as
\[\tilde{P}_{\ell}^m = \sqrt{\frac{(2\ell + 1)(\ell-m)!}{2(\ell+m)!}} P_{\ell m} = \sqrt{\pi} \bar{P}_{\ell m},\]
and satisfy
\[\int_{-1}^{1} \tilde{P}_{\ell m}(x) \tilde{P}_{k m}(x) dx = \delta_{ℓk}\]
SphericalHarmonics.Unnormalized([; csphase = true])
: The polynomials $P_{ℓm}$ that satisfy $P_{ℓm}(1)=\delta_{m0}$
within numerical precision bounds, as well as
\[\int_{-1}^{1} P_{\ell m}(x) P_{k m}(x) dx = \frac{2(\ell+m)!}{(2\ell+1)(\ell-m)!}\delta_{ℓk}\]
In each case setting csphase = false
will lead to an overall factor of $(-1)^m$ being multiplied to the polynomials.
Returns an AbstractVector
that may be indexed using (ℓ,m)
pairs aside from the canonical indexing with Int
s.
The precision of the result may be increased by using arbitrary-precision arguments.
Examples
julia> P = computePlmcostheta(pi/2, lmax = 1)
3-element normalized AssociatedLegendrePolynomials{Float64} for lmax = 1 and cosθ = 6.123e-17:
0.3989422804014327
4.231083042742082e-17
-0.4886025119029199
julia> P[(0,0)]
0.3989422804014327
julia> P = computePlmcostheta(big(pi)/2, lmax = 1)
3-element normalized AssociatedLegendrePolynomials{BigFloat} for lmax = 1 and cosθ = 5.485e-78:
0.3989422804014326779399460599343818684758586311649346576659258296706579258993008
3.789785583114350800838137317730900078444216599640987847808409161681770236721676e-78
-0.4886025119029199215863846228383470045758856081942277021382431574458410003616367
SphericalHarmonics.computePlmx!
— FunctioncomputePlmx!(S::SphericalHarmonicsCache, x[, lmax::Integer]; [norm = SphericalHarmonics.LMnorm()])
Compute an entire set of normalized Associated Legendre Polynomials $\bar{P}_ℓ^m(x)$ using the pre-computed coefficients in S
, and store the result in S
. If lmax
is not provided, the value of lmax
for which coefficients have been computed in S
is used.
The keyword argument norm
may be used to specify the how the polynomials are normalized. See computePlmcostheta
for the possible normalization options.
SphericalHarmonics.computePlmx!
— MethodcomputePlmx!(P::AbstractVector, x, lmax::Integer, coeff::AbstractMatrix; [norm = SphericalHarmonics.LMnorm()])
Compute an entire set of normalized Associated Legendre Polynomials $\bar{P}_ℓ^m(x)$ using the given coefficients, and store in the array P
. The matrix coeff
may be computed using compute_coefficients
.
The argument x
needs to lie in $-1 ≤ x ≤ 1$. The function implicitly assumes that $x = \cos\theta$ where $0 ≤ \theta ≤ π$.
The keyword argument norm
may be used to specify the how the polynomials are normalized. See computePlmcostheta
for the possible normalization options.
computePlmx!(P::AbstractVector, x, lmax::Integer, m::Integer, coeff::AbstractMatrix)
Compute the set of normalized Associated Legendre Polynomials $\bar{P}_ℓ^m(x)$ for for the specified $m$ and all $ℓ$ lying in $|m| ≤ ℓ ≤ ℓ_\mathrm{max}$ .
SphericalHarmonics.computePlmx
— MethodcomputePlmx(x; lmax::Integer, [m::Integer], [norm = SphericalHarmonics.LMnorm()])
computePlmx(x, lmax::Integer[, m::Integer]; [norm = SphericalHarmonics.LMnorm()])
Compute an entire set of normalized Associated Legendre Polynomials $\bar{P}_ℓ^m(x)$ where $0 ≤ ℓ ≤ ℓ_\mathrm{max}$ and $0 ≤ m ≤ ℓ$. If m
is provided then only the polynomials for that azimuthal order are computed.
The argument x
needs to lie in $-1 ≤ x ≤ 1$. The function implicitly assumes that $x = \cos\theta$ where $0 ≤ \theta ≤ π$.
The keyword argument norm
may be used to specify the how the polynomials are normalized. See computePlmcostheta
for the possible normalization options.
SphericalHarmonics.computeYlm
— FunctioncomputeYlm(θ, ϕ; lmax::Integer, [m::Integer], [m_range = SphericalHarmonics.FullRange], [SHType = SphericalHarmonics.ComplexHarmonics()])
computeYlm(θ::SphericalHarmonics.Pole; lmax::Integer, [m_range = SphericalHarmonics.FullRange], [SHType = SphericalHarmonics.ComplexHarmonics()])
Compute an entire set of spherical harmonics $Y_{ℓ,m}(θ,ϕ)$ for $0 ≤ ℓ ≤ ℓ_\mathrm{max}$ and $-ℓ ≤ m ≤ ℓ$, for colatitude $\theta$ and azimuth $\phi$. If $m$ is provided then only the polynomials for the specified $m$ are computed.
Returns an AbstractVector
that may be indexed using (l,m)
pairs aside from the canonical indexing with Int
s.
The optional argument m_range
decides if the spherical harmonics for negative m
values are computed. By default the functions for all values of m
are evaluated. Setting m_range
to SphericalHarmonics.ZeroTo
would result in only functions for m ≥ 0
being evaluated.
The optional argument SHType
may be used to choose between real and complex harmonics. To compute real spherical harmonics, set this to SphericalHarmonics.RealHarmonics()
. The complex harmonics are defined as
\[Y_{\ell,m}\left(\theta,\phi\right)=\begin{cases} \frac{1}{\sqrt{2}}\bar{P}_{\ell}^{m}\left(\cos\theta\right)\exp\left(im\phi\right), & m\geq0,\\ \left(-1\right)^{m}Y_{\ell,-m}^{*}\left(\theta,\phi\right), & m<0. \end{cases}\]
This definition corresponds to Laplace spherical harmonics, whereas the quantum mechanical definition prepends a Condon-Shortley phase on the harmonics.
The real spherical harmonics are defined as
\[Y_{\ell,m}\left(\theta,\phi\right)=\begin{cases} \bar{P}_{\ell}^{\left|m\right|}\left(\cos\theta\right)\sin\left|m\right|\phi, & m<0,\\ \bar{P}_{\ell}^{0}\left(\cos\theta\right)/\sqrt{2}, & m=0,\\ \bar{P}_{\ell}^{m}\left(\cos\theta\right)\cos m\phi, & m>0. \end{cases}\]
The precision of the result may be increased by using arbitrary-precision arguments.
Examples
julia> Y = computeYlm(pi/2, 0, lmax = 1)
4-element SHArray(::Vector{ComplexF64}, (ML(0:1, -1:1),)):
0.2820947917738782 + 0.0im
0.3454941494713355 - 0.0im
2.9918275112863375e-17 + 0.0im
-0.3454941494713355 - 0.0im
julia> Y[(1,-1)] # index using (l,m)
0.3454941494713355 - 0.0im
julia> Y = computeYlm(big(pi)/2, big(0), lmax = big(1)) # Arbitrary precision
4-element SHArray(::Vector{Complex{BigFloat}}, (ML(0:1, -1:1),)):
0.2820947917738781434740397257803862929220253146644994284220428608553212342207478 + 0.0im
0.3454941494713354792652446460318896831393773703262433134867073548945156550201567 - 0.0im
2.679783085063171668225419916118067917387251852939708540164955895366691604430101e-78 + 0.0im
-0.3454941494713354792652446460318896831393773703262433134867073548945156550201567 - 0.0im
julia> computeYlm(SphericalHarmonics.NorthPole(), 0, lmax = 1)
4-element SHArray(::Vector{ComplexF64}, (ML(0:1, -1:1),)):
0.2820947917738782 + 0.0im
-0.0 + 0.0im
0.48860251190292 + 0.0im
0.0 + 0.0im
julia> Y = computeYlm(pi/2, pi/3, lmax = 1, m_range = SphericalHarmonics.ZeroTo, SHType = SphericalHarmonics.RealHarmonics())
3-element SHArray(::Vector{Float64}, (ML(0:1, 0:1),)):
0.2820947917738782
2.9918275112863375e-17
-0.24430125595146002
SphericalHarmonics.computeYlm!
— FunctioncomputeYlm!(S::SphericalHarmonicsCache, θ, ϕ, [lmax::Integer])
Compute an entire set of spherical harmonics $Y_{ℓ,m}(θ,ϕ)$ for $0 ≤ ℓ ≤ ℓ_\mathrm{max}$ using the pre-computed associated Legendre polynomials saved in S
, and store the result in S
. If lmax
is not provided, the value of lmax
for which associated Legendre polynomials have been computed in S
is used.
This function assumes that the associated Legendre polynomials have been pre-computed, and does not perform any check on their values. In general computeYlm!(S::SphericalHarmonicsCache, θ, ϕ, lmax)
should only be called after a preceeding call to computePlmcostheta!(S, θ, lmax)
in order to obtain meaningful results.
SphericalHarmonics.computeYlm!
— FunctioncomputeYlm!(Y::AbstractVector, P::AbstractVector, θ, ϕ; lmax::Integer, [m::Integer], [m_range = SphericalHarmonics.FullRange], [SHType = SphericalHarmonics.ComplexHarmonics()])
computeYlm!(Y::AbstractVector, P::AbstractVector, θ, ϕ, lmax::Integer, [m::Integer, [m_range = SphericalHarmonics.FullRange, [SHType = SphericalHarmonics.ComplexHarmonics()]]])
Compute an entire set of spherical harmonics $Y_{ℓ,m}(θ,ϕ)$ using the precomputed associated Legendre Polynomials $\bar{P}_ℓ^m(\cos θ)$, and store in the array Y
. The array P
may be computed using computePlmcostheta
.
The optional argument m_range
decides if the spherical harmonics for negative m
values are computed. By default the functions for all values of m
are evaluated. Setting m_range
to SphericalHarmonics.ZeroTo
would result in only functions for m ≥ 0
being evaluated. Providing m
would override this, in which case only the polynomials corresponding to the azimuthal order m
would be computed.
The optional argument SHType
may be used to choose between real and complex harmonics. To compute real spherical harmonics, set this to RealHarmonics()
.
This function assumes that the associated Legendre Polynomials have been pre-computed, and does not perform any check on the values of P
.
SphericalHarmonics.computeYlm!
— FunctioncomputeYlm!(Y::AbstractVector, θ, ϕ; lmax::Integer, [m::Integer] [m_range = SphericalHarmonics.FullRange], [SHType = SphericalHarmonics.ComplexHarmonics()])
Compute an entire set of spherical harmonics $Y_{ℓ,m}(θ,ϕ)$ for $0 ≤ ℓ ≤ ℓ_\mathrm{max}$, and store them in the array Y
.
The optional argument m_range
decides if the spherical harmonics for negative m
values are computed. By default the functions for all values of m
are evaluated. Setting m_range
to SphericalHarmonics.ZeroTo
would result in only functions for m ≥ 0
being evaluated. Providing m
would override this, in which case only the polynomials corresponding to the azimuthal order m
would be computed.
The optional argument SHType
may be used to choose between real and complex harmonics. To compute real spherical harmonics, set this to SphericalHarmonics.RealHarmonics()
.
SphericalHarmonics.compute_coefficients
— MethodSphericalHarmonics.compute_coefficients(lmax)
Precompute coefficients $a_ℓ^m$ and $b_ℓ^m$ for all $2 ≤ ℓ ≤ ℓ_\mathrm{max}$ and $0 ≤ m ≤ ℓ-2$.
SphericalHarmonics.compute_coefficients(lmax, m)
Precompute coefficients $a_ℓ^m$ and $b_ℓ^m$ for all $|m| + 2 ≤ ℓ ≤ ℓ_\mathrm{max}$ and the specified $m$.
SphericalHarmonics.sphericalharmonic
— MethodSphericalHarmonics.sphericalharmonic(θ, ϕ; l, m, [SHType = ComplexHarmonics()], [coeff = nothing])
Evaluate the spherical harmonic $Y_{ℓ,m}(θ, ϕ)$. The flag SHType
sets the type of the harmonic computed, and setting this to RealHarmonics()
would evaluate real spherical harmonics. Optionally a precomputed matrix of coefficients returned by compute_coefficients
may be provided.
Example
julia> SphericalHarmonics.sphericalharmonic(π/3, π/3, l = 500, m = 250)
-0.18910100312194328 - 0.32753254516944075im
julia> SphericalHarmonics.sphericalharmonic(π/3, π/3, l = 500, m = 250, SHType = SphericalHarmonics.RealHarmonics())
-0.26742920327340913