xref: /netbsd-src/external/gpl3/gcc.old/dist/libstdc++-v3/include/tr1/legendre_function.tcc (revision 8feb0f0b7eaff0608f8350bbfa3098827b4bb91b)
136ac495dSmrg // Special functions -*- C++ -*-
236ac495dSmrg 
3*8feb0f0bSmrg // Copyright (C) 2006-2020 Free Software Foundation, Inc.
436ac495dSmrg //
536ac495dSmrg // This file is part of the GNU ISO C++ Library.  This library is free
636ac495dSmrg // software; you can redistribute it and/or modify it under the
736ac495dSmrg // terms of the GNU General Public License as published by the
836ac495dSmrg // Free Software Foundation; either version 3, or (at your option)
936ac495dSmrg // any later version.
1036ac495dSmrg //
1136ac495dSmrg // This library is distributed in the hope that it will be useful,
1236ac495dSmrg // but WITHOUT ANY WARRANTY; without even the implied warranty of
1336ac495dSmrg // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
1436ac495dSmrg // GNU General Public License for more details.
1536ac495dSmrg //
1636ac495dSmrg // Under Section 7 of GPL version 3, you are granted additional
1736ac495dSmrg // permissions described in the GCC Runtime Library Exception, version
1836ac495dSmrg // 3.1, as published by the Free Software Foundation.
1936ac495dSmrg 
2036ac495dSmrg // You should have received a copy of the GNU General Public License and
2136ac495dSmrg // a copy of the GCC Runtime Library Exception along with this program;
2236ac495dSmrg // see the files COPYING3 and COPYING.RUNTIME respectively.  If not, see
2336ac495dSmrg // <http://www.gnu.org/licenses/>.
2436ac495dSmrg 
2536ac495dSmrg /** @file tr1/legendre_function.tcc
2636ac495dSmrg  *  This is an internal header file, included by other library headers.
2736ac495dSmrg  *  Do not attempt to use it directly. @headername{tr1/cmath}
2836ac495dSmrg  */
2936ac495dSmrg 
3036ac495dSmrg //
3136ac495dSmrg // ISO C++ 14882 TR1: 5.2  Special functions
3236ac495dSmrg //
3336ac495dSmrg 
3436ac495dSmrg // Written by Edward Smith-Rowland based on:
3536ac495dSmrg //   (1) Handbook of Mathematical Functions,
3636ac495dSmrg //       ed. Milton Abramowitz and Irene A. Stegun,
3736ac495dSmrg //       Dover Publications,
3836ac495dSmrg //       Section 8, pp. 331-341
3936ac495dSmrg //   (2) The Gnu Scientific Library, http://www.gnu.org/software/gsl
4036ac495dSmrg //   (3) Numerical Recipes in C, by W. H. Press, S. A. Teukolsky,
4136ac495dSmrg //       W. T. Vetterling, B. P. Flannery, Cambridge University Press (1992),
4236ac495dSmrg //       2nd ed, pp. 252-254
4336ac495dSmrg 
4436ac495dSmrg #ifndef _GLIBCXX_TR1_LEGENDRE_FUNCTION_TCC
4536ac495dSmrg #define _GLIBCXX_TR1_LEGENDRE_FUNCTION_TCC 1
4636ac495dSmrg 
47c0a68be4Smrg #include <tr1/special_function_util.h>
4836ac495dSmrg 
4936ac495dSmrg namespace std _GLIBCXX_VISIBILITY(default)
5036ac495dSmrg {
51a2dc1f3fSmrg _GLIBCXX_BEGIN_NAMESPACE_VERSION
52a2dc1f3fSmrg 
5336ac495dSmrg #if _GLIBCXX_USE_STD_SPEC_FUNCS
5436ac495dSmrg # define _GLIBCXX_MATH_NS ::std
5536ac495dSmrg #elif defined(_GLIBCXX_TR1_CMATH)
5636ac495dSmrg namespace tr1
5736ac495dSmrg {
5836ac495dSmrg # define _GLIBCXX_MATH_NS ::std::tr1
5936ac495dSmrg #else
6036ac495dSmrg # error do not include this header directly, use <cmath> or <tr1/cmath>
6136ac495dSmrg #endif
6236ac495dSmrg   // [5.2] Special functions
6336ac495dSmrg 
6436ac495dSmrg   // Implementation-space details.
6536ac495dSmrg   namespace __detail
6636ac495dSmrg   {
6736ac495dSmrg     /**
68c0a68be4Smrg      *   @brief  Return the Legendre polynomial by recursion on degree
6936ac495dSmrg      *           @f$ l @f$.
7036ac495dSmrg      *
7136ac495dSmrg      *   The Legendre function of @f$ l @f$ and @f$ x @f$,
7236ac495dSmrg      *   @f$ P_l(x) @f$, is defined by:
7336ac495dSmrg      *   @f[
7436ac495dSmrg      *     P_l(x) = \frac{1}{2^l l!}\frac{d^l}{dx^l}(x^2 - 1)^{l}
7536ac495dSmrg      *   @f]
7636ac495dSmrg      *
77c0a68be4Smrg      *   @param  l  The degree of the Legendre polynomial.  @f$l >= 0@f$.
7836ac495dSmrg      *   @param  x  The argument of the Legendre polynomial.  @f$|x| <= 1@f$.
7936ac495dSmrg      */
8036ac495dSmrg     template<typename _Tp>
8136ac495dSmrg     _Tp
__poly_legendre_p(unsigned int __l,_Tp __x)8236ac495dSmrg     __poly_legendre_p(unsigned int __l, _Tp __x)
8336ac495dSmrg     {
8436ac495dSmrg 
85c0a68be4Smrg       if (__isnan(__x))
8636ac495dSmrg         return std::numeric_limits<_Tp>::quiet_NaN();
8736ac495dSmrg       else if (__x == +_Tp(1))
8836ac495dSmrg         return +_Tp(1);
8936ac495dSmrg       else if (__x == -_Tp(1))
9036ac495dSmrg         return (__l % 2 == 1 ? -_Tp(1) : +_Tp(1));
9136ac495dSmrg       else
9236ac495dSmrg         {
9336ac495dSmrg           _Tp __p_lm2 = _Tp(1);
9436ac495dSmrg           if (__l == 0)
9536ac495dSmrg             return __p_lm2;
9636ac495dSmrg 
9736ac495dSmrg           _Tp __p_lm1 = __x;
9836ac495dSmrg           if (__l == 1)
9936ac495dSmrg             return __p_lm1;
10036ac495dSmrg 
10136ac495dSmrg           _Tp __p_l = 0;
10236ac495dSmrg           for (unsigned int __ll = 2; __ll <= __l; ++__ll)
10336ac495dSmrg             {
10436ac495dSmrg               //  This arrangement is supposed to be better for roundoff
10536ac495dSmrg               //  protection, Arfken, 2nd Ed, Eq 12.17a.
10636ac495dSmrg               __p_l = _Tp(2) * __x * __p_lm1 - __p_lm2
10736ac495dSmrg                     - (__x * __p_lm1 - __p_lm2) / _Tp(__ll);
10836ac495dSmrg               __p_lm2 = __p_lm1;
10936ac495dSmrg               __p_lm1 = __p_l;
11036ac495dSmrg             }
11136ac495dSmrg 
11236ac495dSmrg           return __p_l;
11336ac495dSmrg         }
11436ac495dSmrg     }
11536ac495dSmrg 
11636ac495dSmrg 
11736ac495dSmrg     /**
11836ac495dSmrg      *   @brief  Return the associated Legendre function by recursion
11936ac495dSmrg      *           on @f$ l @f$.
12036ac495dSmrg      *
12136ac495dSmrg      *   The associated Legendre function is derived from the Legendre function
12236ac495dSmrg      *   @f$ P_l(x) @f$ by the Rodrigues formula:
12336ac495dSmrg      *   @f[
12436ac495dSmrg      *     P_l^m(x) = (1 - x^2)^{m/2}\frac{d^m}{dx^m}P_l(x)
12536ac495dSmrg      *   @f]
126c0a68be4Smrg      *   @note @f$ P_l^m(x) = 0 @f$ if @f$ m > l @f$.
12736ac495dSmrg      *
128c0a68be4Smrg      *   @param  l  The degree of the associated Legendre function.
12936ac495dSmrg      *              @f$ l >= 0 @f$.
13036ac495dSmrg      *   @param  m  The order of the associated Legendre function.
13136ac495dSmrg      *   @param  x  The argument of the associated Legendre function.
13236ac495dSmrg      *              @f$ |x| <= 1 @f$.
133c0a68be4Smrg      *   @param  phase  The phase of the associated Legendre function.
134c0a68be4Smrg      *                  Use -1 for the Condon-Shortley phase convention.
13536ac495dSmrg      */
13636ac495dSmrg     template<typename _Tp>
13736ac495dSmrg     _Tp
__assoc_legendre_p(unsigned int __l,unsigned int __m,_Tp __x,_Tp __phase=_Tp (+1))138c0a68be4Smrg     __assoc_legendre_p(unsigned int __l, unsigned int __m, _Tp __x,
139c0a68be4Smrg 		       _Tp __phase = _Tp(+1))
14036ac495dSmrg     {
14136ac495dSmrg 
142c0a68be4Smrg       if (__m > __l)
143c0a68be4Smrg         return _Tp(0);
14436ac495dSmrg       else if (__isnan(__x))
14536ac495dSmrg         return std::numeric_limits<_Tp>::quiet_NaN();
14636ac495dSmrg       else if (__m == 0)
14736ac495dSmrg         return __poly_legendre_p(__l, __x);
14836ac495dSmrg       else
14936ac495dSmrg         {
15036ac495dSmrg           _Tp __p_mm = _Tp(1);
15136ac495dSmrg           if (__m > 0)
15236ac495dSmrg             {
15336ac495dSmrg               //  Two square roots seem more accurate more of the time
15436ac495dSmrg               //  than just one.
15536ac495dSmrg               _Tp __root = std::sqrt(_Tp(1) - __x) * std::sqrt(_Tp(1) + __x);
15636ac495dSmrg               _Tp __fact = _Tp(1);
15736ac495dSmrg               for (unsigned int __i = 1; __i <= __m; ++__i)
15836ac495dSmrg                 {
159c0a68be4Smrg                   __p_mm *= __phase * __fact * __root;
16036ac495dSmrg                   __fact += _Tp(2);
16136ac495dSmrg                 }
16236ac495dSmrg             }
16336ac495dSmrg           if (__l == __m)
16436ac495dSmrg             return __p_mm;
16536ac495dSmrg 
16636ac495dSmrg           _Tp __p_mp1m = _Tp(2 * __m + 1) * __x * __p_mm;
16736ac495dSmrg           if (__l == __m + 1)
16836ac495dSmrg             return __p_mp1m;
16936ac495dSmrg 
17036ac495dSmrg           _Tp __p_lm2m = __p_mm;
17136ac495dSmrg           _Tp __P_lm1m = __p_mp1m;
17236ac495dSmrg           _Tp __p_lm = _Tp(0);
17336ac495dSmrg           for (unsigned int __j = __m + 2; __j <= __l; ++__j)
17436ac495dSmrg             {
17536ac495dSmrg               __p_lm = (_Tp(2 * __j - 1) * __x * __P_lm1m
17636ac495dSmrg                       - _Tp(__j + __m - 1) * __p_lm2m) / _Tp(__j - __m);
17736ac495dSmrg               __p_lm2m = __P_lm1m;
17836ac495dSmrg               __P_lm1m = __p_lm;
17936ac495dSmrg             }
18036ac495dSmrg 
18136ac495dSmrg           return __p_lm;
18236ac495dSmrg         }
18336ac495dSmrg     }
18436ac495dSmrg 
18536ac495dSmrg 
18636ac495dSmrg     /**
18736ac495dSmrg      *   @brief  Return the spherical associated Legendre function.
18836ac495dSmrg      *
18936ac495dSmrg      *   The spherical associated Legendre function of @f$ l @f$, @f$ m @f$,
19036ac495dSmrg      *   and @f$ \theta @f$ is defined as @f$ Y_l^m(\theta,0) @f$ where
19136ac495dSmrg      *   @f[
19236ac495dSmrg      *      Y_l^m(\theta,\phi) = (-1)^m[\frac{(2l+1)}{4\pi}
19336ac495dSmrg      *                                  \frac{(l-m)!}{(l+m)!}]
19436ac495dSmrg      *                     P_l^m(\cos\theta) \exp^{im\phi}
19536ac495dSmrg      *   @f]
19636ac495dSmrg      *   is the spherical harmonic function and @f$ P_l^m(x) @f$ is the
19736ac495dSmrg      *   associated Legendre function.
19836ac495dSmrg      *
19936ac495dSmrg      *   This function differs from the associated Legendre function by
20036ac495dSmrg      *   argument (@f$x = \cos(\theta)@f$) and by a normalization factor
20136ac495dSmrg      *   but this factor is rather large for large @f$ l @f$ and @f$ m @f$
20236ac495dSmrg      *   and so this function is stable for larger differences of @f$ l @f$
20336ac495dSmrg      *   and @f$ m @f$.
204c0a68be4Smrg      *   @note Unlike the case for __assoc_legendre_p the Condon-Shortley
205c0a68be4Smrg      *         phase factor @f$ (-1)^m @f$ is present here.
206c0a68be4Smrg      *   @note @f$ Y_l^m(\theta) = 0 @f$ if @f$ m > l @f$.
20736ac495dSmrg      *
208c0a68be4Smrg      *   @param  l  The degree of the spherical associated Legendre function.
20936ac495dSmrg      *              @f$ l >= 0 @f$.
21036ac495dSmrg      *   @param  m  The order of the spherical associated Legendre function.
21136ac495dSmrg      *   @param  theta  The radian angle argument of the spherical associated
21236ac495dSmrg      *                  Legendre function.
21336ac495dSmrg      */
21436ac495dSmrg     template <typename _Tp>
21536ac495dSmrg     _Tp
__sph_legendre(unsigned int __l,unsigned int __m,_Tp __theta)21636ac495dSmrg     __sph_legendre(unsigned int __l, unsigned int __m, _Tp __theta)
21736ac495dSmrg     {
21836ac495dSmrg       if (__isnan(__theta))
21936ac495dSmrg         return std::numeric_limits<_Tp>::quiet_NaN();
22036ac495dSmrg 
22136ac495dSmrg       const _Tp __x = std::cos(__theta);
22236ac495dSmrg 
223c0a68be4Smrg       if (__m > __l)
224c0a68be4Smrg         return _Tp(0);
22536ac495dSmrg       else if (__m == 0)
22636ac495dSmrg         {
22736ac495dSmrg           _Tp __P = __poly_legendre_p(__l, __x);
22836ac495dSmrg           _Tp __fact = std::sqrt(_Tp(2 * __l + 1)
22936ac495dSmrg                      / (_Tp(4) * __numeric_constants<_Tp>::__pi()));
23036ac495dSmrg           __P *= __fact;
23136ac495dSmrg           return __P;
23236ac495dSmrg         }
23336ac495dSmrg       else if (__x == _Tp(1) || __x == -_Tp(1))
23436ac495dSmrg         {
23536ac495dSmrg           //  m > 0 here
23636ac495dSmrg           return _Tp(0);
23736ac495dSmrg         }
23836ac495dSmrg       else
23936ac495dSmrg         {
24036ac495dSmrg           // m > 0 and |x| < 1 here
24136ac495dSmrg 
24236ac495dSmrg           // Starting value for recursion.
24336ac495dSmrg           // Y_m^m(x) = sqrt( (2m+1)/(4pi m) gamma(m+1/2)/gamma(m) )
24436ac495dSmrg           //             (-1)^m (1-x^2)^(m/2) / pi^(1/4)
24536ac495dSmrg           const _Tp __sgn = ( __m % 2 == 1 ? -_Tp(1) : _Tp(1));
24636ac495dSmrg           const _Tp __y_mp1m_factor = __x * std::sqrt(_Tp(2 * __m + 3));
24736ac495dSmrg #if _GLIBCXX_USE_C99_MATH_TR1
24836ac495dSmrg           const _Tp __lncirc = _GLIBCXX_MATH_NS::log1p(-__x * __x);
24936ac495dSmrg #else
25036ac495dSmrg           const _Tp __lncirc = std::log(_Tp(1) - __x * __x);
25136ac495dSmrg #endif
25236ac495dSmrg           //  Gamma(m+1/2) / Gamma(m)
25336ac495dSmrg #if _GLIBCXX_USE_C99_MATH_TR1
25436ac495dSmrg           const _Tp __lnpoch = _GLIBCXX_MATH_NS::lgamma(_Tp(__m + _Tp(0.5L)))
25536ac495dSmrg                              - _GLIBCXX_MATH_NS::lgamma(_Tp(__m));
25636ac495dSmrg #else
25736ac495dSmrg           const _Tp __lnpoch = __log_gamma(_Tp(__m + _Tp(0.5L)))
25836ac495dSmrg                              - __log_gamma(_Tp(__m));
25936ac495dSmrg #endif
26036ac495dSmrg           const _Tp __lnpre_val =
26136ac495dSmrg                     -_Tp(0.25L) * __numeric_constants<_Tp>::__lnpi()
26236ac495dSmrg                     + _Tp(0.5L) * (__lnpoch + __m * __lncirc);
263c0a68be4Smrg           const _Tp __sr = std::sqrt((_Tp(2) + _Tp(1) / __m)
26436ac495dSmrg                          / (_Tp(4) * __numeric_constants<_Tp>::__pi()));
26536ac495dSmrg           _Tp __y_mm = __sgn * __sr * std::exp(__lnpre_val);
26636ac495dSmrg           _Tp __y_mp1m = __y_mp1m_factor * __y_mm;
26736ac495dSmrg 
26836ac495dSmrg           if (__l == __m)
26936ac495dSmrg             return __y_mm;
27036ac495dSmrg           else if (__l == __m + 1)
27136ac495dSmrg             return __y_mp1m;
27236ac495dSmrg           else
27336ac495dSmrg             {
27436ac495dSmrg               _Tp __y_lm = _Tp(0);
27536ac495dSmrg 
27636ac495dSmrg               // Compute Y_l^m, l > m+1, upward recursion on l.
277*8feb0f0bSmrg               for (unsigned int __ll = __m + 2; __ll <= __l; ++__ll)
27836ac495dSmrg                 {
27936ac495dSmrg                   const _Tp __rat1 = _Tp(__ll - __m) / _Tp(__ll + __m);
28036ac495dSmrg                   const _Tp __rat2 = _Tp(__ll - __m - 1) / _Tp(__ll + __m - 1);
28136ac495dSmrg                   const _Tp __fact1 = std::sqrt(__rat1 * _Tp(2 * __ll + 1)
28236ac495dSmrg                                                        * _Tp(2 * __ll - 1));
28336ac495dSmrg                   const _Tp __fact2 = std::sqrt(__rat1 * __rat2 * _Tp(2 * __ll + 1)
28436ac495dSmrg                                                                 / _Tp(2 * __ll - 3));
28536ac495dSmrg                   __y_lm = (__x * __y_mp1m * __fact1
28636ac495dSmrg                          - (__ll + __m - 1) * __y_mm * __fact2) / _Tp(__ll - __m);
28736ac495dSmrg                   __y_mm = __y_mp1m;
28836ac495dSmrg                   __y_mp1m = __y_lm;
28936ac495dSmrg                 }
29036ac495dSmrg 
29136ac495dSmrg               return __y_lm;
29236ac495dSmrg             }
29336ac495dSmrg         }
29436ac495dSmrg     }
29536ac495dSmrg   } // namespace __detail
29636ac495dSmrg #undef _GLIBCXX_MATH_NS
29736ac495dSmrg #if ! _GLIBCXX_USE_STD_SPEC_FUNCS && defined(_GLIBCXX_TR1_CMATH)
29836ac495dSmrg } // namespace tr1
29936ac495dSmrg #endif
300a2dc1f3fSmrg 
301a2dc1f3fSmrg _GLIBCXX_END_NAMESPACE_VERSION
30236ac495dSmrg }
30336ac495dSmrg 
30436ac495dSmrg #endif // _GLIBCXX_TR1_LEGENDRE_FUNCTION_TCC
305