Next: , Previous: Legendre Polynomials, Up: Legendre Functions and Spherical Harmonics


7.24.2 Associated Legendre Polynomials and Spherical Harmonics

The following functions compute the associated Legendre Polynomials P_l^m(x). Note that this function grows combinatorially with l and can overflow for l larger than about 150. There is no trouble for small m, but overflow occurs when m and l are both large. Rather than allow overflows, these functions refuse to calculate P_l^m(x) and return GSL_EOVRFLW when they can sense that l and m are too big.

If you want to calculate a spherical harmonic, then do not use these functions. Instead use gsl_sf_legendre_sphPlm below, which uses a similar recursion, but with the normalized functions.

— Function: double gsl_sf_legendre_Plm (int l, int m, double x)
— Function: int gsl_sf_legendre_Plm_e (int l, int m, double x, gsl_sf_result * result)

These routines compute the associated Legendre polynomial P_l^m(x) for m >= 0, l >= m, |x| <= 1.

— Function: int gsl_sf_legendre_Plm_array (int lmax, int m, double x, double result_array[])
— Function: int gsl_sf_legendre_Plm_deriv_array (int lmax, int m, double x, double result_array[], double result_deriv_array[])

These functions compute arrays of Legendre polynomials P_l^m(x) and derivatives dP_l^m(x)/dx, for m >= 0, l = |m|, ..., lmax, |x| <= 1.

— Function: double gsl_sf_legendre_sphPlm (int l, int m, double x)
— Function: int gsl_sf_legendre_sphPlm_e (int l, int m, double x, gsl_sf_result * result)

These routines compute the normalized associated Legendre polynomial \sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(x) suitable for use in spherical harmonics. The parameters must satisfy m >= 0, l >= m, |x| <= 1. Theses routines avoid the overflows that occur for the standard normalization of P_l^m(x).

— Function: int gsl_sf_legendre_sphPlm_array (int lmax, int m, double x, double result_array[])
— Function: int gsl_sf_legendre_sphPlm_deriv_array (int lmax, int m, double x, double result_array[], double result_deriv_array[])

These functions compute arrays of normalized associated Legendre functions \sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(x), and derivatives, for m >= 0, l = |m|, ..., lmax, |x| <= 1.0

— Function: int gsl_sf_legendre_array_size (const int lmax, const int m)

This function returns the size of result_array[] needed for the array versions of P_l^m(x), lmax - m + 1. An inline version of this function is used when HAVE_INLINE is defined.