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_prime(int n, T x);
20 template <class T, class ``__Policy``>
21 ``__sf_result`` legendre_p_prime(int n, T x, const ``__Policy``&);
23 template <class T, class ``__Policy``>
24 std::vector<T> legendre_p_zeros(int l, const ``__Policy``&);
27 std::vector<T> legendre_p_zeros(int l);
30 ``__sf_result`` legendre_p(int n, int m, T x);
32 template <class T, class ``__Policy``>
33 ``__sf_result`` legendre_p(int n, int m, T x, const ``__Policy``&);
36 ``__sf_result`` legendre_q(unsigned n, T x);
38 template <class T, class ``__Policy``>
39 ``__sf_result`` legendre_q(unsigned n, T x, const ``__Policy``&);
41 template <class T1, class T2, class T3>
42 ``__sf_result`` legendre_next(unsigned l, T1 x, T2 Pl, T3 Plm1);
44 template <class T1, class T2, class T3>
45 ``__sf_result`` legendre_next(unsigned l, unsigned m, T1 x, T2 Pl, T3 Plm1);
50 The return type of these functions is computed using the __arg_promotion_rules:
51 note than when there is a single template argument the result is the same type
52 as that argument or `double` if the template argument is an integer type.
59 ``__sf_result`` legendre_p(int l, T x);
61 template <class T, class ``__Policy``>
62 ``__sf_result`` legendre_p(int l, T x, const ``__Policy``&);
64 Returns the Legendre Polynomial of the first kind:
68 Requires -1 <= x <= 1, otherwise returns the result of __domain_error.
70 Negative orders are handled via the reflection formula:
72 [:P[sub -l-1](x) = P[sub l](x)]
74 The following graph illustrates the behaviour of the first few
80 ``__sf_result`` legendre_p_prime(int n, T x);
82 template <class T, class ``__Policy``>
83 ``__sf_result`` legendre_p_prime(int n, T x, const ``__Policy``&);
85 Returns the derivatives of the Legendre polynomials.
87 template <class T, class ``__Policy``>
88 std::vector<T> legendre_p_zeros(int l, const ``__Policy``&);
91 std::vector<T> legendre_p_zeros(int l);
93 The zeros of the Legendre polynomials are calculated by Newton's method using an initial guess given by Tricomi with root bracketing provided by Szego.
95 Since the Legendre polynomials are alternatively even and odd, only the non-negative zeros are returned.
96 For the odd Legendre polynomials, the first zero is always zero.
97 The rest of the zeros are returned in increasing order.
99 Note that the argument to the routine is an integer, and the output is a floating-point type.
100 Hence the template argument is mandatory.
101 The time to extract a single root is linear in `l` (this is scaling to evaluate the Legendre polynomials), so recovering all roots is [bigo](`l`[super 2]).
102 Algorithms with linear scaling [@ https://doi.org/10.1137/06067016X exist] for recovering all roots, but requires tooling not currently built into boost.math.
103 This implementation proceeds under the assumption that calculating zeros of these functions will not be a bottleneck for any workflow.
106 ``__sf_result`` legendre_p(int l, int m, T x);
108 template <class T, class ``__Policy``>
109 ``__sf_result`` legendre_p(int l, int m, T x, const ``__Policy``&);
111 Returns the associated Legendre polynomial of the first kind:
113 [equation legendre_1]
115 Requires -1 <= x <= 1, otherwise returns the result of __domain_error.
117 Negative values of /l/ and /m/ are handled via the identity relations:
119 [equation legendre_3]
121 [caution The definition of the associated Legendre polynomial used here
122 includes a leading Condon-Shortley phase term of (-1)[super m]. This
123 matches the definition given by Abramowitz and Stegun (8.6.6) and that
124 used by [@http://mathworld.wolfram.com/LegendrePolynomial.html Mathworld]
125 and [@http://documents.wolfram.com/mathematica/functions/LegendreP
126 Mathematica's LegendreP function]. However, uses in the literature
127 do not always include this phase term, and strangely the specification
128 for the associated Legendre function in the C++ TR1 (assoc_legendre)
129 also omits it, in spite of stating that it uses Abramowitz and Stegun
130 as the final arbiter on these matters.
134 [@http://mathworld.wolfram.com/LegendrePolynomial.html
135 Weisstein, Eric W. "Legendre Polynomial."
136 From MathWorld--A Wolfram Web Resource].
138 Abramowitz, M. and Stegun, I. A. (Eds.). "Legendre Functions" and
139 "Orthogonal Polynomials." Ch. 22 in Chs. 8 and 22 in Handbook of
140 Mathematical Functions with Formulas, Graphs, and Mathematical Tables,
141 9th printing. New York: Dover, pp. 331-339 and 771-802, 1972.
145 ``__sf_result`` legendre_q(unsigned n, T x);
147 template <class T, class ``__Policy``>
148 ``__sf_result`` legendre_q(unsigned n, T x, const ``__Policy``&);
150 Returns the value of the Legendre polynomial that is the second solution
151 to the Legendre differential equation, for example:
153 [equation legendre_2]
155 Requires -1 <= x <= 1, otherwise __domain_error is called.
157 The following graph illustrates the first few Legendre functions of the
162 template <class T1, class T2, class T3>
163 ``__sf_result`` legendre_next(unsigned l, T1 x, T2 Pl, T3 Plm1);
165 Implements the three term recurrence relation for the Legendre
166 polynomials, this function can be used to create a sequence of
167 values evaluated at the same /x/, and for rising /l/. This recurrence
168 relation holds for Legendre Polynomials of both the first and second kinds.
170 [equation legendre_4]
172 For example we could produce a vector of the first 10 polynomial
175 double x = 0.5; // Abscissa value
177 v.push_back(legendre_p(0, x));
178 v.push_back(legendre_p(1, x));
179 for(unsigned l = 1; l < 10; ++l)
180 v.push_back(legendre_next(l, x, v[l], v[l-1]));
181 // Double check values:
182 for(unsigned l = 1; l < 10; ++l)
183 assert(v[l] == legendre_p(l, x));
185 Formally the arguments are:
188 [[l][The degree of the last polynomial calculated.]]
189 [[x][The abscissa value]]
190 [[Pl][The value of the polynomial evaluated at degree /l/.]]
191 [[Plm1][The value of the polynomial evaluated at degree /l-1/.]]
194 template <class T1, class T2, class T3>
195 ``__sf_result`` legendre_next(unsigned l, unsigned m, T1 x, T2 Pl, T3 Plm1);
197 Implements the three term recurrence relation for the Associated Legendre
198 polynomials, this function can be used to create a sequence of
199 values evaluated at the same /x/, and for rising /l/.
201 [equation legendre_5]
203 For example we could produce a vector of the first m+10 polynomial
206 double x = 0.5; // Abscissa value
209 v.push_back(legendre_p(m, m, x));
210 v.push_back(legendre_p(1 + m, m, x));
211 for(unsigned l = 1; l < 10; ++l)
212 v.push_back(legendre_next(l + 10, m, x, v[l], v[l-1]));
213 // Double check values:
214 for(unsigned l = 1; l < 10; ++l)
215 assert(v[l] == legendre_p(10 + l, m, x));
217 Formally the arguments are:
220 [[l][The degree of the last polynomial calculated.]]
221 [[m][The order of the Associated Polynomial.]]
222 [[x][The abscissa value]]
223 [[Pl][The value of the polynomial evaluated at degree /l/.]]
224 [[Plm1][The value of the polynomial evaluated at degree /l-1/.]]
229 The following table shows peak errors (in units of epsilon)
230 for various domains of input arguments.
231 Note that only results for the widest floating point type on the system are
232 given as narrower types have __zero_error.
238 [table_legendre_p_associated_]
240 Note that the worst errors occur when the order increases, values greater than
241 ~120 are very unlikely to produce sensible results, especially in the associated
242 polynomial case when the degree is also large. Further the relative errors
243 are likely to grow arbitrarily large when the function is very close to a root.
247 A mixture of spot tests of values calculated using functions.wolfram.com,
248 and randomly generated test data are
249 used: the test data was computed using
250 [@http://shoup.net/ntl/doc/RR.txt NTL::RR] at 1000-bit precision.
254 These functions are implemented using the stable three term
255 recurrence relations. These relations guarantee low absolute error
256 but cannot guarantee low relative error near one of the roots of the
259 [endsect] [/section:beta_function The Beta Function]
261 Copyright 2006 John Maddock and Paul A. Bristow.
262 Distributed under the Boost Software License, Version 1.0.
263 (See accompanying file LICENSE_1_0.txt or copy at
264 http://www.boost.org/LICENSE_1_0.txt).