1 [section:legendre Legendre (and Associated) Polynomials]
6 #include <boost/math/special_functions/legendre.hpp>
9 namespace boost{ namespace math{
12 ``__sf_result`` legendre_p(int n, T x);
14 template <class T, class ``__Policy``>
15 ``__sf_result`` legendre_p(int n, T x, const ``__Policy``&);
18 ``__sf_result`` legendre_p(int n, int m, T x);
20 template <class T, class ``__Policy``>
21 ``__sf_result`` legendre_p(int n, int m, T x, const ``__Policy``&);
24 ``__sf_result`` legendre_q(unsigned n, T x);
26 template <class T, class ``__Policy``>
27 ``__sf_result`` legendre_q(unsigned n, T x, const ``__Policy``&);
29 template <class T1, class T2, class T3>
30 ``__sf_result`` legendre_next(unsigned l, T1 x, T2 Pl, T3 Plm1);
32 template <class T1, class T2, class T3>
33 ``__sf_result`` legendre_next(unsigned l, unsigned m, T1 x, T2 Pl, T3 Plm1);
38 The return type of these functions is computed using the __arg_promotion_rules:
39 note than when there is a single template argument the result is the same type
40 as that argument or `double` if the template argument is an integer type.
47 ``__sf_result`` legendre_p(int l, T x);
49 template <class T, class ``__Policy``>
50 ``__sf_result`` legendre_p(int l, T x, const ``__Policy``&);
52 Returns the Legendre Polynomial of the first kind:
56 Requires -1 <= x <= 1, otherwise returns the result of __domain_error.
58 Negative orders are handled via the reflection formula:
60 P[sub -l-1](x) = P[sub l](x)
62 The following graph illustrates the behaviour of the first few
68 ``__sf_result`` legendre_p(int l, int m, T x);
70 template <class T, class ``__Policy``>
71 ``__sf_result`` legendre_p(int l, int m, T x, const ``__Policy``&);
73 Returns the associated Legendre polynomial of the first kind:
77 Requires -1 <= x <= 1, otherwise returns the result of __domain_error.
79 Negative values of /l/ and /m/ are handled via the identity relations:
83 [caution The definition of the associated Legendre polynomial used here
84 includes a leading Condon-Shortley phase term of (-1)[super m]. This
85 matches the definition given by Abramowitz and Stegun (8.6.6) and that
86 used by [@http://mathworld.wolfram.com/LegendrePolynomial.html Mathworld]
87 and [@http://documents.wolfram.com/mathematica/functions/LegendreP
88 Mathematica's LegendreP function]. However, uses in the literature
89 do not always include this phase term, and strangely the specification
90 for the associated Legendre function in the C++ TR1 (assoc_legendre)
91 also omits it, in spite of stating that it uses Abramowitz and Stegun
92 as the final arbiter on these matters.
96 [@http://mathworld.wolfram.com/LegendrePolynomial.html
97 Weisstein, Eric W. "Legendre Polynomial."
98 From MathWorld--A Wolfram Web Resource].
100 Abramowitz, M. and Stegun, I. A. (Eds.). "Legendre Functions" and
101 "Orthogonal Polynomials." Ch. 22 in Chs. 8 and 22 in Handbook of
102 Mathematical Functions with Formulas, Graphs, and Mathematical Tables,
103 9th printing. New York: Dover, pp. 331-339 and 771-802, 1972.
107 ``__sf_result`` legendre_q(unsigned n, T x);
109 template <class T, class ``__Policy``>
110 ``__sf_result`` legendre_q(unsigned n, T x, const ``__Policy``&);
112 Returns the value of the Legendre polynomial that is the second solution
113 to the Legendre differential equation, for example:
115 [equation legendre_2]
117 Requires -1 <= x <= 1, otherwise __domain_error is called.
119 The following graph illustrates the first few Legendre functions of the
124 template <class T1, class T2, class T3>
125 ``__sf_result`` legendre_next(unsigned l, T1 x, T2 Pl, T3 Plm1);
127 Implements the three term recurrence relation for the Legendre
128 polynomials, this function can be used to create a sequence of
129 values evaluated at the same /x/, and for rising /l/. This recurrence
130 relation holds for Legendre Polynomials of both the first and second kinds.
132 [equation legendre_4]
134 For example we could produce a vector of the first 10 polynomial
137 double x = 0.5; // Abscissa value
139 v.push_back(legendre_p(0, x));
140 v.push_back(legendre_p(1, x));
141 for(unsigned l = 1; l < 10; ++l)
142 v.push_back(legendre_next(l, x, v[l], v[l-1]));
143 // Double check values:
144 for(unsigned l = 1; l < 10; ++l)
145 assert(v[l] == legendre_p(l, x));
147 Formally the arguments are:
150 [[l][The degree of the last polynomial calculated.]]
151 [[x][The abscissa value]]
152 [[Pl][The value of the polynomial evaluated at degree /l/.]]
153 [[Plm1][The value of the polynomial evaluated at degree /l-1/.]]
156 template <class T1, class T2, class T3>
157 ``__sf_result`` legendre_next(unsigned l, unsigned m, T1 x, T2 Pl, T3 Plm1);
159 Implements the three term recurrence relation for the Associated Legendre
160 polynomials, this function can be used to create a sequence of
161 values evaluated at the same /x/, and for rising /l/.
163 [equation legendre_5]
165 For example we could produce a vector of the first m+10 polynomial
168 double x = 0.5; // Abscissa value
171 v.push_back(legendre_p(m, m, x));
172 v.push_back(legendre_p(1 + m, m, x));
173 for(unsigned l = 1; l < 10; ++l)
174 v.push_back(legendre_next(l + 10, m, x, v[l], v[l-1]));
175 // Double check values:
176 for(unsigned l = 1; l < 10; ++l)
177 assert(v[l] == legendre_p(10 + l, m, x));
179 Formally the arguments are:
182 [[l][The degree of the last polynomial calculated.]]
183 [[m][The order of the Associated Polynomial.]]
184 [[x][The abscissa value]]
185 [[Pl][The value of the polynomial evaluated at degree /l/.]]
186 [[Plm1][The value of the polynomial evaluated at degree /l-1/.]]
191 The following table shows peak errors (in units of epsilon)
192 for various domains of input arguments.
193 Note that only results for the widest floating point type on the system are
194 given as narrower types have __zero_error.
200 [table_legendre_p_associated_]
202 Note that the worst errors occur when the order increases, values greater than
203 ~120 are very unlikely to produce sensible results, especially in the associated
204 polynomial case when the degree is also large. Further the relative errors
205 are likely to grow arbitrarily large when the function is very close to a root.
209 A mixture of spot tests of values calculated using functions.wolfram.com,
210 and randomly generated test data are
211 used: the test data was computed using
212 [@http://shoup.net/ntl/doc/RR.txt NTL::RR] at 1000-bit precision.
216 These functions are implemented using the stable three term
217 recurrence relations. These relations guarantee low absolute error
218 but cannot guarantee low relative error near one of the roots of the
221 [endsect][/section:beta_function The Beta Function]
223 Copyright 2006 John Maddock and Paul A. Bristow.
224 Distributed under the Boost Software License, Version 1.0.
225 (See accompanying file LICENSE_1_0.txt or copy at
226 http://www.boost.org/LICENSE_1_0.txt).