OSDN Git Service

2008-01-31 Ralf Wildenhues <Ralf.Wildenhues@gmx.de>
[pf3gnuchains/gcc-fork.git] / libstdc++-v3 / include / tr1 / poly_laguerre.tcc
1 // Special functions -*- C++ -*-
2
3 // Copyright (C) 2006-2007
4 // Free Software Foundation, Inc.
5 //
6 // This file is part of the GNU ISO C++ Library.  This library is free
7 // software; you can redistribute it and/or modify it under the
8 // terms of the GNU General Public License as published by the
9 // Free Software Foundation; either version 2, or (at your option)
10 // any later version.
11 //
12 // This library is distributed in the hope that it will be useful,
13 // but WITHOUT ANY WARRANTY; without even the implied warranty of
14 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
15 // GNU General Public License for more details.
16 //
17 // You should have received a copy of the GNU General Public License along
18 // with this library; see the file COPYING.  If not, write to the Free
19 // Software Foundation, 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301,
20 // USA.
21 //
22 // As a special exception, you may use this file as part of a free software
23 // library without restriction.  Specifically, if other files instantiate
24 // templates or use macros or inline functions from this file, or you compile
25 // this file and link it with other files to produce an executable, this
26 // file does not by itself cause the resulting executable to be covered by
27 // the GNU General Public License.  This exception does not however
28 // invalidate any other reasons why the executable file might be covered by
29 // the GNU General Public License.
30
31 /** @file tr1/poly_laguerre.tcc
32  *  This is an internal header file, included by other library headers.
33  *  You should not attempt to use it directly.
34  */
35
36 //
37 // ISO C++ 14882 TR1: 5.2  Special functions
38 //
39
40 // Written by Edward Smith-Rowland based on:
41 //   (1) Handbook of Mathematical Functions,
42 //       Ed. Milton Abramowitz and Irene A. Stegun,
43 //       Dover Publications,
44 //       Section 13, pp. 509-510, Section 22 pp. 773-802
45 //   (2) The Gnu Scientific Library, http://www.gnu.org/software/gsl
46
47 #ifndef _GLIBCXX_TR1_POLY_LAGUERRE_TCC
48 #define _GLIBCXX_TR1_POLY_LAGUERRE_TCC 1
49
50 namespace std
51 {
52 namespace tr1
53 {
54
55   // [5.2] Special functions
56
57   /**
58    * @ingroup tr1_math_spec_func
59    * @{
60    */
61
62   //
63   // Implementation-space details.
64   //
65   namespace __detail
66   {
67
68
69     /**
70      *   @brief This routine returns the associated Laguerre polynomial 
71      *          of order @f$ n @f$, degree @f$ \alpha @f$ for large n.
72      *   Abramowitz & Stegun, 13.5.21
73      *
74      *   @param __n The order of the Laguerre function.
75      *   @param __alpha The degree of the Laguerre function.
76      *   @param __x The argument of the Laguerre function.
77      *   @return The value of the Laguerre function of order n,
78      *           degree @f$ \alpha @f$, and argument x.
79      *
80      *  This is from the GNU Scientific Library.
81      */
82     template<typename _Tpa, typename _Tp>
83     _Tp
84     __poly_laguerre_large_n(const unsigned __n, const _Tpa __alpha1,
85                             const _Tp __x)
86     {
87       const _Tp __a = -_Tp(__n);
88       const _Tp __b = _Tp(__alpha1) + _Tp(1);
89       const _Tp __eta = _Tp(2) * __b - _Tp(4) * __a;
90       const _Tp __cos2th = __x / __eta;
91       const _Tp __sin2th = _Tp(1) - __cos2th;
92       const _Tp __th = std::acos(std::sqrt(__cos2th));
93       const _Tp __pre_h = __numeric_constants<_Tp>::__pi_2()
94                         * __numeric_constants<_Tp>::__pi_2()
95                         * __eta * __eta * __cos2th * __sin2th;
96
97 #if _GLIBCXX_USE_C99_MATH_TR1
98       const _Tp __lg_b = std::tr1::lgamma(_Tp(__n) + __b);
99       const _Tp __lnfact = std::tr1::lgamma(_Tp(__n + 1));
100 #else
101       const _Tp __lg_b = __log_gamma(_Tp(__n) + __b);
102       const _Tp __lnfact = __log_gamma(_Tp(__n + 1));
103 #endif
104
105       _Tp __pre_term1 = _Tp(0.5L) * (_Tp(1) - __b)
106                       * std::log(_Tp(0.25L) * __x * __eta);
107       _Tp __pre_term2 = _Tp(0.25L) * std::log(__pre_h);
108       _Tp __lnpre = __lg_b - __lnfact + _Tp(0.5L) * __x
109                       + __pre_term1 - __pre_term2;
110       _Tp __ser_term1 = std::sin(__a * __numeric_constants<_Tp>::__pi());
111       _Tp __ser_term2 = std::sin(_Tp(0.25L) * __eta
112                               * (_Tp(2) * __th
113                                - std::sin(_Tp(2) * __th))
114                                + __numeric_constants<_Tp>::__pi_4());
115       _Tp __ser = __ser_term1 + __ser_term2;
116
117       return std::exp(__lnpre) * __ser;
118     }
119
120
121     /**
122      *  @brief  Evaluate the polynomial based on the confluent hypergeometric
123      *          function in a safe way, with no restriction on the arguments.
124      *
125      *   The associated Laguerre function is defined by
126      *   @f[
127      *       L_n^\alpha(x) = \frac{(\alpha + 1)_n}{n!}
128      *                       _1F_1(-n; \alpha + 1; x)
129      *   @f]
130      *   where @f$ (\alpha)_n @f$ is the Pochhammer symbol and
131      *   @f$ _1F_1(a; c; x) @f$ is the confluent hypergeometric function.
132      *
133      *  This function assumes x != 0.
134      *
135      *  This is from the GNU Scientific Library.
136      */
137     template<typename _Tpa, typename _Tp>
138     _Tp
139     __poly_laguerre_hyperg(const unsigned int __n, const _Tpa __alpha1,
140                            const _Tp __x)
141     {
142       const _Tp __b = _Tp(__alpha1) + _Tp(1);
143       const _Tp __mx = -__x;
144       const _Tp __tc_sgn = (__x < _Tp(0) ? _Tp(1)
145                          : ((__n % 2 == 1) ? -_Tp(1) : _Tp(1)));
146       //  Get |x|^n/n!
147       _Tp __tc = _Tp(1);
148       const _Tp __ax = std::abs(__x);
149       for (unsigned int __k = 1; __k <= __n; ++__k)
150         __tc *= (__ax / __k);
151
152       _Tp __term = __tc * __tc_sgn;
153       _Tp __sum = __term;
154       for (int __k = int(__n) - 1; __k >= 0; --__k)
155         {
156           __term *= ((__b + _Tp(__k)) / _Tp(int(__n) - __k))
157                   * _Tp(__k + 1) / __mx;
158           __sum += __term;
159         }
160
161       return __sum;
162     }
163
164
165     /**
166      *   @brief This routine returns the associated Laguerre polynomial 
167      *          of order @f$ n @f$, degree @f$ \alpha @f$: @f$ L_n^\alpha(x) @f$
168      *          by recursion.
169      *
170      *   The associated Laguerre function is defined by
171      *   @f[
172      *       L_n^\alpha(x) = \frac{(\alpha + 1)_n}{n!}
173      *                       _1F_1(-n; \alpha + 1; x)
174      *   @f]
175      *   where @f$ (\alpha)_n @f$ is the Pochhammer symbol and
176      *   @f$ _1F_1(a; c; x) @f$ is the confluent hypergeometric function.
177      *
178      *   The associated Laguerre polynomial is defined for integral
179      *   @f$ \alpha = m @f$ by:
180      *   @f[
181      *       L_n^m(x) = (-1)^m \frac{d^m}{dx^m} L_{n + m}(x)
182      *   @f]
183      *   where the Laguerre polynomial is defined by:
184      *   @f[
185      *       L_n(x) = \frac{e^x}{n!} \frac{d^n}{dx^n} (x^ne^{-x})
186      *   @f]
187      *
188      *   @param __n The order of the Laguerre function.
189      *   @param __alpha The degree of the Laguerre function.
190      *   @param __x The argument of the Laguerre function.
191      *   @return The value of the Laguerre function of order n,
192      *           degree @f$ \alpha @f$, and argument x.
193      */
194     template<typename _Tpa, typename _Tp>
195     _Tp
196     __poly_laguerre_recursion(const unsigned int __n,
197                               const _Tpa __alpha1, const _Tp __x)
198     {
199       //   Compute l_0.
200       _Tp __l_0 = _Tp(1);
201       if  (__n == 0)
202         return __l_0;
203
204       //  Compute l_1^alpha.
205       _Tp __l_1 = -__x + _Tp(1) + _Tp(__alpha1);
206       if  (__n == 1)
207         return __l_1;
208
209       //  Compute l_n^alpha by recursion on n.
210       _Tp __l_n2 = __l_0;
211       _Tp __l_n1 = __l_1;
212       _Tp __l_n = _Tp(0);
213       for  (unsigned int __nn = 2; __nn <= __n; ++__nn)
214         {
215             __l_n = (_Tp(2 * __nn - 1) + _Tp(__alpha1) - __x)
216                   * __l_n1 / _Tp(__nn)
217                   - (_Tp(__nn - 1) + _Tp(__alpha1)) * __l_n2 / _Tp(__nn);
218             __l_n2 = __l_n1;
219             __l_n1 = __l_n;
220         }
221
222       return __l_n;
223     }
224
225
226     /**
227      *   @brief This routine returns the associated Laguerre polynomial
228      *          of order n, degree @f$ \alpha @f$: @f$ L_n^alpha(x) @f$.
229      *
230      *   The associated Laguerre function is defined by
231      *   @f[
232      *       L_n^\alpha(x) = \frac{(\alpha + 1)_n}{n!}
233      *                       _1F_1(-n; \alpha + 1; x)
234      *   @f]
235      *   where @f$ (\alpha)_n @f$ is the Pochhammer symbol and
236      *   @f$ _1F_1(a; c; x) @f$ is the confluent hypergeometric function.
237      *
238      *   The associated Laguerre polynomial is defined for integral
239      *   @f$ \alpha = m @f$ by:
240      *   @f[
241      *       L_n^m(x) = (-1)^m \frac{d^m}{dx^m} L_{n + m}(x)
242      *   @f]
243      *   where the Laguerre polynomial is defined by:
244      *   @f[
245      *       L_n(x) = \frac{e^x}{n!} \frac{d^n}{dx^n} (x^ne^{-x})
246      *   @f]
247      *
248      *   @param __n The order of the Laguerre function.
249      *   @param __alpha The degree of the Laguerre function.
250      *   @param __x The argument of the Laguerre function.
251      *   @return The value of the Laguerre function of order n,
252      *           degree @f$ \alpha @f$, and argument x.
253      */
254     template<typename _Tpa, typename _Tp>
255     inline _Tp
256     __poly_laguerre(const unsigned int __n, const _Tpa __alpha1,
257                     const _Tp __x)
258     {
259       if (__x < _Tp(0))
260         std::__throw_domain_error(__N("Negative argument "
261                                       "in __poly_laguerre."));
262       //  Return NaN on NaN input.
263       else if (__isnan(__x))
264         return std::numeric_limits<_Tp>::quiet_NaN();
265       else if (__n == 0)
266         return _Tp(1);
267       else if (__n == 1)
268         return _Tp(1) + _Tp(__alpha1) - __x;
269       else if (__x == _Tp(0))
270         {
271           _Tp __prod = _Tp(__alpha1) + _Tp(1);
272           for (unsigned int __k = 2; __k <= __n; ++__k)
273             __prod *= (_Tp(__alpha1) + _Tp(__k)) / _Tp(__k);
274           return __prod;
275         }
276       else if (__n > 10000000 && _Tp(__alpha1) > -_Tp(1)
277             && __x < _Tp(2) * (_Tp(__alpha1) + _Tp(1)) + _Tp(4 * __n))
278         return __poly_laguerre_large_n(__n, __alpha1, __x);
279       else if (_Tp(__alpha1) >= _Tp(0)
280            || (__x > _Tp(0) && _Tp(__alpha1) < -_Tp(__n + 1)))
281         return __poly_laguerre_recursion(__n, __alpha1, __x);
282       else
283         return __poly_laguerre_hyperg(__n, __alpha1, __x);
284     }
285
286
287     /**
288      *   @brief This routine returns the associated Laguerre polynomial
289      *          of order n, degree m: @f$ L_n^m @f$.
290      *
291      *   The associated Laguerre polynomial is defined for integral
292      *   @f$ \alpha = m @f$ by:
293      *   @f[
294      *       L_n^m(x) = (-1)^m \frac{d^m}{dx^m} L_{n + m}(x)
295      *   @f]
296      *   where the Laguerre polynomial is defined by:
297      *   @f[
298      *       L_n(x) = \frac{e^x}{n!} \frac{d^n}{dx^n} (x^ne^{-x})
299      *   @f]
300      *
301      *   @param __n The order of the Laguerre polynomial.
302      *   @param __m The degree of the Laguerre polynomial.
303      *   @param __x The argument of the Laguerre polynomial.
304      *   @return The value of the associated Laguerre polynomial of order n,
305      *           degree m, and argument x.
306      */
307     template<typename _Tp>
308     inline _Tp
309     __assoc_laguerre(const unsigned int __n, const unsigned int __m,
310                      const _Tp __x)
311     {
312       return __poly_laguerre<unsigned int, _Tp>(__n, __m, __x);
313     }
314
315
316     /**
317      *   @brief This routine returns the associated Laguerre polynomial
318      *          of order n: @f$ L_n(x) @f$.
319      *
320      *   The Laguerre polynomial is defined by:
321      *   @f[
322      *       L_n(x) = \frac{e^x}{n!} \frac{d^n}{dx^n} (x^ne^{-x})
323      *   @f]
324      *
325      *   @param __n The order of the Laguerre polynomial.
326      *   @param __x The argument of the Laguerre polynomial.
327      *   @return The value of the Laguerre polynomial of order n
328      *           and argument x.
329      */
330     template<typename _Tp>
331     inline _Tp
332     __laguerre(const unsigned int __n, const _Tp __x)
333     {
334       return __poly_laguerre<unsigned int, _Tp>(__n, 0, __x);
335     }
336
337   } // namespace std::tr1::__detail
338
339   /* @} */ // group tr1_math_spec_func
340
341 }
342 }
343
344 #endif // _GLIBCXX_TR1_POLY_LAGUERRE_TCC