1 /* specfunc/bessel.h 2 * 3 * Copyright (C) 1996, 1997, 1998, 1999, 2000 Gerard Jungman 4 * 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. 9 * 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. 14 * 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. 18 */ 19 20 /* Author: G. Jungman */ 21 22 #ifndef _BESSEL_H_ 23 #define _BESSEL_H_ 24 25 #include "gsl_sf_result.h" 26 27 28 /* Taylor expansion for J_nu(x) or I_nu(x) 29 * sign = -1 ==> Jnu 30 * sign = +1 ==> Inu 31 */ 32 int gsl_sf_bessel_IJ_taylor_e(const double nu, const double x, 33 const int sign, 34 const int kmax, 35 const double threshold, 36 gsl_sf_result * result 37 ); 38 39 int gsl_sf_bessel_Jnu_asympx_e(const double nu, const double x, gsl_sf_result * result); 40 int gsl_sf_bessel_Ynu_asympx_e(const double nu, const double x, gsl_sf_result * result); 41 42 int gsl_sf_bessel_Inu_scaled_asympx_e(const double nu, const double x, gsl_sf_result * result); 43 int gsl_sf_bessel_Knu_scaled_asympx_e(const double nu, const double x, gsl_sf_result * result); 44 45 int gsl_sf_bessel_Inu_scaled_asymp_unif_e(const double nu, const double x, gsl_sf_result * result); 46 int gsl_sf_bessel_Knu_scaled_asymp_unif_e(const double nu, const double x, gsl_sf_result * result); 47 48 49 /* ratio = J_{nu+1}(x) / J_nu(x) 50 * sgn = sgn(J_nu(x)) 51 */ 52 int 53 gsl_sf_bessel_J_CF1(const double nu, const double x, double * ratio, double * sgn); 54 55 56 /* ratio = I_{nu+1}(x) / I_nu(x) 57 */ 58 int 59 gsl_sf_bessel_I_CF1_ser(const double nu, const double x, double * ratio); 60 61 62 /* Evaluate the Steed method continued fraction CF2 for 63 * 64 * (J' + i Y')/(J + i Y) := P + i Q 65 */ 66 int 67 gsl_sf_bessel_JY_steed_CF2(const double nu, const double x, 68 double * P, double * Q); 69 70 71 int 72 gsl_sf_bessel_JY_mu_restricted(const double mu, const double x, 73 gsl_sf_result * Jmu, gsl_sf_result * Jmup1, 74 gsl_sf_result * Ymu, gsl_sf_result * Ymup1); 75 76 77 int 78 gsl_sf_bessel_K_scaled_steed_temme_CF2(const double nu, const double x, 79 double * K_nu, double * K_nup1, 80 double * Kp_nu); 81 82 83 /* These are of use in calculating the oscillating 84 * Bessel functions. 85 * cos(y - pi/4 + eps) 86 * sin(y - pi/4 + eps) 87 */ 88 int gsl_sf_bessel_cos_pi4_e(double y, double eps, gsl_sf_result * result); 89 int gsl_sf_bessel_sin_pi4_e(double y, double eps, gsl_sf_result * result); 90 91 92 #endif /* !_BESSEL_H_ */ 93