Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 10:00:59

0001 /* specfunc/gsl_sf_bessel.h
0002  * 
0003  * Copyright (C) 1996, 1997, 1998, 1999, 2000 Gerard Jungman
0004  * 
0005  * This program is free software; you can redistribute it and/or modify
0006  * it under the terms of the GNU General Public License as published by
0007  * the Free Software Foundation; either version 3 of the License, or (at
0008  * your option) any later version.
0009  * 
0010  * This program is distributed in the hope that it will be useful, but
0011  * WITHOUT ANY WARRANTY; without even the implied warranty of
0012  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
0013  * General Public License for more details.
0014  * 
0015  * You should have received a copy of the GNU General Public License
0016  * along with this program; if not, write to the Free Software
0017  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
0018  */
0019 
0020 /* Author:  G. Jungman */
0021 
0022 #ifndef __GSL_SF_BESSEL_H__
0023 #define __GSL_SF_BESSEL_H__
0024 
0025 #include <stdlib.h>
0026 #include <gsl/gsl_mode.h>
0027 #include <gsl/gsl_precision.h>
0028 #include <gsl/gsl_sf_result.h>
0029 
0030 #undef __BEGIN_DECLS
0031 #undef __END_DECLS
0032 #ifdef __cplusplus
0033 # define __BEGIN_DECLS extern "C" {
0034 # define __END_DECLS }
0035 #else
0036 # define __BEGIN_DECLS /* empty */
0037 # define __END_DECLS /* empty */
0038 #endif
0039 
0040 __BEGIN_DECLS
0041 
0042 
0043 /* Regular Bessel Function J_0(x)
0044  *
0045  * exceptions: none
0046  */
0047 int gsl_sf_bessel_J0_e(const double x,  gsl_sf_result * result);
0048 double gsl_sf_bessel_J0(const double x);
0049 
0050 
0051 /* Regular Bessel Function J_1(x)
0052  *
0053  * exceptions: GSL_EUNDRFLW
0054  */
0055 int gsl_sf_bessel_J1_e(const double x, gsl_sf_result * result);
0056 double gsl_sf_bessel_J1(const double x);
0057 
0058 
0059 /* Regular Bessel Function J_n(x)
0060  *
0061  * exceptions: GSL_EUNDRFLW
0062  */
0063 int gsl_sf_bessel_Jn_e(int n, double x, gsl_sf_result * result);
0064 double gsl_sf_bessel_Jn(const int n, const double x);
0065 
0066 
0067 /* Regular Bessel Function J_n(x),  nmin <= n <= nmax
0068  *
0069  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0070  */
0071 int gsl_sf_bessel_Jn_array(int nmin, int nmax, double x, double * result_array);
0072 
0073 
0074 /* Irregular Bessel function Y_0(x)
0075  *
0076  * x > 0.0
0077  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0078  */
0079 int gsl_sf_bessel_Y0_e(const double x, gsl_sf_result * result);
0080 double gsl_sf_bessel_Y0(const double x);
0081 
0082 
0083 /* Irregular Bessel function Y_1(x)
0084  *
0085  * x > 0.0
0086  * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
0087  */
0088 int gsl_sf_bessel_Y1_e(const double x, gsl_sf_result * result);
0089 double gsl_sf_bessel_Y1(const double x);
0090 
0091 
0092 /* Irregular Bessel function Y_n(x)
0093  *
0094  * x > 0.0
0095  * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
0096  */
0097 int gsl_sf_bessel_Yn_e(int n,const double x, gsl_sf_result * result);
0098 double gsl_sf_bessel_Yn(const int n,const double x);
0099 
0100 
0101 /* Irregular Bessel function Y_n(x), nmin <= n <= nmax
0102  *
0103  * x > 0.0
0104  * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
0105  */
0106 int gsl_sf_bessel_Yn_array(const int nmin, const int nmax, const double x, double * result_array);
0107 
0108 
0109 /* Regular modified Bessel function I_0(x)
0110  *
0111  * exceptions: GSL_EOVRFLW
0112  */
0113 int gsl_sf_bessel_I0_e(const double x, gsl_sf_result * result);
0114 double gsl_sf_bessel_I0(const double x);
0115 
0116 
0117 /* Regular modified Bessel function I_1(x)
0118  *
0119  * exceptions: GSL_EOVRFLW, GSL_EUNDRFLW
0120  */
0121 int gsl_sf_bessel_I1_e(const double x, gsl_sf_result * result);
0122 double gsl_sf_bessel_I1(const double x);
0123 
0124 
0125 /* Regular modified Bessel function I_n(x)
0126  *
0127  * exceptions: GSL_EOVRFLW, GSL_EUNDRFLW
0128  */
0129 int gsl_sf_bessel_In_e(const int n, const double x, gsl_sf_result * result);
0130 double gsl_sf_bessel_In(const int n, const double x);
0131 
0132 
0133 /* Regular modified Bessel function  I_n(x) for n=nmin,...,nmax
0134  *
0135  * nmin >=0, nmax >= nmin
0136  * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
0137  */
0138 int gsl_sf_bessel_In_array(const int nmin, const int nmax, const double x, double * result_array);
0139 
0140 
0141 /* Scaled regular modified Bessel function
0142  *  exp(-|x|) I_0(x)
0143  *
0144  * exceptions: none
0145  */
0146 int gsl_sf_bessel_I0_scaled_e(const double x, gsl_sf_result * result);
0147 double gsl_sf_bessel_I0_scaled(const double x);
0148 
0149 
0150 /* Scaled regular modified Bessel function
0151  *  exp(-|x|) I_1(x)
0152  *
0153  * exceptions: GSL_EUNDRFLW
0154  */
0155 int gsl_sf_bessel_I1_scaled_e(const double x, gsl_sf_result * result);
0156 double gsl_sf_bessel_I1_scaled(const double x);
0157 
0158 
0159 /* Scaled regular modified Bessel function
0160  *  exp(-|x|) I_n(x)
0161  *
0162  * exceptions: GSL_EUNDRFLW
0163  */
0164 int gsl_sf_bessel_In_scaled_e(int n, const double x, gsl_sf_result * result);
0165 double gsl_sf_bessel_In_scaled(const int n, const double x);
0166 
0167 
0168 /* Scaled regular modified Bessel function
0169  *  exp(-|x|) I_n(x)  for n=nmin,...,nmax
0170  *
0171  * nmin >=0, nmax >= nmin
0172  * exceptions: GSL_EUNDRFLW
0173  */
0174 int gsl_sf_bessel_In_scaled_array(const int nmin, const int nmax, const double x, double * result_array);
0175 
0176 
0177 /* Irregular modified Bessel function K_0(x)
0178  *
0179  * x > 0.0
0180  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0181  */
0182 int gsl_sf_bessel_K0_e(const double x, gsl_sf_result * result);
0183 double gsl_sf_bessel_K0(const double x);
0184 
0185 
0186 /* Irregular modified Bessel function K_1(x)
0187  *
0188  * x > 0.0
0189  * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
0190  */
0191 int gsl_sf_bessel_K1_e(const double x, gsl_sf_result * result);
0192 double gsl_sf_bessel_K1(const double x);
0193 
0194 
0195 /* Irregular modified Bessel function K_n(x)
0196  *
0197  * x > 0.0
0198  * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
0199  */
0200 int gsl_sf_bessel_Kn_e(const int n, const double x, gsl_sf_result * result);
0201 double gsl_sf_bessel_Kn(const int n, const double x);
0202 
0203 
0204 /* Irregular modified Bessel function  K_n(x)  for n=nmin,...,nmax
0205  *
0206  * x > 0.0, nmin >=0, nmax >= nmin
0207  * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
0208  */
0209 int gsl_sf_bessel_Kn_array(const int nmin, const int nmax, const double x, double * result_array);
0210 
0211 
0212 /* Scaled irregular modified Bessel function
0213  *  exp(x) K_0(x)
0214  *
0215  * x > 0.0
0216  * exceptions: GSL_EDOM
0217  */
0218 int gsl_sf_bessel_K0_scaled_e(const double x, gsl_sf_result * result);
0219 double gsl_sf_bessel_K0_scaled(const double x);
0220 
0221 
0222 /* Scaled irregular modified Bessel function
0223  *  exp(x) K_1(x)
0224  *
0225  * x > 0.0
0226  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0227  */
0228 int gsl_sf_bessel_K1_scaled_e(const double x, gsl_sf_result * result); 
0229 double gsl_sf_bessel_K1_scaled(const double x);
0230 
0231 
0232 /* Scaled irregular modified Bessel function
0233  *  exp(x) K_n(x)
0234  *
0235  * x > 0.0
0236  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0237  */
0238 int gsl_sf_bessel_Kn_scaled_e(int n, const double x, gsl_sf_result * result);
0239 double gsl_sf_bessel_Kn_scaled(const int n, const double x);
0240 
0241 
0242 /* Scaled irregular modified Bessel function  exp(x) K_n(x)  for n=nmin,...,nmax
0243  *
0244  * x > 0.0, nmin >=0, nmax >= nmin
0245  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0246  */
0247 int gsl_sf_bessel_Kn_scaled_array(const int nmin, const int nmax, const double x, double * result_array);
0248 
0249 
0250 /* Regular spherical Bessel function j_0(x) = sin(x)/x
0251  *
0252  * exceptions: none
0253  */
0254 int gsl_sf_bessel_j0_e(const double x, gsl_sf_result * result);
0255 double gsl_sf_bessel_j0(const double x);
0256 
0257 
0258 /* Regular spherical Bessel function j_1(x) = (sin(x)/x - cos(x))/x
0259  *
0260  * exceptions: GSL_EUNDRFLW
0261  */
0262 int gsl_sf_bessel_j1_e(const double x, gsl_sf_result * result);
0263 double gsl_sf_bessel_j1(const double x);
0264 
0265 
0266 /* Regular spherical Bessel function j_2(x) = ((3/x^2 - 1)sin(x) - 3cos(x)/x)/x
0267  *
0268  * exceptions: GSL_EUNDRFLW
0269  */
0270 int gsl_sf_bessel_j2_e(const double x, gsl_sf_result * result);
0271 double gsl_sf_bessel_j2(const double x);
0272 
0273 
0274 /* Regular spherical Bessel function j_l(x)
0275  *
0276  * l >= 0, x >= 0.0
0277  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0278  */
0279 int gsl_sf_bessel_jl_e(const int l, const double x, gsl_sf_result * result);
0280 double gsl_sf_bessel_jl(const int l, const double x);
0281 
0282 
0283 /* Regular spherical Bessel function j_l(x) for l=0,1,...,lmax
0284  *
0285  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0286  */
0287 int gsl_sf_bessel_jl_array(const int lmax, const double x, double * result_array);
0288 
0289 
0290 /* Regular spherical Bessel function j_l(x) for l=0,1,...,lmax
0291  * Uses Steed's method.
0292  *
0293  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0294  */
0295 int gsl_sf_bessel_jl_steed_array(const int lmax, const double x, double * jl_x_array);
0296 
0297 
0298 /* Irregular spherical Bessel function y_0(x)
0299  *
0300  * exceptions: none
0301  */
0302 int gsl_sf_bessel_y0_e(const double x, gsl_sf_result * result);
0303 double gsl_sf_bessel_y0(const double x);
0304 
0305 
0306 /* Irregular spherical Bessel function y_1(x)
0307  *
0308  * exceptions: GSL_EUNDRFLW
0309  */
0310 int gsl_sf_bessel_y1_e(const double x, gsl_sf_result * result);
0311 double gsl_sf_bessel_y1(const double x);
0312 
0313 
0314 /* Irregular spherical Bessel function y_2(x)
0315  *
0316  * exceptions: GSL_EUNDRFLW
0317  */
0318 int gsl_sf_bessel_y2_e(const double x, gsl_sf_result * result);
0319 double gsl_sf_bessel_y2(const double x);
0320 
0321 
0322 /* Irregular spherical Bessel function y_l(x)
0323  *
0324  * exceptions: GSL_EUNDRFLW
0325  */
0326 int gsl_sf_bessel_yl_e(int l, const double x, gsl_sf_result * result);
0327 double gsl_sf_bessel_yl(const int l, const double x);
0328 
0329 
0330 /* Irregular spherical Bessel function y_l(x) for l=0,1,...,lmax
0331  *
0332  * exceptions: GSL_EUNDRFLW
0333  */
0334 int gsl_sf_bessel_yl_array(const int lmax, const double x, double * result_array);
0335 
0336 
0337 /* Regular scaled modified spherical Bessel function
0338  *
0339  * Exp[-|x|] i_0(x)
0340  *
0341  * exceptions: none
0342  */
0343 int gsl_sf_bessel_i0_scaled_e(const double x, gsl_sf_result * result);
0344 double gsl_sf_bessel_i0_scaled(const double x);
0345 
0346 
0347 /* Regular scaled modified spherical Bessel function
0348  *
0349  * Exp[-|x|] i_1(x)
0350  *
0351  * exceptions: GSL_EUNDRFLW
0352  */
0353 int gsl_sf_bessel_i1_scaled_e(const double x, gsl_sf_result * result);
0354 double gsl_sf_bessel_i1_scaled(const double x);
0355 
0356 
0357 /* Regular scaled modified spherical Bessel function
0358  *
0359  * Exp[-|x|] i_2(x)
0360  *
0361  * exceptions: GSL_EUNDRFLW
0362  */
0363 int gsl_sf_bessel_i2_scaled_e(const double x, gsl_sf_result * result);
0364 double gsl_sf_bessel_i2_scaled(const double x);
0365 
0366 
0367 /* Regular scaled modified spherical Bessel functions
0368  *
0369  * Exp[-|x|] i_l(x)
0370  *
0371  * i_l(x) = Sqrt[Pi/(2x)] BesselI[l+1/2,x]
0372  *
0373  * l >= 0
0374  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0375  */
0376 int gsl_sf_bessel_il_scaled_e(const int l, double x, gsl_sf_result * result);
0377 double gsl_sf_bessel_il_scaled(const int l, const double x);
0378 
0379 
0380 /* Regular scaled modified spherical Bessel functions
0381  *
0382  * Exp[-|x|] i_l(x)
0383  * for l=0,1,...,lmax
0384  *
0385  * exceptions: GSL_EUNDRFLW
0386  */
0387 int gsl_sf_bessel_il_scaled_array(const int lmax, const double x, double * result_array);
0388 
0389 
0390 /* Irregular scaled modified spherical Bessel function
0391  * Exp[x] k_0(x)
0392  *
0393  * x > 0.0
0394  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0395  */
0396 int gsl_sf_bessel_k0_scaled_e(const double x, gsl_sf_result * result);
0397 double gsl_sf_bessel_k0_scaled(const double x);
0398 
0399 
0400 /* Irregular modified spherical Bessel function
0401  * Exp[x] k_1(x)
0402  *
0403  * x > 0.0
0404  * exceptions: GSL_EDOM, GSL_EUNDRFLW, GSL_EOVRFLW
0405  */
0406 int gsl_sf_bessel_k1_scaled_e(const double x, gsl_sf_result * result);
0407 double gsl_sf_bessel_k1_scaled(const double x);
0408 
0409 
0410 /* Irregular modified spherical Bessel function
0411  * Exp[x] k_2(x)
0412  *
0413  * x > 0.0
0414  * exceptions: GSL_EDOM, GSL_EUNDRFLW, GSL_EOVRFLW
0415  */
0416 int gsl_sf_bessel_k2_scaled_e(const double x, gsl_sf_result * result);
0417 double gsl_sf_bessel_k2_scaled(const double x);
0418 
0419 
0420 /* Irregular modified spherical Bessel function
0421  * Exp[x] k_l[x]
0422  *
0423  * k_l(x) = Sqrt[Pi/(2x)] BesselK[l+1/2,x]
0424  *
0425  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0426  */
0427 int gsl_sf_bessel_kl_scaled_e(int l, const double x, gsl_sf_result * result);
0428 double gsl_sf_bessel_kl_scaled(const int l, const double x);
0429 
0430 
0431 /* Irregular scaled modified spherical Bessel function
0432  * Exp[x] k_l(x)
0433  *
0434  * for l=0,1,...,lmax
0435  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0436  */
0437 int gsl_sf_bessel_kl_scaled_array(const int lmax, const double x, double * result_array);
0438 
0439 
0440 /* Regular cylindrical Bessel function J_nu(x)
0441  *
0442  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0443  */
0444 int gsl_sf_bessel_Jnu_e(const double nu, const double x, gsl_sf_result * result);
0445 double gsl_sf_bessel_Jnu(const double nu, const double x);
0446 
0447 
0448 /* Irregular cylindrical Bessel function Y_nu(x)
0449  *
0450  * exceptions:  
0451  */
0452 int gsl_sf_bessel_Ynu_e(double nu, double x, gsl_sf_result * result);
0453 double gsl_sf_bessel_Ynu(const double nu, const double x);
0454 
0455 
0456 /* Regular cylindrical Bessel function J_nu(x)
0457  * evaluated at a series of x values. The array
0458  * contains the x values. They are assumed to be
0459  * strictly ordered and positive. The array is
0460  * over-written with the values of J_nu(x_i).
0461  *
0462  * exceptions: GSL_EDOM, GSL_EINVAL
0463  */
0464 int gsl_sf_bessel_sequence_Jnu_e(double nu, gsl_mode_t mode, size_t size, double * v);
0465 
0466 
0467 /* Scaled modified cylindrical Bessel functions
0468  *
0469  * Exp[-|x|] BesselI[nu, x]
0470  * x >= 0, nu >= 0
0471  *
0472  * exceptions: GSL_EDOM
0473  */
0474 int gsl_sf_bessel_Inu_scaled_e(double nu, double x, gsl_sf_result * result);
0475 double gsl_sf_bessel_Inu_scaled(double nu, double x);
0476 
0477 
0478 /* Modified cylindrical Bessel functions
0479  *
0480  * BesselI[nu, x]
0481  * x >= 0, nu >= 0
0482  *
0483  * exceptions: GSL_EDOM, GSL_EOVRFLW
0484  */
0485 int gsl_sf_bessel_Inu_e(double nu, double x, gsl_sf_result * result);
0486 double gsl_sf_bessel_Inu(double nu, double x);
0487 
0488 
0489 /* Scaled modified cylindrical Bessel functions
0490  *
0491  * Exp[+|x|] BesselK[nu, x]
0492  * x > 0, nu >= 0
0493  *
0494  * exceptions: GSL_EDOM
0495  */
0496 int gsl_sf_bessel_Knu_scaled_e(const double nu, const double x, gsl_sf_result * result);
0497 double gsl_sf_bessel_Knu_scaled(const double nu, const double x);
0498 
0499 int gsl_sf_bessel_Knu_scaled_e10_e(const double nu, const double x, gsl_sf_result_e10 * result);
0500 
0501 /* Modified cylindrical Bessel functions
0502  *
0503  * BesselK[nu, x]
0504  * x > 0, nu >= 0
0505  *
0506  * exceptions: GSL_EDOM, GSL_EUNDRFLW
0507  */
0508 int gsl_sf_bessel_Knu_e(const double nu, const double x, gsl_sf_result * result);
0509 double gsl_sf_bessel_Knu(const double nu, const double x);
0510 
0511 
0512 /* Logarithm of modified cylindrical Bessel functions.
0513  *
0514  * Log[BesselK[nu, x]]
0515  * x > 0, nu >= 0
0516  *
0517  * exceptions: GSL_EDOM
0518  */
0519 int gsl_sf_bessel_lnKnu_e(const double nu, const double x, gsl_sf_result * result);
0520 double gsl_sf_bessel_lnKnu(const double nu, const double x);
0521 
0522 
0523 /* s'th positive zero of the Bessel function J_0(x).
0524  *
0525  * exceptions: 
0526  */
0527 int gsl_sf_bessel_zero_J0_e(unsigned int s, gsl_sf_result * result);
0528 double gsl_sf_bessel_zero_J0(unsigned int s);
0529 
0530 
0531 /* s'th positive zero of the Bessel function J_1(x).
0532  *
0533  * exceptions: 
0534  */
0535 int gsl_sf_bessel_zero_J1_e(unsigned int s, gsl_sf_result * result);
0536 double gsl_sf_bessel_zero_J1(unsigned int s);
0537 
0538 
0539 /* s'th positive zero of the Bessel function J_nu(x).
0540  *
0541  * exceptions: 
0542  */
0543 int gsl_sf_bessel_zero_Jnu_e(double nu, unsigned int s, gsl_sf_result * result);
0544 double gsl_sf_bessel_zero_Jnu(double nu, unsigned int s);
0545 
0546 
0547 __END_DECLS
0548 
0549 #endif /* __GSL_SF_BESSEL_H__ */