1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122
|
<html lang="en">
<head>
<title>Associated Legendre Polynomials and Spherical Harmonics - GNU Scientific Library -- Reference Manual</title>
<meta http-equiv="Content-Type" content="text/html">
<meta name="description" content="GNU Scientific Library -- Reference Manual">
<meta name="generator" content="makeinfo 4.11">
<link title="Top" rel="start" href="index.html#Top">
<link rel="up" href="Legendre-Functions-and-Spherical-Harmonics.html" title="Legendre Functions and Spherical Harmonics">
<link rel="prev" href="Legendre-Polynomials.html" title="Legendre Polynomials">
<link rel="next" href="Conical-Functions.html" title="Conical Functions">
<link href="http://www.gnu.org/software/texinfo/" rel="generator-home" title="Texinfo Homepage">
<!--
Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010 The GSL Team.
Permission is granted to copy, distribute and/or modify this document
under the terms of the GNU Free Documentation License, Version 1.3 or
any later version published by the Free Software Foundation; with the
Invariant Sections being ``GNU General Public License'' and ``Free Software
Needs Free Documentation'', the Front-Cover text being ``A GNU Manual'',
and with the Back-Cover Text being (a) (see below). A copy of the
license is included in the section entitled ``GNU Free Documentation
License''.
(a) The Back-Cover Text is: ``You have the freedom to copy and modify this
GNU Manual.''-->
<meta http-equiv="Content-Style-Type" content="text/css">
<style type="text/css"><!--
pre.display { font-family:inherit }
pre.format { font-family:inherit }
pre.smalldisplay { font-family:inherit; font-size:smaller }
pre.smallformat { font-family:inherit; font-size:smaller }
pre.smallexample { font-size:smaller }
pre.smalllisp { font-size:smaller }
span.sc { font-variant:small-caps }
span.roman { font-family:serif; font-weight:normal; }
span.sansserif { font-family:sans-serif; font-weight:normal; }
--></style>
</head>
<body>
<div class="node">
<p>
<a name="Associated-Legendre-Polynomials-and-Spherical-Harmonics"></a>
Next: <a rel="next" accesskey="n" href="Conical-Functions.html">Conical Functions</a>,
Previous: <a rel="previous" accesskey="p" href="Legendre-Polynomials.html">Legendre Polynomials</a>,
Up: <a rel="up" accesskey="u" href="Legendre-Functions-and-Spherical-Harmonics.html">Legendre Functions and Spherical Harmonics</a>
<hr>
</div>
<h4 class="subsection">7.24.2 Associated Legendre Polynomials and Spherical Harmonics</h4>
<p>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 <code>GSL_EOVRFLW</code> when
they can sense that l and m are too big.
<p>If you want to calculate a spherical harmonic, then <em>do not</em> use
these functions. Instead use <code>gsl_sf_legendre_sphPlm</code> below,
which uses a similar recursion, but with the normalized functions.
<div class="defun">
— Function: double <b>gsl_sf_legendre_Plm</b> (<var>int l, int m, double x</var>)<var><a name="index-gsl_005fsf_005flegendre_005fPlm-719"></a></var><br>
— Function: int <b>gsl_sf_legendre_Plm_e</b> (<var>int l, int m, double x, gsl_sf_result * result</var>)<var><a name="index-gsl_005fsf_005flegendre_005fPlm_005fe-720"></a></var><br>
<blockquote><p>These routines compute the associated Legendre polynomial
P_l^m(x) for <!-- {$m \ge 0$} -->
m >= 0, <!-- {$l \ge m$} -->
l >= m, <!-- {$|x| \le 1$} -->
|x| <= 1.
<!-- Exceptional Return Values: GSL_EDOM, GSL_EOVRFLW -->
</p></blockquote></div>
<div class="defun">
— Function: int <b>gsl_sf_legendre_Plm_array</b> (<var>int lmax, int m, double x, double result_array</var>[])<var><a name="index-gsl_005fsf_005flegendre_005fPlm_005farray-721"></a></var><br>
— Function: int <b>gsl_sf_legendre_Plm_deriv_array</b> (<var>int lmax, int m, double x, double result_array</var>[]<var>, double result_deriv_array</var>[])<var><a name="index-gsl_005fsf_005flegendre_005fPlm_005fderiv_005farray-722"></a></var><br>
<blockquote><p>These functions compute arrays of Legendre polynomials
P_l^m(x) and derivatives dP_l^m(x)/dx,
for <!-- {$m \ge 0$} -->
m >= 0, <!-- {$l = |m|, \dots, lmax$} -->
l = |m|, ..., lmax, <!-- {$|x| \le 1$} -->
|x| <= 1.
<!-- Exceptional Return Values: GSL_EDOM, GSL_EOVRFLW -->
</p></blockquote></div>
<div class="defun">
— Function: double <b>gsl_sf_legendre_sphPlm</b> (<var>int l, int m, double x</var>)<var><a name="index-gsl_005fsf_005flegendre_005fsphPlm-723"></a></var><br>
— Function: int <b>gsl_sf_legendre_sphPlm_e</b> (<var>int l, int m, double x, gsl_sf_result * result</var>)<var><a name="index-gsl_005fsf_005flegendre_005fsphPlm_005fe-724"></a></var><br>
<blockquote><p>These routines compute the normalized associated Legendre polynomial
<!-- {$\sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(x)$} -->
\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 \ge 0$} -->
m >= 0, <!-- {$l \ge m$} -->
l >= m, <!-- {$|x| \le 1$} -->
|x| <= 1. Theses routines avoid the overflows
that occur for the standard normalization of P_l^m(x).
<!-- Exceptional Return Values: GSL_EDOM -->
</p></blockquote></div>
<div class="defun">
— Function: int <b>gsl_sf_legendre_sphPlm_array</b> (<var>int lmax, int m, double x, double result_array</var>[])<var><a name="index-gsl_005fsf_005flegendre_005fsphPlm_005farray-725"></a></var><br>
— Function: int <b>gsl_sf_legendre_sphPlm_deriv_array</b> (<var>int lmax, int m, double x, double result_array</var>[]<var>, double result_deriv_array</var>[])<var><a name="index-gsl_005fsf_005flegendre_005fsphPlm_005fderiv_005farray-726"></a></var><br>
<blockquote><p>These functions compute arrays of normalized associated Legendre functions
<!-- {$\sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(x)$} -->
\sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(x),
and derivatives,
for <!-- {$m \ge 0$} -->
m >= 0, <!-- {$l = |m|, \dots, lmax$} -->
l = |m|, ..., lmax, <!-- {$|x| \le 1$} -->
|x| <= 1.0
<!-- Exceptional Return Values: GSL_EDOM -->
</p></blockquote></div>
<div class="defun">
— Function: int <b>gsl_sf_legendre_array_size</b> (<var>const int lmax, const int m</var>)<var><a name="index-gsl_005fsf_005flegendre_005farray_005fsize-727"></a></var><br>
<blockquote><p>This function returns the size of <var>result_array</var>[] needed for the array
versions of P_l^m(x), <var>lmax</var> - <var>m</var> + 1. An inline version of this function is used when <code>HAVE_INLINE</code> is defined.
<!-- Exceptional Return Values: none -->
</p></blockquote></div>
<hr>The GNU Scientific Library - a free numerical library licensed under the GNU GPL<br>Back to the <a href="/software/gsl/">GNU Scientific Library Homepage</a></body></html>
|