1 /* specfunc/gsl_sf_bessel.h
3 * Copyright (C) 1996, 1997, 1998, 1999, 2000 Gerard Jungman
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 3 of the License, or (at
8 * your option) any later version.
10 * This program is distributed in the hope that it will be useful, but
11 * WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 * General Public License for more details.
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
20 /* Author: G. Jungman */
22 #ifndef __GSL_SF_BESSEL_H__
23 #define __GSL_SF_BESSEL_H__
26 #include <gsl/gsl_mode.h>
27 #include <gsl/gsl_precision.h>
28 #include <gsl/gsl_sf_result.h>
33 # define __BEGIN_DECLS extern "C" {
34 # define __END_DECLS }
36 # define __BEGIN_DECLS /* empty */
37 # define __END_DECLS /* empty */
43 /* Regular Bessel Function J_0(x)
47 int gsl_sf_bessel_J0_e(const double x, gsl_sf_result * result);
48 double gsl_sf_bessel_J0(const double x);
51 /* Regular Bessel Function J_1(x)
53 * exceptions: GSL_EUNDRFLW
55 int gsl_sf_bessel_J1_e(const double x, gsl_sf_result * result);
56 double gsl_sf_bessel_J1(const double x);
59 /* Regular Bessel Function J_n(x)
61 * exceptions: GSL_EUNDRFLW
63 int gsl_sf_bessel_Jn_e(int n, double x, gsl_sf_result * result);
64 double gsl_sf_bessel_Jn(const int n, const double x);
67 /* Regular Bessel Function J_n(x), nmin <= n <= nmax
69 * exceptions: GSL_EDOM, GSL_EUNDRFLW
71 int gsl_sf_bessel_Jn_array(int nmin, int nmax, double x, double * result_array);
74 /* Irregular Bessel function Y_0(x)
77 * exceptions: GSL_EDOM, GSL_EUNDRFLW
79 int gsl_sf_bessel_Y0_e(const double x, gsl_sf_result * result);
80 double gsl_sf_bessel_Y0(const double x);
83 /* Irregular Bessel function Y_1(x)
86 * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
88 int gsl_sf_bessel_Y1_e(const double x, gsl_sf_result * result);
89 double gsl_sf_bessel_Y1(const double x);
92 /* Irregular Bessel function Y_n(x)
95 * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
97 int gsl_sf_bessel_Yn_e(int n,const double x, gsl_sf_result * result);
98 double gsl_sf_bessel_Yn(const int n,const double x);
101 /* Irregular Bessel function Y_n(x), nmin <= n <= nmax
104 * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
106 int gsl_sf_bessel_Yn_array(const int nmin, const int nmax, const double x, double * result_array);
109 /* Regular modified Bessel function I_0(x)
111 * exceptions: GSL_EOVRFLW
113 int gsl_sf_bessel_I0_e(const double x, gsl_sf_result * result);
114 double gsl_sf_bessel_I0(const double x);
117 /* Regular modified Bessel function I_1(x)
119 * exceptions: GSL_EOVRFLW, GSL_EUNDRFLW
121 int gsl_sf_bessel_I1_e(const double x, gsl_sf_result * result);
122 double gsl_sf_bessel_I1(const double x);
125 /* Regular modified Bessel function I_n(x)
127 * exceptions: GSL_EOVRFLW, GSL_EUNDRFLW
129 int gsl_sf_bessel_In_e(const int n, const double x, gsl_sf_result * result);
130 double gsl_sf_bessel_In(const int n, const double x);
133 /* Regular modified Bessel function I_n(x) for n=nmin,...,nmax
135 * nmin >=0, nmax >= nmin
136 * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
138 int gsl_sf_bessel_In_array(const int nmin, const int nmax, const double x, double * result_array);
141 /* Scaled regular modified Bessel function
146 int gsl_sf_bessel_I0_scaled_e(const double x, gsl_sf_result * result);
147 double gsl_sf_bessel_I0_scaled(const double x);
150 /* Scaled regular modified Bessel function
153 * exceptions: GSL_EUNDRFLW
155 int gsl_sf_bessel_I1_scaled_e(const double x, gsl_sf_result * result);
156 double gsl_sf_bessel_I1_scaled(const double x);
159 /* Scaled regular modified Bessel function
162 * exceptions: GSL_EUNDRFLW
164 int gsl_sf_bessel_In_scaled_e(int n, const double x, gsl_sf_result * result);
165 double gsl_sf_bessel_In_scaled(const int n, const double x);
168 /* Scaled regular modified Bessel function
169 * exp(-|x|) I_n(x) for n=nmin,...,nmax
171 * nmin >=0, nmax >= nmin
172 * exceptions: GSL_EUNDRFLW
174 int gsl_sf_bessel_In_scaled_array(const int nmin, const int nmax, const double x, double * result_array);
177 /* Irregular modified Bessel function K_0(x)
180 * exceptions: GSL_EDOM, GSL_EUNDRFLW
182 int gsl_sf_bessel_K0_e(const double x, gsl_sf_result * result);
183 double gsl_sf_bessel_K0(const double x);
186 /* Irregular modified Bessel function K_1(x)
189 * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
191 int gsl_sf_bessel_K1_e(const double x, gsl_sf_result * result);
192 double gsl_sf_bessel_K1(const double x);
195 /* Irregular modified Bessel function K_n(x)
198 * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
200 int gsl_sf_bessel_Kn_e(const int n, const double x, gsl_sf_result * result);
201 double gsl_sf_bessel_Kn(const int n, const double x);
204 /* Irregular modified Bessel function K_n(x) for n=nmin,...,nmax
206 * x > 0.0, nmin >=0, nmax >= nmin
207 * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
209 int gsl_sf_bessel_Kn_array(const int nmin, const int nmax, const double x, double * result_array);
212 /* Scaled irregular modified Bessel function
216 * exceptions: GSL_EDOM
218 int gsl_sf_bessel_K0_scaled_e(const double x, gsl_sf_result * result);
219 double gsl_sf_bessel_K0_scaled(const double x);
222 /* Scaled irregular modified Bessel function
226 * exceptions: GSL_EDOM, GSL_EUNDRFLW
228 int gsl_sf_bessel_K1_scaled_e(const double x, gsl_sf_result * result);
229 double gsl_sf_bessel_K1_scaled(const double x);
232 /* Scaled irregular modified Bessel function
236 * exceptions: GSL_EDOM, GSL_EUNDRFLW
238 int gsl_sf_bessel_Kn_scaled_e(int n, const double x, gsl_sf_result * result);
239 double gsl_sf_bessel_Kn_scaled(const int n, const double x);
242 /* Scaled irregular modified Bessel function exp(x) K_n(x) for n=nmin,...,nmax
244 * x > 0.0, nmin >=0, nmax >= nmin
245 * exceptions: GSL_EDOM, GSL_EUNDRFLW
247 int gsl_sf_bessel_Kn_scaled_array(const int nmin, const int nmax, const double x, double * result_array);
250 /* Regular spherical Bessel function j_0(x) = sin(x)/x
254 int gsl_sf_bessel_j0_e(const double x, gsl_sf_result * result);
255 double gsl_sf_bessel_j0(const double x);
258 /* Regular spherical Bessel function j_1(x) = (sin(x)/x - cos(x))/x
260 * exceptions: GSL_EUNDRFLW
262 int gsl_sf_bessel_j1_e(const double x, gsl_sf_result * result);
263 double gsl_sf_bessel_j1(const double x);
266 /* Regular spherical Bessel function j_2(x) = ((3/x^2 - 1)sin(x) - 3cos(x)/x)/x
268 * exceptions: GSL_EUNDRFLW
270 int gsl_sf_bessel_j2_e(const double x, gsl_sf_result * result);
271 double gsl_sf_bessel_j2(const double x);
274 /* Regular spherical Bessel function j_l(x)
277 * exceptions: GSL_EDOM, GSL_EUNDRFLW
279 int gsl_sf_bessel_jl_e(const int l, const double x, gsl_sf_result * result);
280 double gsl_sf_bessel_jl(const int l, const double x);
283 /* Regular spherical Bessel function j_l(x) for l=0,1,...,lmax
285 * exceptions: GSL_EDOM, GSL_EUNDRFLW
287 int gsl_sf_bessel_jl_array(const int lmax, const double x, double * result_array);
290 /* Regular spherical Bessel function j_l(x) for l=0,1,...,lmax
291 * Uses Steed's method.
293 * exceptions: GSL_EDOM, GSL_EUNDRFLW
295 int gsl_sf_bessel_jl_steed_array(const int lmax, const double x, double * jl_x_array);
298 /* Irregular spherical Bessel function y_0(x)
302 int gsl_sf_bessel_y0_e(const double x, gsl_sf_result * result);
303 double gsl_sf_bessel_y0(const double x);
306 /* Irregular spherical Bessel function y_1(x)
308 * exceptions: GSL_EUNDRFLW
310 int gsl_sf_bessel_y1_e(const double x, gsl_sf_result * result);
311 double gsl_sf_bessel_y1(const double x);
314 /* Irregular spherical Bessel function y_2(x)
316 * exceptions: GSL_EUNDRFLW
318 int gsl_sf_bessel_y2_e(const double x, gsl_sf_result * result);
319 double gsl_sf_bessel_y2(const double x);
322 /* Irregular spherical Bessel function y_l(x)
324 * exceptions: GSL_EUNDRFLW
326 int gsl_sf_bessel_yl_e(int l, const double x, gsl_sf_result * result);
327 double gsl_sf_bessel_yl(const int l, const double x);
330 /* Irregular spherical Bessel function y_l(x) for l=0,1,...,lmax
332 * exceptions: GSL_EUNDRFLW
334 int gsl_sf_bessel_yl_array(const int lmax, const double x, double * result_array);
337 /* Regular scaled modified spherical Bessel function
343 int gsl_sf_bessel_i0_scaled_e(const double x, gsl_sf_result * result);
344 double gsl_sf_bessel_i0_scaled(const double x);
347 /* Regular scaled modified spherical Bessel function
351 * exceptions: GSL_EUNDRFLW
353 int gsl_sf_bessel_i1_scaled_e(const double x, gsl_sf_result * result);
354 double gsl_sf_bessel_i1_scaled(const double x);
357 /* Regular scaled modified spherical Bessel function
361 * exceptions: GSL_EUNDRFLW
363 int gsl_sf_bessel_i2_scaled_e(const double x, gsl_sf_result * result);
364 double gsl_sf_bessel_i2_scaled(const double x);
367 /* Regular scaled modified spherical Bessel functions
371 * i_l(x) = Sqrt[Pi/(2x)] BesselI[l+1/2,x]
374 * exceptions: GSL_EDOM, GSL_EUNDRFLW
376 int gsl_sf_bessel_il_scaled_e(const int l, double x, gsl_sf_result * result);
377 double gsl_sf_bessel_il_scaled(const int l, const double x);
380 /* Regular scaled modified spherical Bessel functions
385 * exceptions: GSL_EUNDRFLW
387 int gsl_sf_bessel_il_scaled_array(const int lmax, const double x, double * result_array);
390 /* Irregular scaled modified spherical Bessel function
394 * exceptions: GSL_EDOM, GSL_EUNDRFLW
396 int gsl_sf_bessel_k0_scaled_e(const double x, gsl_sf_result * result);
397 double gsl_sf_bessel_k0_scaled(const double x);
400 /* Irregular modified spherical Bessel function
404 * exceptions: GSL_EDOM, GSL_EUNDRFLW, GSL_EOVRFLW
406 int gsl_sf_bessel_k1_scaled_e(const double x, gsl_sf_result * result);
407 double gsl_sf_bessel_k1_scaled(const double x);
410 /* Irregular modified spherical Bessel function
414 * exceptions: GSL_EDOM, GSL_EUNDRFLW, GSL_EOVRFLW
416 int gsl_sf_bessel_k2_scaled_e(const double x, gsl_sf_result * result);
417 double gsl_sf_bessel_k2_scaled(const double x);
420 /* Irregular modified spherical Bessel function
423 * k_l(x) = Sqrt[Pi/(2x)] BesselK[l+1/2,x]
425 * exceptions: GSL_EDOM, GSL_EUNDRFLW
427 int gsl_sf_bessel_kl_scaled_e(int l, const double x, gsl_sf_result * result);
428 double gsl_sf_bessel_kl_scaled(const int l, const double x);
431 /* Irregular scaled modified spherical Bessel function
435 * exceptions: GSL_EDOM, GSL_EUNDRFLW
437 int gsl_sf_bessel_kl_scaled_array(const int lmax, const double x, double * result_array);
440 /* Regular cylindrical Bessel function J_nu(x)
442 * exceptions: GSL_EDOM, GSL_EUNDRFLW
444 int gsl_sf_bessel_Jnu_e(const double nu, const double x, gsl_sf_result * result);
445 double gsl_sf_bessel_Jnu(const double nu, const double x);
448 /* Irregular cylindrical Bessel function Y_nu(x)
452 int gsl_sf_bessel_Ynu_e(double nu, double x, gsl_sf_result * result);
453 double gsl_sf_bessel_Ynu(const double nu, const double x);
456 /* Regular cylindrical Bessel function J_nu(x)
457 * evaluated at a series of x values. The array
458 * contains the x values. They are assumed to be
459 * strictly ordered and positive. The array is
460 * over-written with the values of J_nu(x_i).
462 * exceptions: GSL_EDOM, GSL_EINVAL
464 int gsl_sf_bessel_sequence_Jnu_e(double nu, gsl_mode_t mode, size_t size, double * v);
467 /* Scaled modified cylindrical Bessel functions
469 * Exp[-|x|] BesselI[nu, x]
472 * exceptions: GSL_EDOM
474 int gsl_sf_bessel_Inu_scaled_e(double nu, double x, gsl_sf_result * result);
475 double gsl_sf_bessel_Inu_scaled(double nu, double x);
478 /* Modified cylindrical Bessel functions
483 * exceptions: GSL_EDOM, GSL_EOVRFLW
485 int gsl_sf_bessel_Inu_e(double nu, double x, gsl_sf_result * result);
486 double gsl_sf_bessel_Inu(double nu, double x);
489 /* Scaled modified cylindrical Bessel functions
491 * Exp[+|x|] BesselK[nu, x]
494 * exceptions: GSL_EDOM
496 int gsl_sf_bessel_Knu_scaled_e(const double nu, const double x, gsl_sf_result * result);
497 double gsl_sf_bessel_Knu_scaled(const double nu, const double x);
500 /* Modified cylindrical Bessel functions
505 * exceptions: GSL_EDOM, GSL_EUNDRFLW
507 int gsl_sf_bessel_Knu_e(const double nu, const double x, gsl_sf_result * result);
508 double gsl_sf_bessel_Knu(const double nu, const double x);
511 /* Logarithm of modified cylindrical Bessel functions.
513 * Log[BesselK[nu, x]]
516 * exceptions: GSL_EDOM
518 int gsl_sf_bessel_lnKnu_e(const double nu, const double x, gsl_sf_result * result);
519 double gsl_sf_bessel_lnKnu(const double nu, const double x);
522 /* s'th positive zero of the Bessel function J_0(x).
526 int gsl_sf_bessel_zero_J0_e(unsigned int s, gsl_sf_result * result);
527 double gsl_sf_bessel_zero_J0(unsigned int s);
530 /* s'th positive zero of the Bessel function J_1(x).
534 int gsl_sf_bessel_zero_J1_e(unsigned int s, gsl_sf_result * result);
535 double gsl_sf_bessel_zero_J1(unsigned int s);
538 /* s'th positive zero of the Bessel function J_nu(x).
542 int gsl_sf_bessel_zero_Jnu_e(double nu, unsigned int s, gsl_sf_result * result);
543 double gsl_sf_bessel_zero_Jnu(double nu, unsigned int s);
548 #endif /* __GSL_SF_BESSEL_H__ */