1 /* specfunc/bessel_J0.c
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 */
23 #include <gsl/gsl_math.h>
24 #include <gsl/gsl_errno.h>
25 #include <gsl/gsl_mode.h>
27 #include "bessel_amp_phase.h"
28 #include <gsl/gsl_sf_trig.h>
29 #include <gsl/gsl_sf_bessel.h>
31 #include "cheb_eval.c"
33 /*-*-*-*-*-*-*-*-*-*-*-* Private Section *-*-*-*-*-*-*-*-*-*-*-*/
36 /* based on SLATEC besj0, 1977 version, w. fullerton */
38 /* chebyshev expansions for Bessel functions
40 series for bj0 on the interval 0. to 1.60000d+01
41 with weighted error 7.47e-18
42 log weighted error 17.13
43 significant figures required 16.98
44 decimal places required 17.68
48 static double bj0_data[13] = {
50 -0.665223007764405132,
52 -0.0332527231700357697,
53 0.0023114179304694015,
54 -0.0000991127741995080,
55 0.0000028916708643998,
56 -0.0000000612108586630,
57 0.0000000009838650793,
58 -0.0000000000124235515,
59 0.0000000000001265433,
60 -0.0000000000000010619,
61 0.0000000000000000074,
63 static cheb_series bj0_cs = {
71 /*-*-*-*-*-*-*-*-*-*-*-* Functions with Error Codes *-*-*-*-*-*-*-*-*-*-*-*/
73 int gsl_sf_bessel_J0_e(const double x, gsl_sf_result * result)
77 /* CHECK_POINTER(result) */
79 if(y < 2.0*GSL_SQRT_DBL_EPSILON) {
85 return cheb_eval_e(&bj0_cs, 0.125*y*y - 1.0, result);
88 const double z = 32.0/(y*y) - 1.0;
92 const int stat_ca = cheb_eval_e(&_gsl_sf_bessel_amp_phase_bm0_cs, z, &ca);
93 const int stat_ct = cheb_eval_e(&_gsl_sf_bessel_amp_phase_bth0_cs, z, &ct);
94 const int stat_cp = gsl_sf_bessel_cos_pi4_e(y, ct.val/y, &cp);
95 const double sqrty = sqrt(y);
96 const double ampl = (0.75 + ca.val) / sqrty;
97 result->val = ampl * cp.val;
98 result->err = fabs(cp.val) * ca.err/sqrty + fabs(ampl) * cp.err;
99 result->err += GSL_DBL_EPSILON * fabs(result->val);
100 return GSL_ERROR_SELECT_3(stat_ca, stat_ct, stat_cp);
104 /*-*-*-*-*-*-*-*-*-* Functions w/ Natural Prototypes *-*-*-*-*-*-*-*-*-*-*/
108 double gsl_sf_bessel_J0(const double x)
110 EVAL_RESULT(gsl_sf_bessel_J0_e(x, &result));