legendre.qbk 8.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265
  1. [section:legendre Legendre (and Associated) Polynomials]
  2. [h4 Synopsis]
  3. ``
  4. #include <boost/math/special_functions/legendre.hpp>
  5. ``
  6. namespace boost{ namespace math{
  7. template <class T>
  8. ``__sf_result`` legendre_p(int n, T x);
  9. template <class T, class ``__Policy``>
  10. ``__sf_result`` legendre_p(int n, T x, const ``__Policy``&);
  11. template <class T>
  12. ``__sf_result`` legendre_p_prime(int n, T x);
  13. template <class T, class ``__Policy``>
  14. ``__sf_result`` legendre_p_prime(int n, T x, const ``__Policy``&);
  15. template <class T, class ``__Policy``>
  16. std::vector<T> legendre_p_zeros(int l, const ``__Policy``&);
  17. template <class T>
  18. std::vector<T> legendre_p_zeros(int l);
  19. template <class T>
  20. ``__sf_result`` legendre_p(int n, int m, T x);
  21. template <class T, class ``__Policy``>
  22. ``__sf_result`` legendre_p(int n, int m, T x, const ``__Policy``&);
  23. template <class T>
  24. ``__sf_result`` legendre_q(unsigned n, T x);
  25. template <class T, class ``__Policy``>
  26. ``__sf_result`` legendre_q(unsigned n, T x, const ``__Policy``&);
  27. template <class T1, class T2, class T3>
  28. ``__sf_result`` legendre_next(unsigned l, T1 x, T2 Pl, T3 Plm1);
  29. template <class T1, class T2, class T3>
  30. ``__sf_result`` legendre_next(unsigned l, unsigned m, T1 x, T2 Pl, T3 Plm1);
  31. }} // namespaces
  32. The return type of these functions is computed using the __arg_promotion_rules:
  33. note than when there is a single template argument the result is the same type
  34. as that argument or `double` if the template argument is an integer type.
  35. [optional_policy]
  36. [h4 Description]
  37. template <class T>
  38. ``__sf_result`` legendre_p(int l, T x);
  39. template <class T, class ``__Policy``>
  40. ``__sf_result`` legendre_p(int l, T x, const ``__Policy``&);
  41. Returns the Legendre Polynomial of the first kind:
  42. [equation legendre_0]
  43. Requires -1 <= x <= 1, otherwise returns the result of __domain_error.
  44. Negative orders are handled via the reflection formula:
  45. [:P[sub -l-1](x) = P[sub l](x)]
  46. The following graph illustrates the behaviour of the first few
  47. Legendre Polynomials:
  48. [graph legendre_p]
  49. template <class T>
  50. ``__sf_result`` legendre_p_prime(int n, T x);
  51. template <class T, class ``__Policy``>
  52. ``__sf_result`` legendre_p_prime(int n, T x, const ``__Policy``&);
  53. Returns the derivatives of the Legendre polynomials.
  54. template <class T, class ``__Policy``>
  55. std::vector<T> legendre_p_zeros(int l, const ``__Policy``&);
  56. template <class T>
  57. std::vector<T> legendre_p_zeros(int l);
  58. 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.
  59. Since the Legendre polynomials are alternatively even and odd, only the non-negative zeros are returned.
  60. For the odd Legendre polynomials, the first zero is always zero.
  61. The rest of the zeros are returned in increasing order.
  62. Note that the argument to the routine is an integer, and the output is a floating-point type.
  63. Hence the template argument is mandatory.
  64. 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]).
  65. Algorithms with linear scaling [@ https://doi.org/10.1137/06067016X exist] for recovering all roots, but requires tooling not currently built into boost.math.
  66. This implementation proceeds under the assumption that calculating zeros of these functions will not be a bottleneck for any workflow.
  67. template <class T>
  68. ``__sf_result`` legendre_p(int l, int m, T x);
  69. template <class T, class ``__Policy``>
  70. ``__sf_result`` legendre_p(int l, int m, T x, const ``__Policy``&);
  71. Returns the associated Legendre polynomial of the first kind:
  72. [equation legendre_1]
  73. Requires -1 <= x <= 1, otherwise returns the result of __domain_error.
  74. Negative values of /l/ and /m/ are handled via the identity relations:
  75. [equation legendre_3]
  76. [caution The definition of the associated Legendre polynomial used here
  77. includes a leading Condon-Shortley phase term of (-1)[super m]. This
  78. matches the definition given by Abramowitz and Stegun (8.6.6) and that
  79. used by [@http://mathworld.wolfram.com/LegendrePolynomial.html Mathworld]
  80. and [@http://documents.wolfram.com/mathematica/functions/LegendreP
  81. Mathematica's LegendreP function]. However, uses in the literature
  82. do not always include this phase term, and strangely the specification
  83. for the associated Legendre function in the C++ TR1 (assoc_legendre)
  84. also omits it, in spite of stating that it uses Abramowitz and Stegun
  85. as the final arbiter on these matters.
  86. See:
  87. [@http://mathworld.wolfram.com/LegendrePolynomial.html
  88. Weisstein, Eric W. "Legendre Polynomial."
  89. From MathWorld--A Wolfram Web Resource].
  90. Abramowitz, M. and Stegun, I. A. (Eds.). "Legendre Functions" and
  91. "Orthogonal Polynomials." Ch. 22 in Chs. 8 and 22 in Handbook of
  92. Mathematical Functions with Formulas, Graphs, and Mathematical Tables,
  93. 9th printing. New York: Dover, pp. 331-339 and 771-802, 1972.
  94. ]
  95. template <class T>
  96. ``__sf_result`` legendre_q(unsigned n, T x);
  97. template <class T, class ``__Policy``>
  98. ``__sf_result`` legendre_q(unsigned n, T x, const ``__Policy``&);
  99. Returns the value of the Legendre polynomial that is the second solution
  100. to the Legendre differential equation, for example:
  101. [equation legendre_2]
  102. Requires -1 <= x <= 1, otherwise __domain_error is called.
  103. The following graph illustrates the first few Legendre functions of the
  104. second kind:
  105. [graph legendre_q]
  106. template <class T1, class T2, class T3>
  107. ``__sf_result`` legendre_next(unsigned l, T1 x, T2 Pl, T3 Plm1);
  108. Implements the three term recurrence relation for the Legendre
  109. polynomials, this function can be used to create a sequence of
  110. values evaluated at the same /x/, and for rising /l/. This recurrence
  111. relation holds for Legendre Polynomials of both the first and second kinds.
  112. [equation legendre_4]
  113. For example we could produce a vector of the first 10 polynomial
  114. values using:
  115. double x = 0.5; // Abscissa value
  116. vector<double> v;
  117. v.push_back(legendre_p(0, x));
  118. v.push_back(legendre_p(1, x));
  119. for(unsigned l = 1; l < 10; ++l)
  120. v.push_back(legendre_next(l, x, v[l], v[l-1]));
  121. // Double check values:
  122. for(unsigned l = 1; l < 10; ++l)
  123. assert(v[l] == legendre_p(l, x));
  124. Formally the arguments are:
  125. [variablelist
  126. [[l][The degree of the last polynomial calculated.]]
  127. [[x][The abscissa value]]
  128. [[Pl][The value of the polynomial evaluated at degree /l/.]]
  129. [[Plm1][The value of the polynomial evaluated at degree /l-1/.]]
  130. ]
  131. template <class T1, class T2, class T3>
  132. ``__sf_result`` legendre_next(unsigned l, unsigned m, T1 x, T2 Pl, T3 Plm1);
  133. Implements the three term recurrence relation for the Associated Legendre
  134. polynomials, this function can be used to create a sequence of
  135. values evaluated at the same /x/, and for rising /l/.
  136. [equation legendre_5]
  137. For example we could produce a vector of the first m+10 polynomial
  138. values using:
  139. double x = 0.5; // Abscissa value
  140. int m = 10; // order
  141. vector<double> v;
  142. v.push_back(legendre_p(m, m, x));
  143. v.push_back(legendre_p(1 + m, m, x));
  144. for(unsigned l = 1; l < 10; ++l)
  145. v.push_back(legendre_next(l + 10, m, x, v[l], v[l-1]));
  146. // Double check values:
  147. for(unsigned l = 1; l < 10; ++l)
  148. assert(v[l] == legendre_p(10 + l, m, x));
  149. Formally the arguments are:
  150. [variablelist
  151. [[l][The degree of the last polynomial calculated.]]
  152. [[m][The order of the Associated Polynomial.]]
  153. [[x][The abscissa value]]
  154. [[Pl][The value of the polynomial evaluated at degree /l/.]]
  155. [[Plm1][The value of the polynomial evaluated at degree /l-1/.]]
  156. ]
  157. [h4 Accuracy]
  158. The following table shows peak errors (in units of epsilon)
  159. for various domains of input arguments.
  160. Note that only results for the widest floating point type on the system are
  161. given as narrower types have __zero_error.
  162. [table_legendre_p]
  163. [table_legendre_q]
  164. [table_legendre_p_associated_]
  165. Note that the worst errors occur when the order increases, values greater than
  166. ~120 are very unlikely to produce sensible results, especially in the associated
  167. polynomial case when the degree is also large. Further the relative errors
  168. are likely to grow arbitrarily large when the function is very close to a root.
  169. [h4 Testing]
  170. A mixture of spot tests of values calculated using functions.wolfram.com,
  171. and randomly generated test data are
  172. used: the test data was computed using
  173. [@http://shoup.net/ntl/doc/RR.txt NTL::RR] at 1000-bit precision.
  174. [h4 Implementation]
  175. These functions are implemented using the stable three term
  176. recurrence relations. These relations guarantee low absolute error
  177. but cannot guarantee low relative error near one of the roots of the
  178. polynomials.
  179. [endsect] [/section:beta_function The Beta Function]
  180. [/
  181. Copyright 2006 John Maddock and Paul A. Bristow.
  182. Distributed under the Boost Software License, Version 1.0.
  183. (See accompanying file LICENSE_1_0.txt or copy at
  184. http://www.boost.org/LICENSE_1_0.txt).
  185. ]