Boost C++ Libraries

...one of the most highly regarded and expertly designed C++ library projects in the world. Herb Sutter and Andrei Alexandrescu, C++ Coding Standards

This is the documentation for an old version of Boost. Click here to view this page for the latest version.
PrevUpHomeNext

Spherical Harmonics

Synopsis
#include <boost/math/special_functions/spherical_harmonic.hpp>
namespace boost{ namespace math{

template <class T1, class T2>
std::complex<calculated-result-type> spherical_harmonic(unsigned n, int m, T1 theta, T2 phi);

template <class T1, class T2, class Policy>
std::complex<calculated-result-type> spherical_harmonic(unsigned n, int m, T1 theta, T2 phi, const Policy&);

template <class T1, class T2>
calculated-result-type spherical_harmonic_r(unsigned n, int m, T1 theta, T2 phi);

template <class T1, class T2, class Policy>
calculated-result-type spherical_harmonic_r(unsigned n, int m, T1 theta, T2 phi, const Policy&);

template <class T1, class T2>
calculated-result-type spherical_harmonic_i(unsigned n, int m, T1 theta, T2 phi);

template <class T1, class T2, class Policy>
calculated-result-type spherical_harmonic_i(unsigned n, int m, T1 theta, T2 phi, const Policy&);

}} // namespaces
Description

The return type of these functions is computed using the result type calculation rules when T1 and T2 are different types.

The final Policy argument is optional and can be used to control the behaviour of the function: how it handles errors, what level of precision to use etc. Refer to the policy documentation for more details.

template <class T1, class T2>
std::complex<calculated-result-type> spherical_harmonic(unsigned n, int m, T1 theta, T2 phi);

template <class T1, class T2, class Policy>
std::complex<calculated-result-type> spherical_harmonic(unsigned n, int m, T1 theta, T2 phi, const Policy&);

Returns the value of the Spherical Harmonic Ynm(theta, phi):

The spherical harmonics Ynm(theta, phi) are the angular portion of the solution to Laplace's equation in spherical coordinates where azimuthal symmetry is not present.

[Caution] Caution

Care must be taken in correctly identifying the arguments to this function: θ is taken as the polar (colatitudinal) coordinate with θ in [0, π], and φas the azimuthal (longitudinal) coordinate with φin [0,2π). This is the convention used in Physics, and matches the definition used by Mathematica in the function SpericalHarmonicY, but is opposite to the usual mathematical conventions.

Some other sources include an additional Condon-Shortley phase term of (-1)m in the definition of this function: note however that our definition of the associated Legendre polynomial already includes this term.

This implementation returns zero for m > n

For θ outside [0, π] and φ outside [0, 2π] this implementation follows the convention used by Mathematica: the function is periodic with period π in θ and 2π in φ. Please note that this is not the behaviour one would get from a casual application of the function's definition. Cautious users should keep θ and φ to the range [0, π] and [0, 2π] respectively.

See: Weisstein, Eric W. "Spherical Harmonic." From MathWorld--A Wolfram Web Resource.

template <class T1, class T2>
calculated-result-type spherical_harmonic_r(unsigned n, int m, T1 theta, T2 phi);

template <class T1, class T2, class Policy>
calculated-result-type spherical_harmonic_r(unsigned n, int m, T1 theta, T2 phi, const Policy&);

Returns the real part of Ynm(theta, phi):

template <class T1, class T2>
calculated-result-type spherical_harmonic_i(unsigned n, int m, T1 theta, T2 phi);

template <class T1, class T2, class Policy>
calculated-result-type spherical_harmonic_i(unsigned n, int m, T1 theta, T2 phi, const Policy&);

Returns the imaginary part of Ynm(theta, phi):

Accuracy

The following table shows peak errors for various domains of input arguments. Note that only results for the widest floating point type on the system are given as narrower types have effectively zero error. Peak errors are the same for both the real and imaginary parts, as the error is dominated by calculation of the associated Legendre polynomials: especially near the roots of the associated Legendre function.

All values are in units of epsilon.

Table 8.38. Error rates for spherical_harmonic_r

GNU C++ version 7.1.0
linux
double

GNU C++ version 7.1.0
linux
long double

Sun compiler version 0x5150
Sun Solaris
long double

Microsoft Visual C++ version 14.1
Win32
double

Spherical Harmonics

Max = 1.58ε (Mean = 0.0707ε)

Max = 2.89e+03ε (Mean = 108ε)

Max = 1.03e+04ε (Mean = 327ε)

Max = 2.27e+04ε (Mean = 725ε)


Table 8.39. Error rates for spherical_harmonic_i

GNU C++ version 7.1.0
linux
double

GNU C++ version 7.1.0
linux
long double

Sun compiler version 0x5150
Sun Solaris
long double

Microsoft Visual C++ version 14.1
Win32
double

Spherical Harmonics

Max = 1.36ε (Mean = 0.0765ε)

Max = 2.89e+03ε (Mean = 108ε)

Max = 1.03e+04ε (Mean = 327ε)

Max = 2.27e+04ε (Mean = 725ε)


Note that the worst errors occur when the degree increases, values greater than ~120 are very unlikely to produce sensible results, especially when the order is also large. Further the relative errors are likely to grow arbitrarily large when the function is very close to a root.

Testing

A mixture of spot tests of values calculated using functions.wolfram.com, and randomly generated test data are used: the test data was computed using NTL::RR at 1000-bit precision.

Implementation

These functions are implemented fairly naively using the formulae given above. Some extra care is taken to prevent roundoff error when converting from polar coordinates (so for example the 1-x2 term used by the associated Legendre functions is calculated without roundoff error using x = cos(theta), and 1-x2 = sin2(theta)). The limiting factor in the error rates for these functions is the need to calculate values near the roots of the associated Legendre functions.


PrevUpHomeNext