2 [section:mbessel Modified Bessel Functions of the First and Second Kinds]
6 `#include <boost/math/special_functions/bessel.hpp>`
8 template <class T1, class T2>
9 ``__sf_result`` cyl_bessel_i(T1 v, T2 x);
11 template <class T1, class T2, class ``__Policy``>
12 ``__sf_result`` cyl_bessel_i(T1 v, T2 x, const ``__Policy``&);
14 template <class T1, class T2>
15 ``__sf_result`` cyl_bessel_k(T1 v, T2 x);
17 template <class T1, class T2, class ``__Policy``>
18 ``__sf_result`` cyl_bessel_k(T1 v, T2 x, const ``__Policy``&);
23 The functions __cyl_bessel_i and __cyl_bessel_k return the result of the
24 modified Bessel functions of the first and second kind respectively:
26 cyl_bessel_i(v, x) = I[sub v](x)
28 cyl_bessel_k(v, x) = K[sub v](x)
36 The return type of these functions is computed using the __arg_promotion_rules
37 when T1 and T2 are different types. The functions are also optimised for the
38 relatively common case that T1 is an integer.
42 The functions return the result of __domain_error whenever the result is
43 undefined or complex. For __cyl_bessel_j this occurs when `x < 0` and v is not
44 an integer, or when `x == 0` and `v != 0`. For __cyl_neumann this occurs
47 The following graph illustrates the exponential behaviour of I[sub v].
51 The following graph illustrates the exponential decay of K[sub v].
57 There are two sets of test values: spot values calculated using
58 [@http://functions.wolfram.com functions.wolfram.com],
59 and a much larger set of tests computed using
60 a simplified version of this implementation
61 (with all the special case handling removed).
65 The following tables show how the accuracy of these functions
66 varies on various platforms, along with comparison to other libraries.
67 Note that only results for the widest floating-point type on the
68 system are given, as narrower types have __zero_error. All values
69 are relative errors in units of epsilon. Note that our test suite
70 includes some fairly extreme inputs which results in most of the worst
71 problem cases in other libraries:
73 [table_cyl_bessel_i_integer_orders_]
77 [table_cyl_bessel_k_integer_orders_]
83 The following are handled as special cases first:
85 When computing I[sub v][space] for ['x < 0], then [nu][space] must be an integer
86 or a domain error occurs. If [nu][space] is an integer, then the function is
87 odd if [nu][space] is odd and even if [nu][space] is even, and we can reflect to
90 For I[sub v][space] with v equal to 0, 1 or 0.5 are handled as special cases.
92 The 0 and 1 cases use polynomial approximations on
93 finite and infinite intervals. The approximating forms
95 [@http://www.advanpix.com/2015/11/11/rational-approximations-for-the-modified-bessel-function-of-the-first-kind-i0-computations-double-precision/
96 "Rational Approximations for the Modified Bessel Function of the First Kind - I[sub 0](x) for Computations with Double Precision"]
97 by Pavel Holoborodko, extended by us to deal with up to 128-bit precision (with different approximations for each target precision).
117 The 0.5 case is a simple trigonometric function:
119 I[sub 0.5](x) = sqrt(2 / [pi]x) * sinh(x)
121 For K[sub v][space] with /v/ an integer, the result is calculated using the
126 starting from K[sub 0][space] and K[sub 1][space] which are calculated
127 using rational the approximations above. These rational approximations are
128 accurate to around 19 digits, and are therefore only used when T has
129 no more than 64 binary digits of precision.
131 When /x/ is small compared to /v/, I[sub v]x[space] is best computed directly from the series:
135 In the general case, we first normalize [nu][space] to \[[^0, [inf]])
136 with the help of the reflection formulae:
142 Let [mu][space] = [nu] - floor([nu] + 1/2), then [mu][space] is the fractional part of
143 [nu][space] such that |[mu]| <= 1/2 (we need this for convergence later). The idea is to
144 calculate K[sub [mu]](x) and K[sub [mu]+1](x), and use them to obtain
145 I[sub [nu]](x) and K[sub [nu]](x).
147 The algorithm is proposed by Temme in
148 N.M. Temme, ['On the numerical evaluation of the modified bessel function
149 of the third kind], Journal of Computational Physics, vol 19, 324 (1975),
150 which needs two continued fractions as well as the Wronskian:
158 The continued fractions are computed using the modified Lentz's method
159 (W.J. Lentz, ['Generating Bessel functions in Mie scattering calculations
160 using continued fractions], Applied Optics, vol 15, 668 (1976)).
161 Their convergence rates depend on ['x], therefore we need
162 different strategies for large ['x] and small ['x].
164 ['x > v], CF1 needs O(['x]) iterations to converge, CF2 converges rapidly.
166 ['x <= v], CF1 converges rapidly, CF2 fails to converge when ['x] [^->] 0.
168 When ['x] is large (['x] > 2), both continued fractions converge (CF1
169 may be slow for really large ['x]). K[sub [mu]][space] and K[sub [mu]+1][space]
178 ['S] is also a series that is summed along with CF2, see
179 I.J. Thompson and A.R. Barnett, ['Modified Bessel functions I_v and K_v
180 of real order and complex argument to selected accuracy], Computer Physics
181 Communications, vol 47, 245 (1987).
183 When ['x] is small (['x] <= 2), CF2 convergence may fail (but CF1
184 works very well). The solution here is Temme's series:
192 f[sub k][space] and h[sub k][space]
193 are also computed by recursions (involving gamma functions), but the
194 formulas are a little complicated, readers are referred to
195 N.M. Temme, ['On the numerical evaluation of the modified Bessel function
196 of the third kind], Journal of Computational Physics, vol 19, 324 (1975).
197 Note: Temme's series converge only for |[mu]| <= 1/2.
199 K[sub [nu]](x) is then calculated from the forward
200 recurrence, as is K[sub [nu]+1](x). With these two values and
201 f[sub [nu]], the Wronskian yields I[sub [nu]](x) directly.
206 Copyright 2006 John Maddock, Paul A. Bristow and Xiaogang Zhang.
207 Distributed under the Boost Software License, Version 1.0.
208 (See accompanying file LICENSE_1_0.txt or copy at
209 http://www.boost.org/LICENSE_1_0.txt).