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_sf_pow_int.h>
26 #include <gsl/gsl_sf_gamma.h>
27 #include <gsl/gsl_sf_trig.h>
28 #include <gsl/gsl_sf_bessel.h>
33 #include "bessel_olver.h"
35 /*-*-*-*-*-*-*-*-*-*-*-* Private Section *-*-*-*-*-*-*-*-*-*-*-*/
37 /* [Abramowitz+Stegun, 10.1.3]
38 * with lmax=15, precision ~ 15D for x < 3
40 * checked OK [GJ] Wed May 13 15:41:25 MDT 1998
42 static int bessel_yl_small_x(int l, const double x, gsl_sf_result * result)
44 gsl_sf_result num_fact;
45 double den = gsl_sf_pow_int(x, l+1);
46 int stat_df = gsl_sf_doublefact_e(2*l-1, &num_fact);
48 if(stat_df != GSL_SUCCESS || den == 0.0) {
49 OVERFLOW_ERROR(result);
59 for(i=1; i<=lmax; i++) {
60 t_coeff /= i*(2*(i-l) - 1);
62 delta = t_power*t_coeff;
64 if(fabs(delta/sum) < 0.5*GSL_DBL_EPSILON) break;
66 result->val = -num_fact.val/den * sum;
67 result->err = GSL_DBL_EPSILON * fabs(result->val);
74 /*-*-*-*-*-*-*-*-*-*-*-* Functions with Error Codes *-*-*-*-*-*-*-*-*-*-*-*/
77 int gsl_sf_bessel_y0_e(const double x, gsl_sf_result * result)
79 /* CHECK_POINTER(result) */
84 else if(1.0/GSL_DBL_MAX > 0.0 && x < 1.0/GSL_DBL_MAX) {
85 OVERFLOW_ERROR(result);
88 gsl_sf_result cos_result;
89 const int stat = gsl_sf_cos_e(x, &cos_result);
90 result->val = -cos_result.val/x;
91 result->err = fabs(cos_result.err/x);
92 result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
98 int gsl_sf_bessel_y1_e(const double x, gsl_sf_result * result)
100 /* CHECK_POINTER(result) */
103 DOMAIN_ERROR(result);
105 else if(x < 1.0/GSL_SQRT_DBL_MAX) {
106 OVERFLOW_ERROR(result);
109 const double y = x*x;
110 const double c1 = 1.0/2.0;
111 const double c2 = -1.0/8.0;
112 const double c3 = 1.0/144.0;
113 const double c4 = -1.0/5760.0;
114 const double c5 = 1.0/403200.0;
115 const double c6 = -1.0/43545600.0;
116 const double sum = 1.0 + y*(c1 + y*(c2 + y*(c3 + y*(c4 + y*(c5 + y*c6)))));
117 result->val = -sum/y;
118 result->err = GSL_DBL_EPSILON * fabs(result->val);
122 gsl_sf_result cos_result;
123 gsl_sf_result sin_result;
124 const int stat_cos = gsl_sf_cos_e(x, &cos_result);
125 const int stat_sin = gsl_sf_sin_e(x, &sin_result);
126 const double cx = cos_result.val;
127 const double sx = sin_result.val;
128 result->val = -(cx/x + sx)/x;
129 result->err = (fabs(cos_result.err/x) + sin_result.err)/fabs(x);
130 result->err += GSL_DBL_EPSILON * (fabs(sx/x) + fabs(cx/(x*x)));
131 return GSL_ERROR_SELECT_2(stat_cos, stat_sin);
136 int gsl_sf_bessel_y2_e(const double x, gsl_sf_result * result)
138 /* CHECK_POINTER(result) */
141 DOMAIN_ERROR(result);
143 else if(x < 1.0/GSL_ROOT3_DBL_MAX) {
144 OVERFLOW_ERROR(result);
147 const double y = x*x;
148 const double c1 = 1.0/6.0;
149 const double c2 = 1.0/24.0;
150 const double c3 = -1.0/144.0;
151 const double c4 = 1.0/3456.0;
152 const double c5 = -1.0/172800.0;
153 const double c6 = 1.0/14515200.0;
154 const double c7 = -1.0/1828915200.0;
155 const double sum = 1.0 + y*(c1 + y*(c2 + y*(c3 + y*(c4 + y*(c5 + y*(c6 + y*c7))))));
156 result->val = -3.0/(x*x*x) * sum;
157 result->err = GSL_DBL_EPSILON * fabs(result->val);
161 gsl_sf_result cos_result;
162 gsl_sf_result sin_result;
163 const int stat_cos = gsl_sf_cos_e(x, &cos_result);
164 const int stat_sin = gsl_sf_sin_e(x, &sin_result);
165 const double sx = sin_result.val;
166 const double cx = cos_result.val;
167 const double a = 3.0/(x*x);
168 result->val = (1.0 - a)/x * cx - a * sx;
169 result->err = cos_result.err * fabs((1.0 - a)/x) + sin_result.err * fabs(a);
170 result->err += GSL_DBL_EPSILON * (fabs(cx/x) + fabs(sx/(x*x)));
171 return GSL_ERROR_SELECT_2(stat_cos, stat_sin);
176 int gsl_sf_bessel_yl_e(int l, const double x, gsl_sf_result * result)
178 /* CHECK_POINTER(result) */
180 if(l < 0 || x <= 0.0) {
181 DOMAIN_ERROR(result);
184 return gsl_sf_bessel_y0_e(x, result);
187 return gsl_sf_bessel_y1_e(x, result);
190 return gsl_sf_bessel_y2_e(x, result);
193 return bessel_yl_small_x(l, x, result);
195 else if(GSL_ROOT3_DBL_EPSILON * x > (l*l + l + 1.0)) {
196 int status = gsl_sf_bessel_Ynu_asympx_e(l + 0.5, x, result);
197 double pre = sqrt((0.5*M_PI)/x);
203 int status = gsl_sf_bessel_Ynu_asymp_Olver_e(l + 0.5, x, result);
204 double pre = sqrt((0.5*M_PI)/x);
213 int stat_1 = gsl_sf_bessel_y1_e(x, &r_by);
214 int stat_0 = gsl_sf_bessel_y0_e(x, &r_bym);
215 double bym = r_bym.val;
216 double by = r_by.val;
220 byp = (2*j+1)/x*by - bym;
225 result->err = fabs(result->val) * (GSL_DBL_EPSILON + fabs(r_by.err/r_by.val) + fabs(r_bym.err/r_bym.val));
227 return GSL_ERROR_SELECT_2(stat_1, stat_0);
232 int gsl_sf_bessel_yl_array(const int lmax, const double x, double * result_array)
234 /* CHECK_POINTER(result_array) */
236 if(lmax < 0 || x <= 0.0) {
237 GSL_ERROR ("error", GSL_EDOM);
238 } else if (lmax == 0) {
239 gsl_sf_result result;
240 int stat = gsl_sf_bessel_y0_e(x, &result);
241 result_array[0] = result.val;
244 gsl_sf_result r_yell;
245 gsl_sf_result r_yellm1;
246 int stat_1 = gsl_sf_bessel_y1_e(x, &r_yell);
247 int stat_0 = gsl_sf_bessel_y0_e(x, &r_yellm1);
249 double yell = r_yell.val;
250 double yellm1 = r_yellm1.val;
253 result_array[0] = yellm1;
254 result_array[1] = yell;
256 for(ell = 1; ell < lmax; ell++) {
257 yellp1 = (2*ell+1)/x * yell - yellm1;
258 result_array[ell+1] = yellp1;
263 return GSL_ERROR_SELECT_2(stat_0, stat_1);
268 /*-*-*-*-*-*-*-*-*-* Functions w/ Natural Prototypes *-*-*-*-*-*-*-*-*-*-*/
272 double gsl_sf_bessel_y0(const double x)
274 EVAL_RESULT(gsl_sf_bessel_y0_e(x, &result));
277 double gsl_sf_bessel_y1(const double x)
279 EVAL_RESULT(gsl_sf_bessel_y1_e(x, &result));
282 double gsl_sf_bessel_y2(const double x)
284 EVAL_RESULT(gsl_sf_bessel_y2_e(x, &result));
287 double gsl_sf_bessel_yl(const int l, const double x)
289 EVAL_RESULT(gsl_sf_bessel_yl_e(l, x, &result));