123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549 |
- /* specfunc/gsl_sf_bessel.h
- *
- * Copyright (C) 1996, 1997, 1998, 1999, 2000 Gerard Jungman
- *
- * This program is free software; you can redistribute it and/or modify
- * it under the terms of the GNU General Public License as published by
- * the Free Software Foundation; either version 3 of the License, or (at
- * your option) any later version.
- *
- * This program is distributed in the hope that it will be useful, but
- * WITHOUT ANY WARRANTY; without even the implied warranty of
- * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
- * General Public License for more details.
- *
- * You should have received a copy of the GNU General Public License
- * along with this program; if not, write to the Free Software
- * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
- */
- /* Author: G. Jungman */
- #ifndef __GSL_SF_BESSEL_H__
- #define __GSL_SF_BESSEL_H__
- #include <stdlib.h>
- #include "gsl_mode.h"
- #include "gsl_precision.h"
- #include "gsl_sf_result.h"
- #undef __BEGIN_DECLS
- #undef __END_DECLS
- #ifdef __cplusplus
- # define __BEGIN_DECLS extern "C" {
- # define __END_DECLS }
- #else
- # define __BEGIN_DECLS /* empty */
- # define __END_DECLS /* empty */
- #endif
- __BEGIN_DECLS
- /* Regular Bessel Function J_0(x)
- *
- * exceptions: none
- */
- int gsl_sf_bessel_J0_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_J0(const double x);
- /* Regular Bessel Function J_1(x)
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_J1_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_J1(const double x);
- /* Regular Bessel Function J_n(x)
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_Jn_e(int n, double x, gsl_sf_result * result);
- double gsl_sf_bessel_Jn(const int n, const double x);
- /* Regular Bessel Function J_n(x), nmin <= n <= nmax
- *
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_Jn_array(int nmin, int nmax, double x, double * result_array);
- /* Irregular Bessel function Y_0(x)
- *
- * x > 0.0
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_Y0_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_Y0(const double x);
- /* Irregular Bessel function Y_1(x)
- *
- * x > 0.0
- * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_Y1_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_Y1(const double x);
- /* Irregular Bessel function Y_n(x)
- *
- * x > 0.0
- * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_Yn_e(int n,const double x, gsl_sf_result * result);
- double gsl_sf_bessel_Yn(const int n,const double x);
- /* Irregular Bessel function Y_n(x), nmin <= n <= nmax
- *
- * x > 0.0
- * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_Yn_array(const int nmin, const int nmax, const double x, double * result_array);
- /* Regular modified Bessel function I_0(x)
- *
- * exceptions: GSL_EOVRFLW
- */
- int gsl_sf_bessel_I0_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_I0(const double x);
- /* Regular modified Bessel function I_1(x)
- *
- * exceptions: GSL_EOVRFLW, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_I1_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_I1(const double x);
- /* Regular modified Bessel function I_n(x)
- *
- * exceptions: GSL_EOVRFLW, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_In_e(const int n, const double x, gsl_sf_result * result);
- double gsl_sf_bessel_In(const int n, const double x);
- /* Regular modified Bessel function I_n(x) for n=nmin,...,nmax
- *
- * nmin >=0, nmax >= nmin
- * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_In_array(const int nmin, const int nmax, const double x, double * result_array);
- /* Scaled regular modified Bessel function
- * exp(-|x|) I_0(x)
- *
- * exceptions: none
- */
- int gsl_sf_bessel_I0_scaled_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_I0_scaled(const double x);
- /* Scaled regular modified Bessel function
- * exp(-|x|) I_1(x)
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_I1_scaled_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_I1_scaled(const double x);
- /* Scaled regular modified Bessel function
- * exp(-|x|) I_n(x)
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_In_scaled_e(int n, const double x, gsl_sf_result * result);
- double gsl_sf_bessel_In_scaled(const int n, const double x);
- /* Scaled regular modified Bessel function
- * exp(-|x|) I_n(x) for n=nmin,...,nmax
- *
- * nmin >=0, nmax >= nmin
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_In_scaled_array(const int nmin, const int nmax, const double x, double * result_array);
- /* Irregular modified Bessel function K_0(x)
- *
- * x > 0.0
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_K0_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_K0(const double x);
- /* Irregular modified Bessel function K_1(x)
- *
- * x > 0.0
- * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_K1_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_K1(const double x);
- /* Irregular modified Bessel function K_n(x)
- *
- * x > 0.0
- * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_Kn_e(const int n, const double x, gsl_sf_result * result);
- double gsl_sf_bessel_Kn(const int n, const double x);
- /* Irregular modified Bessel function K_n(x) for n=nmin,...,nmax
- *
- * x > 0.0, nmin >=0, nmax >= nmin
- * exceptions: GSL_EDOM, GSL_EOVRFLW, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_Kn_array(const int nmin, const int nmax, const double x, double * result_array);
- /* Scaled irregular modified Bessel function
- * exp(x) K_0(x)
- *
- * x > 0.0
- * exceptions: GSL_EDOM
- */
- int gsl_sf_bessel_K0_scaled_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_K0_scaled(const double x);
- /* Scaled irregular modified Bessel function
- * exp(x) K_1(x)
- *
- * x > 0.0
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_K1_scaled_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_K1_scaled(const double x);
- /* Scaled irregular modified Bessel function
- * exp(x) K_n(x)
- *
- * x > 0.0
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_Kn_scaled_e(int n, const double x, gsl_sf_result * result);
- double gsl_sf_bessel_Kn_scaled(const int n, const double x);
- /* Scaled irregular modified Bessel function exp(x) K_n(x) for n=nmin,...,nmax
- *
- * x > 0.0, nmin >=0, nmax >= nmin
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_Kn_scaled_array(const int nmin, const int nmax, const double x, double * result_array);
- /* Regular spherical Bessel function j_0(x) = sin(x)/x
- *
- * exceptions: none
- */
- int gsl_sf_bessel_j0_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_j0(const double x);
- /* Regular spherical Bessel function j_1(x) = (sin(x)/x - cos(x))/x
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_j1_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_j1(const double x);
- /* Regular spherical Bessel function j_2(x) = ((3/x^2 - 1)sin(x) - 3cos(x)/x)/x
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_j2_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_j2(const double x);
- /* Regular spherical Bessel function j_l(x)
- *
- * l >= 0, x >= 0.0
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_jl_e(const int l, const double x, gsl_sf_result * result);
- double gsl_sf_bessel_jl(const int l, const double x);
- /* Regular spherical Bessel function j_l(x) for l=0,1,...,lmax
- *
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_jl_array(const int lmax, const double x, double * result_array);
- /* Regular spherical Bessel function j_l(x) for l=0,1,...,lmax
- * Uses Steed's method.
- *
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_jl_steed_array(const int lmax, const double x, double * jl_x_array);
- /* Irregular spherical Bessel function y_0(x)
- *
- * exceptions: none
- */
- int gsl_sf_bessel_y0_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_y0(const double x);
- /* Irregular spherical Bessel function y_1(x)
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_y1_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_y1(const double x);
- /* Irregular spherical Bessel function y_2(x)
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_y2_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_y2(const double x);
- /* Irregular spherical Bessel function y_l(x)
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_yl_e(int l, const double x, gsl_sf_result * result);
- double gsl_sf_bessel_yl(const int l, const double x);
- /* Irregular spherical Bessel function y_l(x) for l=0,1,...,lmax
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_yl_array(const int lmax, const double x, double * result_array);
- /* Regular scaled modified spherical Bessel function
- *
- * Exp[-|x|] i_0(x)
- *
- * exceptions: none
- */
- int gsl_sf_bessel_i0_scaled_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_i0_scaled(const double x);
- /* Regular scaled modified spherical Bessel function
- *
- * Exp[-|x|] i_1(x)
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_i1_scaled_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_i1_scaled(const double x);
- /* Regular scaled modified spherical Bessel function
- *
- * Exp[-|x|] i_2(x)
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_i2_scaled_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_i2_scaled(const double x);
- /* Regular scaled modified spherical Bessel functions
- *
- * Exp[-|x|] i_l(x)
- *
- * i_l(x) = Sqrt[Pi/(2x)] BesselI[l+1/2,x]
- *
- * l >= 0
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_il_scaled_e(const int l, double x, gsl_sf_result * result);
- double gsl_sf_bessel_il_scaled(const int l, const double x);
- /* Regular scaled modified spherical Bessel functions
- *
- * Exp[-|x|] i_l(x)
- * for l=0,1,...,lmax
- *
- * exceptions: GSL_EUNDRFLW
- */
- int gsl_sf_bessel_il_scaled_array(const int lmax, const double x, double * result_array);
- /* Irregular scaled modified spherical Bessel function
- * Exp[x] k_0(x)
- *
- * x > 0.0
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_k0_scaled_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_k0_scaled(const double x);
- /* Irregular modified spherical Bessel function
- * Exp[x] k_1(x)
- *
- * x > 0.0
- * exceptions: GSL_EDOM, GSL_EUNDRFLW, GSL_EOVRFLW
- */
- int gsl_sf_bessel_k1_scaled_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_k1_scaled(const double x);
- /* Irregular modified spherical Bessel function
- * Exp[x] k_2(x)
- *
- * x > 0.0
- * exceptions: GSL_EDOM, GSL_EUNDRFLW, GSL_EOVRFLW
- */
- int gsl_sf_bessel_k2_scaled_e(const double x, gsl_sf_result * result);
- double gsl_sf_bessel_k2_scaled(const double x);
- /* Irregular modified spherical Bessel function
- * Exp[x] k_l[x]
- *
- * k_l(x) = Sqrt[Pi/(2x)] BesselK[l+1/2,x]
- *
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_kl_scaled_e(int l, const double x, gsl_sf_result * result);
- double gsl_sf_bessel_kl_scaled(const int l, const double x);
- /* Irregular scaled modified spherical Bessel function
- * Exp[x] k_l(x)
- *
- * for l=0,1,...,lmax
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_kl_scaled_array(const int lmax, const double x, double * result_array);
- /* Regular cylindrical Bessel function J_nu(x)
- *
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_Jnu_e(const double nu, const double x, gsl_sf_result * result);
- double gsl_sf_bessel_Jnu(const double nu, const double x);
- /* Irregular cylindrical Bessel function Y_nu(x)
- *
- * exceptions:
- */
- int gsl_sf_bessel_Ynu_e(double nu, double x, gsl_sf_result * result);
- double gsl_sf_bessel_Ynu(const double nu, const double x);
- /* Regular cylindrical Bessel function J_nu(x)
- * evaluated at a series of x values. The array
- * contains the x values. They are assumed to be
- * strictly ordered and positive. The array is
- * over-written with the values of J_nu(x_i).
- *
- * exceptions: GSL_EDOM, GSL_EINVAL
- */
- int gsl_sf_bessel_sequence_Jnu_e(double nu, gsl_mode_t mode, size_t size, double * v);
- /* Scaled modified cylindrical Bessel functions
- *
- * Exp[-|x|] BesselI[nu, x]
- * x >= 0, nu >= 0
- *
- * exceptions: GSL_EDOM
- */
- int gsl_sf_bessel_Inu_scaled_e(double nu, double x, gsl_sf_result * result);
- double gsl_sf_bessel_Inu_scaled(double nu, double x);
- /* Modified cylindrical Bessel functions
- *
- * BesselI[nu, x]
- * x >= 0, nu >= 0
- *
- * exceptions: GSL_EDOM, GSL_EOVRFLW
- */
- int gsl_sf_bessel_Inu_e(double nu, double x, gsl_sf_result * result);
- double gsl_sf_bessel_Inu(double nu, double x);
- /* Scaled modified cylindrical Bessel functions
- *
- * Exp[+|x|] BesselK[nu, x]
- * x > 0, nu >= 0
- *
- * exceptions: GSL_EDOM
- */
- int gsl_sf_bessel_Knu_scaled_e(const double nu, const double x, gsl_sf_result * result);
- double gsl_sf_bessel_Knu_scaled(const double nu, const double x);
- /* Modified cylindrical Bessel functions
- *
- * BesselK[nu, x]
- * x > 0, nu >= 0
- *
- * exceptions: GSL_EDOM, GSL_EUNDRFLW
- */
- int gsl_sf_bessel_Knu_e(const double nu, const double x, gsl_sf_result * result);
- double gsl_sf_bessel_Knu(const double nu, const double x);
- /* Logarithm of modified cylindrical Bessel functions.
- *
- * Log[BesselK[nu, x]]
- * x > 0, nu >= 0
- *
- * exceptions: GSL_EDOM
- */
- int gsl_sf_bessel_lnKnu_e(const double nu, const double x, gsl_sf_result * result);
- double gsl_sf_bessel_lnKnu(const double nu, const double x);
- /* s'th positive zero of the Bessel function J_0(x).
- *
- * exceptions:
- */
- int gsl_sf_bessel_zero_J0_e(unsigned int s, gsl_sf_result * result);
- double gsl_sf_bessel_zero_J0(unsigned int s);
- /* s'th positive zero of the Bessel function J_1(x).
- *
- * exceptions:
- */
- int gsl_sf_bessel_zero_J1_e(unsigned int s, gsl_sf_result * result);
- double gsl_sf_bessel_zero_J1(unsigned int s);
- /* s'th positive zero of the Bessel function J_nu(x).
- *
- * exceptions:
- */
- int gsl_sf_bessel_zero_Jnu_e(double nu, unsigned int s, gsl_sf_result * result);
- double gsl_sf_bessel_zero_Jnu(double nu, unsigned int s);
- __END_DECLS
- #endif /* __GSL_SF_BESSEL_H__ */
|