libstdc++
bessel_function.tcc
Go to the documentation of this file.
1 // Special functions -*- C++ -*-
2 
3 // Copyright (C) 2006, 2007, 2008, 2009
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 3, 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 // Under Section 7 of GPL version 3, you are granted additional
18 // permissions described in the GCC Runtime Library Exception, version
19 // 3.1, as published by the Free Software Foundation.
20 
21 // You should have received a copy of the GNU General Public License and
22 // a copy of the GCC Runtime Library Exception along with this program;
23 // see the files COPYING3 and COPYING.RUNTIME respectively. If not, see
24 // <http://www.gnu.org/licenses/>.
25 
26 /** @file tr1/bessel_function.tcc
27  * This is an internal header file, included by other library headers.
28  * You should not attempt to use it directly.
29  */
30 
31 //
32 // ISO C++ 14882 TR1: 5.2 Special functions
33 //
34 
35 // Written by Edward Smith-Rowland.
36 //
37 // References:
38 // (1) Handbook of Mathematical Functions,
39 // ed. Milton Abramowitz and Irene A. Stegun,
40 // Dover Publications,
41 // Section 9, pp. 355-434, Section 10 pp. 435-478
42 // (2) The Gnu Scientific Library, http://www.gnu.org/software/gsl
43 // (3) Numerical Recipes in C, by W. H. Press, S. A. Teukolsky,
44 // W. T. Vetterling, B. P. Flannery, Cambridge University Press (1992),
45 // 2nd ed, pp. 240-245
46 
47 #ifndef _GLIBCXX_TR1_BESSEL_FUNCTION_TCC
48 #define _GLIBCXX_TR1_BESSEL_FUNCTION_TCC 1
49 
50 #include "special_function_util.h"
51 
52 namespace std
53 {
54 namespace tr1
55 {
56 
57  // [5.2] Special functions
58 
59  // Implementation-space details.
60  namespace __detail
61  {
62 
63  /**
64  * @brief Compute the gamma functions required by the Temme series
65  * expansions of @f$ N_\nu(x) @f$ and @f$ K_\nu(x) @f$.
66  * @f[
67  * \Gamma_1 = \frac{1}{2\mu}
68  * [\frac{1}{\Gamma(1 - \mu)} - \frac{1}{\Gamma(1 + \mu)}]
69  * @f]
70  * and
71  * @f[
72  * \Gamma_2 = \frac{1}{2}
73  * [\frac{1}{\Gamma(1 - \mu)} + \frac{1}{\Gamma(1 + \mu)}]
74  * @f]
75  * where @f$ -1/2 <= \mu <= 1/2 @f$ is @f$ \mu = \nu - N @f$ and @f$ N @f$.
76  * is the nearest integer to @f$ \nu @f$.
77  * The values of \f$ \Gamma(1 + \mu) \f$ and \f$ \Gamma(1 - \mu) \f$
78  * are returned as well.
79  *
80  * The accuracy requirements on this are exquisite.
81  *
82  * @param __mu The input parameter of the gamma functions.
83  * @param __gam1 The output function \f$ \Gamma_1(\mu) \f$
84  * @param __gam2 The output function \f$ \Gamma_2(\mu) \f$
85  * @param __gampl The output function \f$ \Gamma(1 + \mu) \f$
86  * @param __gammi The output function \f$ \Gamma(1 - \mu) \f$
87  */
88  template <typename _Tp>
89  void
90  __gamma_temme(const _Tp __mu,
91  _Tp & __gam1, _Tp & __gam2, _Tp & __gampl, _Tp & __gammi)
92  {
93 #if _GLIBCXX_USE_C99_MATH_TR1
94  __gampl = _Tp(1) / std::tr1::tgamma(_Tp(1) + __mu);
95  __gammi = _Tp(1) / std::tr1::tgamma(_Tp(1) - __mu);
96 #else
97  __gampl = _Tp(1) / __gamma(_Tp(1) + __mu);
98  __gammi = _Tp(1) / __gamma(_Tp(1) - __mu);
99 #endif
100 
102  __gam1 = -_Tp(__numeric_constants<_Tp>::__gamma_e());
103  else
104  __gam1 = (__gammi - __gampl) / (_Tp(2) * __mu);
105 
106  __gam2 = (__gammi + __gampl) / (_Tp(2));
107 
108  return;
109  }
110 
111 
112  /**
113  * @brief Compute the Bessel @f$ J_\nu(x) @f$ and Neumann
114  * @f$ N_\nu(x) @f$ functions and their first derivatives
115  * @f$ J'_\nu(x) @f$ and @f$ N'_\nu(x) @f$ respectively.
116  * These four functions are computed together for numerical
117  * stability.
118  *
119  * @param __nu The order of the Bessel functions.
120  * @param __x The argument of the Bessel functions.
121  * @param __Jnu The output Bessel function of the first kind.
122  * @param __Nnu The output Neumann function (Bessel function of the second kind).
123  * @param __Jpnu The output derivative of the Bessel function of the first kind.
124  * @param __Npnu The output derivative of the Neumann function.
125  */
126  template <typename _Tp>
127  void
128  __bessel_jn(const _Tp __nu, const _Tp __x,
129  _Tp & __Jnu, _Tp & __Nnu, _Tp & __Jpnu, _Tp & __Npnu)
130  {
131  if (__x == _Tp(0))
132  {
133  if (__nu == _Tp(0))
134  {
135  __Jnu = _Tp(1);
136  __Jpnu = _Tp(0);
137  }
138  else if (__nu == _Tp(1))
139  {
140  __Jnu = _Tp(0);
141  __Jpnu = _Tp(0.5L);
142  }
143  else
144  {
145  __Jnu = _Tp(0);
146  __Jpnu = _Tp(0);
147  }
150  return;
151  }
152 
153  const _Tp __eps = std::numeric_limits<_Tp>::epsilon();
154  // When the multiplier is N i.e.
155  // fp_min = N * min()
156  // Then J_0 and N_0 tank at x = 8 * N (J_0 = 0 and N_0 = nan)!
157  //const _Tp __fp_min = _Tp(20) * std::numeric_limits<_Tp>::min();
158  const _Tp __fp_min = std::sqrt(std::numeric_limits<_Tp>::min());
159  const int __max_iter = 15000;
160  const _Tp __x_min = _Tp(2);
161 
162  const int __nl = (__x < __x_min
163  ? static_cast<int>(__nu + _Tp(0.5L))
164  : std::max(0, static_cast<int>(__nu - __x + _Tp(1.5L))));
165 
166  const _Tp __mu = __nu - __nl;
167  const _Tp __mu2 = __mu * __mu;
168  const _Tp __xi = _Tp(1) / __x;
169  const _Tp __xi2 = _Tp(2) * __xi;
170  _Tp __w = __xi2 / __numeric_constants<_Tp>::__pi();
171  int __isign = 1;
172  _Tp __h = __nu * __xi;
173  if (__h < __fp_min)
174  __h = __fp_min;
175  _Tp __b = __xi2 * __nu;
176  _Tp __d = _Tp(0);
177  _Tp __c = __h;
178  int __i;
179  for (__i = 1; __i <= __max_iter; ++__i)
180  {
181  __b += __xi2;
182  __d = __b - __d;
183  if (std::abs(__d) < __fp_min)
184  __d = __fp_min;
185  __c = __b - _Tp(1) / __c;
186  if (std::abs(__c) < __fp_min)
187  __c = __fp_min;
188  __d = _Tp(1) / __d;
189  const _Tp __del = __c * __d;
190  __h *= __del;
191  if (__d < _Tp(0))
192  __isign = -__isign;
193  if (std::abs(__del - _Tp(1)) < __eps)
194  break;
195  }
196  if (__i > __max_iter)
197  std::__throw_runtime_error(__N("Argument x too large in __bessel_jn; "
198  "try asymptotic expansion."));
199  _Tp __Jnul = __isign * __fp_min;
200  _Tp __Jpnul = __h * __Jnul;
201  _Tp __Jnul1 = __Jnul;
202  _Tp __Jpnu1 = __Jpnul;
203  _Tp __fact = __nu * __xi;
204  for ( int __l = __nl; __l >= 1; --__l )
205  {
206  const _Tp __Jnutemp = __fact * __Jnul + __Jpnul;
207  __fact -= __xi;
208  __Jpnul = __fact * __Jnutemp - __Jnul;
209  __Jnul = __Jnutemp;
210  }
211  if (__Jnul == _Tp(0))
212  __Jnul = __eps;
213  _Tp __f= __Jpnul / __Jnul;
214  _Tp __Nmu, __Nnu1, __Npmu, __Jmu;
215  if (__x < __x_min)
216  {
217  const _Tp __x2 = __x / _Tp(2);
218  const _Tp __pimu = __numeric_constants<_Tp>::__pi() * __mu;
219  _Tp __fact = (std::abs(__pimu) < __eps
220  ? _Tp(1) : __pimu / std::sin(__pimu));
221  _Tp __d = -std::log(__x2);
222  _Tp __e = __mu * __d;
223  _Tp __fact2 = (std::abs(__e) < __eps
224  ? _Tp(1) : std::sinh(__e) / __e);
225  _Tp __gam1, __gam2, __gampl, __gammi;
226  __gamma_temme(__mu, __gam1, __gam2, __gampl, __gammi);
227  _Tp __ff = (_Tp(2) / __numeric_constants<_Tp>::__pi())
228  * __fact * (__gam1 * std::cosh(__e) + __gam2 * __fact2 * __d);
229  __e = std::exp(__e);
230  _Tp __p = __e / (__numeric_constants<_Tp>::__pi() * __gampl);
231  _Tp __q = _Tp(1) / (__e * __numeric_constants<_Tp>::__pi() * __gammi);
232  const _Tp __pimu2 = __pimu / _Tp(2);
233  _Tp __fact3 = (std::abs(__pimu2) < __eps
234  ? _Tp(1) : std::sin(__pimu2) / __pimu2 );
235  _Tp __r = __numeric_constants<_Tp>::__pi() * __pimu2 * __fact3 * __fact3;
236  _Tp __c = _Tp(1);
237  __d = -__x2 * __x2;
238  _Tp __sum = __ff + __r * __q;
239  _Tp __sum1 = __p;
240  for (__i = 1; __i <= __max_iter; ++__i)
241  {
242  __ff = (__i * __ff + __p + __q) / (__i * __i - __mu2);
243  __c *= __d / _Tp(__i);
244  __p /= _Tp(__i) - __mu;
245  __q /= _Tp(__i) + __mu;
246  const _Tp __del = __c * (__ff + __r * __q);
247  __sum += __del;
248  const _Tp __del1 = __c * __p - __i * __del;
249  __sum1 += __del1;
250  if ( std::abs(__del) < __eps * (_Tp(1) + std::abs(__sum)) )
251  break;
252  }
253  if ( __i > __max_iter )
254  std::__throw_runtime_error(__N("Bessel y series failed to converge "
255  "in __bessel_jn."));
256  __Nmu = -__sum;
257  __Nnu1 = -__sum1 * __xi2;
258  __Npmu = __mu * __xi * __Nmu - __Nnu1;
259  __Jmu = __w / (__Npmu - __f * __Nmu);
260  }
261  else
262  {
263  _Tp __a = _Tp(0.25L) - __mu2;
264  _Tp __q = _Tp(1);
265  _Tp __p = -__xi / _Tp(2);
266  _Tp __br = _Tp(2) * __x;
267  _Tp __bi = _Tp(2);
268  _Tp __fact = __a * __xi / (__p * __p + __q * __q);
269  _Tp __cr = __br + __q * __fact;
270  _Tp __ci = __bi + __p * __fact;
271  _Tp __den = __br * __br + __bi * __bi;
272  _Tp __dr = __br / __den;
273  _Tp __di = -__bi / __den;
274  _Tp __dlr = __cr * __dr - __ci * __di;
275  _Tp __dli = __cr * __di + __ci * __dr;
276  _Tp __temp = __p * __dlr - __q * __dli;
277  __q = __p * __dli + __q * __dlr;
278  __p = __temp;
279  int __i;
280  for (__i = 2; __i <= __max_iter; ++__i)
281  {
282  __a += _Tp(2 * (__i - 1));
283  __bi += _Tp(2);
284  __dr = __a * __dr + __br;
285  __di = __a * __di + __bi;
286  if (std::abs(__dr) + std::abs(__di) < __fp_min)
287  __dr = __fp_min;
288  __fact = __a / (__cr * __cr + __ci * __ci);
289  __cr = __br + __cr * __fact;
290  __ci = __bi - __ci * __fact;
291  if (std::abs(__cr) + std::abs(__ci) < __fp_min)
292  __cr = __fp_min;
293  __den = __dr * __dr + __di * __di;
294  __dr /= __den;
295  __di /= -__den;
296  __dlr = __cr * __dr - __ci * __di;
297  __dli = __cr * __di + __ci * __dr;
298  __temp = __p * __dlr - __q * __dli;
299  __q = __p * __dli + __q * __dlr;
300  __p = __temp;
301  if (std::abs(__dlr - _Tp(1)) + std::abs(__dli) < __eps)
302  break;
303  }
304  if (__i > __max_iter)
305  std::__throw_runtime_error(__N("Lentz's method failed "
306  "in __bessel_jn."));
307  const _Tp __gam = (__p - __f) / __q;
308  __Jmu = std::sqrt(__w / ((__p - __f) * __gam + __q));
309 #if _GLIBCXX_USE_C99_MATH_TR1
310  __Jmu = std::tr1::copysign(__Jmu, __Jnul);
311 #else
312  if (__Jmu * __Jnul < _Tp(0))
313  __Jmu = -__Jmu;
314 #endif
315  __Nmu = __gam * __Jmu;
316  __Npmu = (__p + __q / __gam) * __Nmu;
317  __Nnu1 = __mu * __xi * __Nmu - __Npmu;
318  }
319  __fact = __Jmu / __Jnul;
320  __Jnu = __fact * __Jnul1;
321  __Jpnu = __fact * __Jpnu1;
322  for (__i = 1; __i <= __nl; ++__i)
323  {
324  const _Tp __Nnutemp = (__mu + __i) * __xi2 * __Nnu1 - __Nmu;
325  __Nmu = __Nnu1;
326  __Nnu1 = __Nnutemp;
327  }
328  __Nnu = __Nmu;
329  __Npnu = __nu * __xi * __Nmu - __Nnu1;
330 
331  return;
332  }
333 
334 
335  /**
336  * @brief This routine computes the asymptotic cylindrical Bessel
337  * and Neumann functions of order nu: \f$ J_{\nu} \f$,
338  * \f$ N_{\nu} \f$.
339  *
340  * References:
341  * (1) Handbook of Mathematical Functions,
342  * ed. Milton Abramowitz and Irene A. Stegun,
343  * Dover Publications,
344  * Section 9 p. 364, Equations 9.2.5-9.2.10
345  *
346  * @param __nu The order of the Bessel functions.
347  * @param __x The argument of the Bessel functions.
348  * @param __Jnu The output Bessel function of the first kind.
349  * @param __Nnu The output Neumann function (Bessel function of the second kind).
350  */
351  template <typename _Tp>
352  void
353  __cyl_bessel_jn_asymp(const _Tp __nu, const _Tp __x,
354  _Tp & __Jnu, _Tp & __Nnu)
355  {
356  const _Tp __coef = std::sqrt(_Tp(2)
357  / (__numeric_constants<_Tp>::__pi() * __x));
358  const _Tp __mu = _Tp(4) * __nu * __nu;
359  const _Tp __mum1 = __mu - _Tp(1);
360  const _Tp __mum9 = __mu - _Tp(9);
361  const _Tp __mum25 = __mu - _Tp(25);
362  const _Tp __mum49 = __mu - _Tp(49);
363  const _Tp __xx = _Tp(64) * __x * __x;
364  const _Tp __P = _Tp(1) - __mum1 * __mum9 / (_Tp(2) * __xx)
365  * (_Tp(1) - __mum25 * __mum49 / (_Tp(12) * __xx));
366  const _Tp __Q = __mum1 / (_Tp(8) * __x)
367  * (_Tp(1) - __mum9 * __mum25 / (_Tp(6) * __xx));
368 
369  const _Tp __chi = __x - (__nu + _Tp(0.5L))
371  const _Tp __c = std::cos(__chi);
372  const _Tp __s = std::sin(__chi);
373 
374  __Jnu = __coef * (__c * __P - __s * __Q);
375  __Nnu = __coef * (__s * __P + __c * __Q);
376 
377  return;
378  }
379 
380 
381  /**
382  * @brief This routine returns the cylindrical Bessel functions
383  * of order \f$ \nu \f$: \f$ J_{\nu} \f$ or \f$ I_{\nu} \f$
384  * by series expansion.
385  *
386  * The modified cylindrical Bessel function is:
387  * @f[
388  * Z_{\nu}(x) = \sum_{k=0}^{\infty}
389  * \frac{\sigma^k (x/2)^{\nu + 2k}}{k!\Gamma(\nu+k+1)}
390  * @f]
391  * where \f$ \sigma = +1 \f$ or\f$ -1 \f$ for
392  * \f$ Z = I \f$ or \f$ J \f$ respectively.
393  *
394  * See Abramowitz & Stegun, 9.1.10
395  * Abramowitz & Stegun, 9.6.7
396  * (1) Handbook of Mathematical Functions,
397  * ed. Milton Abramowitz and Irene A. Stegun,
398  * Dover Publications,
399  * Equation 9.1.10 p. 360 and Equation 9.6.10 p. 375
400  *
401  * @param __nu The order of the Bessel function.
402  * @param __x The argument of the Bessel function.
403  * @param __sgn The sign of the alternate terms
404  * -1 for the Bessel function of the first kind.
405  * +1 for the modified Bessel function of the first kind.
406  * @return The output Bessel function.
407  */
408  template <typename _Tp>
409  _Tp
410  __cyl_bessel_ij_series(const _Tp __nu, const _Tp __x, const _Tp __sgn,
411  const unsigned int __max_iter)
412  {
413 
414  const _Tp __x2 = __x / _Tp(2);
415  _Tp __fact = __nu * std::log(__x2);
416 #if _GLIBCXX_USE_C99_MATH_TR1
417  __fact -= std::tr1::lgamma(__nu + _Tp(1));
418 #else
419  __fact -= __log_gamma(__nu + _Tp(1));
420 #endif
421  __fact = std::exp(__fact);
422  const _Tp __xx4 = __sgn * __x2 * __x2;
423  _Tp __Jn = _Tp(1);
424  _Tp __term = _Tp(1);
425 
426  for (unsigned int __i = 1; __i < __max_iter; ++__i)
427  {
428  __term *= __xx4 / (_Tp(__i) * (__nu + _Tp(__i)));
429  __Jn += __term;
430  if (std::abs(__term / __Jn) < std::numeric_limits<_Tp>::epsilon())
431  break;
432  }
433 
434  return __fact * __Jn;
435  }
436 
437 
438  /**
439  * @brief Return the Bessel function of order \f$ \nu \f$:
440  * \f$ J_{\nu}(x) \f$.
441  *
442  * The cylindrical Bessel function is:
443  * @f[
444  * J_{\nu}(x) = \sum_{k=0}^{\infty}
445  * \frac{(-1)^k (x/2)^{\nu + 2k}}{k!\Gamma(\nu+k+1)}
446  * @f]
447  *
448  * @param __nu The order of the Bessel function.
449  * @param __x The argument of the Bessel function.
450  * @return The output Bessel function.
451  */
452  template<typename _Tp>
453  _Tp
454  __cyl_bessel_j(const _Tp __nu, const _Tp __x)
455  {
456  if (__nu < _Tp(0) || __x < _Tp(0))
457  std::__throw_domain_error(__N("Bad argument "
458  "in __cyl_bessel_j."));
459  else if (__isnan(__nu) || __isnan(__x))
461  else if (__x * __x < _Tp(10) * (__nu + _Tp(1)))
462  return __cyl_bessel_ij_series(__nu, __x, -_Tp(1), 200);
463  else if (__x > _Tp(1000))
464  {
465  _Tp __J_nu, __N_nu;
466  __cyl_bessel_jn_asymp(__nu, __x, __J_nu, __N_nu);
467  return __J_nu;
468  }
469  else
470  {
471  _Tp __J_nu, __N_nu, __Jp_nu, __Np_nu;
472  __bessel_jn(__nu, __x, __J_nu, __N_nu, __Jp_nu, __Np_nu);
473  return __J_nu;
474  }
475  }
476 
477 
478  /**
479  * @brief Return the Neumann function of order \f$ \nu \f$:
480  * \f$ N_{\nu}(x) \f$.
481  *
482  * The Neumann function is defined by:
483  * @f[
484  * N_{\nu}(x) = \frac{J_{\nu}(x) \cos \nu\pi - J_{-\nu}(x)}
485  * {\sin \nu\pi}
486  * @f]
487  * where for integral \f$ \nu = n \f$ a limit is taken:
488  * \f$ lim_{\nu \to n} \f$.
489  *
490  * @param __nu The order of the Neumann function.
491  * @param __x The argument of the Neumann function.
492  * @return The output Neumann function.
493  */
494  template<typename _Tp>
495  _Tp
496  __cyl_neumann_n(const _Tp __nu, const _Tp __x)
497  {
498  if (__nu < _Tp(0) || __x < _Tp(0))
499  std::__throw_domain_error(__N("Bad argument "
500  "in __cyl_neumann_n."));
501  else if (__isnan(__nu) || __isnan(__x))
503  else if (__x > _Tp(1000))
504  {
505  _Tp __J_nu, __N_nu;
506  __cyl_bessel_jn_asymp(__nu, __x, __J_nu, __N_nu);
507  return __N_nu;
508  }
509  else
510  {
511  _Tp __J_nu, __N_nu, __Jp_nu, __Np_nu;
512  __bessel_jn(__nu, __x, __J_nu, __N_nu, __Jp_nu, __Np_nu);
513  return __N_nu;
514  }
515  }
516 
517 
518  /**
519  * @brief Compute the spherical Bessel @f$ j_n(x) @f$
520  * and Neumann @f$ n_n(x) @f$ functions and their first
521  * derivatives @f$ j'_n(x) @f$ and @f$ n'_n(x) @f$
522  * respectively.
523  *
524  * @param __n The order of the spherical Bessel function.
525  * @param __x The argument of the spherical Bessel function.
526  * @param __j_n The output spherical Bessel function.
527  * @param __n_n The output spherical Neumann function.
528  * @param __jp_n The output derivative of the spherical Bessel function.
529  * @param __np_n The output derivative of the spherical Neumann function.
530  */
531  template <typename _Tp>
532  void
533  __sph_bessel_jn(const unsigned int __n, const _Tp __x,
534  _Tp & __j_n, _Tp & __n_n, _Tp & __jp_n, _Tp & __np_n)
535  {
536  const _Tp __nu = _Tp(__n) + _Tp(0.5L);
537 
538  _Tp __J_nu, __N_nu, __Jp_nu, __Np_nu;
539  __bessel_jn(__nu, __x, __J_nu, __N_nu, __Jp_nu, __Np_nu);
540 
541  const _Tp __factor = __numeric_constants<_Tp>::__sqrtpio2()
542  / std::sqrt(__x);
543 
544  __j_n = __factor * __J_nu;
545  __n_n = __factor * __N_nu;
546  __jp_n = __factor * __Jp_nu - __j_n / (_Tp(2) * __x);
547  __np_n = __factor * __Np_nu - __n_n / (_Tp(2) * __x);
548 
549  return;
550  }
551 
552 
553  /**
554  * @brief Return the spherical Bessel function
555  * @f$ j_n(x) @f$ of order n.
556  *
557  * The spherical Bessel function is defined by:
558  * @f[
559  * j_n(x) = \left( \frac{\pi}{2x} \right) ^{1/2} J_{n+1/2}(x)
560  * @f]
561  *
562  * @param __n The order of the spherical Bessel function.
563  * @param __x The argument of the spherical Bessel function.
564  * @return The output spherical Bessel function.
565  */
566  template <typename _Tp>
567  _Tp
568  __sph_bessel(const unsigned int __n, const _Tp __x)
569  {
570  if (__x < _Tp(0))
571  std::__throw_domain_error(__N("Bad argument "
572  "in __sph_bessel."));
573  else if (__isnan(__x))
575  else if (__x == _Tp(0))
576  {
577  if (__n == 0)
578  return _Tp(1);
579  else
580  return _Tp(0);
581  }
582  else
583  {
584  _Tp __j_n, __n_n, __jp_n, __np_n;
585  __sph_bessel_jn(__n, __x, __j_n, __n_n, __jp_n, __np_n);
586  return __j_n;
587  }
588  }
589 
590 
591  /**
592  * @brief Return the spherical Neumann function
593  * @f$ n_n(x) @f$.
594  *
595  * The spherical Neumann function is defined by:
596  * @f[
597  * n_n(x) = \left( \frac{\pi}{2x} \right) ^{1/2} N_{n+1/2}(x)
598  * @f]
599  *
600  * @param __n The order of the spherical Neumann function.
601  * @param __x The argument of the spherical Neumann function.
602  * @return The output spherical Neumann function.
603  */
604  template <typename _Tp>
605  _Tp
606  __sph_neumann(const unsigned int __n, const _Tp __x)
607  {
608  if (__x < _Tp(0))
609  std::__throw_domain_error(__N("Bad argument "
610  "in __sph_neumann."));
611  else if (__isnan(__x))
613  else if (__x == _Tp(0))
615  else
616  {
617  _Tp __j_n, __n_n, __jp_n, __np_n;
618  __sph_bessel_jn(__n, __x, __j_n, __n_n, __jp_n, __np_n);
619  return __n_n;
620  }
621  }
622 
623  } // namespace std::tr1::__detail
624 }
625 }
626 
627 #endif // _GLIBCXX_TR1_BESSEL_FUNCTION_TCC
_Tp __cyl_bessel_ij_series(const _Tp __nu, const _Tp __x, const _Tp __sgn, const unsigned int __max_iter)
This routine returns the cylindrical Bessel functions of order : or by series expansion.
complex< _Tp > log(const complex< _Tp > &)
Return complex natural logarithm of z.
Definition: complex:782
A structure for numeric constants.
void __cyl_bessel_jn_asymp(const _Tp __nu, const _Tp __x, _Tp &__Jnu, _Tp &__Nnu)
This routine computes the asymptotic cylindrical Bessel and Neumann functions of order nu: ...
void __sph_bessel_jn(const unsigned int __n, const _Tp __x, _Tp &__j_n, _Tp &__n_n, _Tp &__jp_n, _Tp &__np_n)
Compute the spherical Bessel and Neumann functions and their first derivatives and respectively...
complex< _Tp > cos(const complex< _Tp > &)
Return complex cosine of z.
Definition: complex:699
complex< _Tp > sqrt(const complex< _Tp > &)
Return complex square root of z.
Definition: complex:891
void __gamma_temme(const _Tp __mu, _Tp &__gam1, _Tp &__gam2, _Tp &__gampl, _Tp &__gammi)
Compute the gamma functions required by the Temme series expansions of and . and where is and ...
_Tp __cyl_neumann_n(const _Tp __nu, const _Tp __x)
Return the Neumann function of order : .
_Tp __cyl_bessel_j(const _Tp __nu, const _Tp __x)
Return the Bessel function of order : .
const _Tp & max(const _Tp &, const _Tp &)
This does what you think it does.
Definition: stl_algobase.h:209
_Tp __sph_bessel(const unsigned int __n, const _Tp __x)
Return the spherical Bessel function of order n.
Properties of fundamental types.
Definition: limits:278
complex< _Tp > sinh(const complex< _Tp > &)
Return complex hyperbolic sine of z.
Definition: complex:847
complex< _Tp > cosh(const complex< _Tp > &)
Return complex hyperbolic cosine of z.
Definition: complex:729
static _Tp epsilon()
Definition: limits:287
_Tp __log_gamma(const _Tp __x)
Return . This will return values even for . To recover the sign of for any argument use __log_gamma_...
Definition: gamma.tcc:221
complex< _Tp > sin(const complex< _Tp > &)
Return complex sine of z.
Definition: complex:817
void __bessel_jn(const _Tp __nu, const _Tp __x, _Tp &__Jnu, _Tp &__Nnu, _Tp &__Jpnu, _Tp &__Npnu)
Compute the Bessel and Neumann functions and their first derivatives and respectively. These four functions are computed together for numerical stability.
_Tp abs(const complex< _Tp > &)
Return magnitude of z.
Definition: complex:594
_Tp __sph_neumann(const unsigned int __n, const _Tp __x)
Return the spherical Neumann function .
_Tp __gamma(const _Tp __x)
Return .
Definition: gamma.tcc:333
static _Tp quiet_NaN()
Definition: limits:293
complex< _Tp > exp(const complex< _Tp > &)
Return complex base e exponential of z.
Definition: complex:755
static _Tp infinity()
Definition: limits:291