Blame specfunc/bessel_j.c

Packit 67cb25
/* specfunc/bessel_j.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 1996,1997,1998,1999,2000,2001,2002,2003 Gerard Jungman
Packit 67cb25
 * 
Packit 67cb25
 * This program is free software; you can redistribute it and/or modify
Packit 67cb25
 * it under the terms of the GNU General Public License as published by
Packit 67cb25
 * the Free Software Foundation; either version 3 of the License, or (at
Packit 67cb25
 * your option) any later version.
Packit 67cb25
 * 
Packit 67cb25
 * This program is distributed in the hope that it will be useful, but
Packit 67cb25
 * WITHOUT ANY WARRANTY; without even the implied warranty of
Packit 67cb25
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
Packit 67cb25
 * General Public License for more details.
Packit 67cb25
 * 
Packit 67cb25
 * You should have received a copy of the GNU General Public License
Packit 67cb25
 * along with this program; if not, write to the Free Software
Packit 67cb25
 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
/* Author:  G. Jungman */
Packit 67cb25
Packit 67cb25
#include <config.h>
Packit 67cb25
#include <gsl/gsl_math.h>
Packit 67cb25
#include <gsl/gsl_errno.h>
Packit 67cb25
#include <gsl/gsl_sf_pow_int.h>
Packit 67cb25
#include <gsl/gsl_sf_trig.h>
Packit 67cb25
#include <gsl/gsl_sf_bessel.h>
Packit 67cb25
Packit 67cb25
#include "error.h"
Packit 67cb25
Packit 67cb25
#include "bessel.h"
Packit 67cb25
#include "bessel_olver.h"
Packit 67cb25
Packit 67cb25
/*-*-*-*-*-*-*-*-*-*-*-* Functions with Error Codes *-*-*-*-*-*-*-*-*-*-*-*/
Packit 67cb25
Packit 67cb25
int gsl_sf_bessel_j0_e(const double x, gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  double ax = fabs(x);
Packit 67cb25
Packit 67cb25
  /* CHECK_POINTER(result) */
Packit 67cb25
Packit 67cb25
  if(ax < 0.5) {
Packit 67cb25
    const double y = x*x;
Packit 67cb25
    const double c1 = -1.0/6.0;
Packit 67cb25
    const double c2 =  1.0/120.0;
Packit 67cb25
    const double c3 = -1.0/5040.0;
Packit 67cb25
    const double c4 =  1.0/362880.0;
Packit 67cb25
    const double c5 = -1.0/39916800.0;
Packit 67cb25
    const double c6 =  1.0/6227020800.0;
Packit 67cb25
    result->val = 1.0 + y*(c1 + y*(c2 + y*(c3 + y*(c4 + y*(c5 + y*c6)))));
Packit 67cb25
    result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    result->val = sin(x) / x;
Packit 67cb25
    result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int gsl_sf_bessel_j1_e(const double x, gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  double ax = fabs(x);
Packit 67cb25
Packit 67cb25
  /* CHECK_POINTER(result) */
Packit 67cb25
Packit 67cb25
  if(x == 0.0) {
Packit 67cb25
    result->val = 0.0;
Packit 67cb25
    result->err = 0.0;
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
  else if(ax < 3.1*GSL_DBL_MIN) {
Packit 67cb25
    UNDERFLOW_ERROR(result);
Packit 67cb25
  }
Packit 67cb25
  else if(ax < 0.25) {
Packit 67cb25
    const double y = x*x;
Packit 67cb25
    const double c1 = -1.0/10.0;
Packit 67cb25
    const double c2 =  1.0/280.0;
Packit 67cb25
    const double c3 = -1.0/15120.0;
Packit 67cb25
    const double c4 =  1.0/1330560.0;
Packit 67cb25
    const double c5 = -1.0/172972800.0;
Packit 67cb25
    const double sum = 1.0 + y*(c1 + y*(c2 + y*(c3 + y*(c4 + y*c5))));
Packit 67cb25
    result->val = x/3.0 * sum;
Packit 67cb25
    result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    const double cos_x = cos(x);
Packit 67cb25
    const double sin_x = sin(x);
Packit 67cb25
    result->val  = (sin_x/x - cos_x)/x;
Packit 67cb25
    result->err  = 2.0 * GSL_DBL_EPSILON * (fabs(sin_x/(x*x)) + fabs(cos_x/x));
Packit 67cb25
    result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int gsl_sf_bessel_j2_e(const double x, gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  double ax = fabs(x);
Packit 67cb25
Packit 67cb25
  /* CHECK_POINTER(result) */
Packit 67cb25
  
Packit 67cb25
  if(x == 0.0) {
Packit 67cb25
    result->val = 0.0;
Packit 67cb25
    result->err = 0.0;
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
  else if(ax < 4.0*GSL_SQRT_DBL_MIN) {
Packit 67cb25
    UNDERFLOW_ERROR(result);
Packit 67cb25
  }
Packit 67cb25
  else if(ax < 1.3) {
Packit 67cb25
    const double y  = x*x;
Packit 67cb25
    const double c1 = -1.0/14.0;
Packit 67cb25
    const double c2 =  1.0/504.0;
Packit 67cb25
    const double c3 = -1.0/33264.0;
Packit 67cb25
    const double c4 =  1.0/3459456.0;
Packit 67cb25
    const double c5 = -1.0/518918400;
Packit 67cb25
    const double c6 =  1.0/105859353600.0;
Packit 67cb25
    const double c7 = -1.0/28158588057600.0;
Packit 67cb25
    const double c8 =  1.0/9461285587353600.0;
Packit 67cb25
    const double c9 = -1.0/3916972233164390400.0;
Packit 67cb25
    const double sum = 1.0+y*(c1+y*(c2+y*(c3+y*(c4+y*(c5+y*(c6+y*(c7+y*(c8+y*c9))))))));
Packit 67cb25
    result->val = y/15.0 * sum;
Packit 67cb25
    result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    /* bug #45730: switch from gsl_sf_{cos,sin} to cos/sin to fix large inputs */
Packit 67cb25
#if 0
Packit 67cb25
    gsl_sf_result cos_result;
Packit 67cb25
    gsl_sf_result sin_result;
Packit 67cb25
    const int stat_cos = gsl_sf_cos_e(x, &cos_result);
Packit 67cb25
    const int stat_sin = gsl_sf_sin_e(x, &sin_result);
Packit 67cb25
    const double cos_x = cos_result.val;
Packit 67cb25
    const double sin_x = sin_result.val;
Packit 67cb25
    const double f = (3.0/(x*x) - 1.0);
Packit 67cb25
    result->val  = (f * sin_x - 3.0*cos_x/x)/x;
Packit 67cb25
    result->err  = fabs(f * sin_result.err/x) + fabs((3.0*cos_result.err/x)/x);
Packit 67cb25
    result->err += 2.0 * GSL_DBL_EPSILON * (fabs(f*sin_x/x) + 3.0*fabs(cos_x/(x*x)));
Packit 67cb25
    result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    return GSL_ERROR_SELECT_2(stat_cos, stat_sin);
Packit 67cb25
#else
Packit 67cb25
    const double cos_x = cos(x);
Packit 67cb25
    const double sin_x = sin(x);
Packit 67cb25
    const double f = (3.0/(x*x) - 1.0);
Packit 67cb25
    result->val  = (f * sin_x - 3.0*cos_x/x)/x;
Packit 67cb25
    result->err = 2.0 * GSL_DBL_EPSILON * (fabs(f*sin_x/x) + 3.0*fabs(cos_x/(x*x)));
Packit 67cb25
    result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    /*return GSL_ERROR_SELECT_2(stat_cos, stat_sin);*/
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
#endif
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_sf_bessel_jl_e(const int l, const double x, gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  if(l < 0 || x < 0.0) {
Packit 67cb25
    DOMAIN_ERROR(result);
Packit 67cb25
  }
Packit 67cb25
  else if(x == 0.0) {
Packit 67cb25
    result->val = ( l > 0 ? 0.0 : 1.0 );
Packit 67cb25
    result->err = 0.0;
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
  else if(l == 0) {
Packit 67cb25
    return gsl_sf_bessel_j0_e(x, result);
Packit 67cb25
  }
Packit 67cb25
  else if(l == 1) {
Packit 67cb25
    return gsl_sf_bessel_j1_e(x, result);
Packit 67cb25
  }
Packit 67cb25
  else if(l == 2) {
Packit 67cb25
    return gsl_sf_bessel_j2_e(x, result);
Packit 67cb25
  }
Packit 67cb25
  else if(x*x < 10.0*(l+0.5)/M_E) {
Packit 67cb25
    gsl_sf_result b;
Packit 67cb25
    int status = gsl_sf_bessel_IJ_taylor_e(l+0.5, x, -1, 50, GSL_DBL_EPSILON, &b);
Packit 67cb25
    double pre   = sqrt((0.5*M_PI)/x);
Packit 67cb25
    result->val  = pre * b.val;
Packit 67cb25
    result->err  = pre * b.err;
Packit 67cb25
    result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    return status;
Packit 67cb25
  }
Packit 67cb25
  else if(GSL_ROOT4_DBL_EPSILON * x > (l*l + l + 1.0)) {
Packit 67cb25
    gsl_sf_result b;
Packit 67cb25
    int status = gsl_sf_bessel_Jnu_asympx_e(l + 0.5, x, &b);
Packit 67cb25
    double pre = sqrt((0.5*M_PI)/x);
Packit 67cb25
    result->val = pre * b.val;
Packit 67cb25
    result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val) + pre * b.err;
Packit 67cb25
    return status;
Packit 67cb25
  }
Packit 67cb25
  else if(l > 1.0/GSL_ROOT6_DBL_EPSILON) {
Packit 67cb25
    gsl_sf_result b;
Packit 67cb25
    int status = gsl_sf_bessel_Jnu_asymp_Olver_e(l + 0.5, x, &b);
Packit 67cb25
    double pre = sqrt((0.5*M_PI)/x);
Packit 67cb25
    result->val = pre * b.val;
Packit 67cb25
    result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val) + pre * b.err;
Packit 67cb25
    return status;
Packit 67cb25
  }
Packit 67cb25
  else if(x > 1000.0 && x > l*l)
Packit 67cb25
  {
Packit 67cb25
    /* We need this path to avoid feeding large x to CF1 below; */
Packit 67cb25
    gsl_sf_result b;
Packit 67cb25
    int status = gsl_sf_bessel_Jnu_asympx_e(l + 0.5, x, &b);
Packit 67cb25
    double pre = sqrt((0.5*M_PI)/x);
Packit 67cb25
    result->val = pre * b.val;
Packit 67cb25
    result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val) + pre * b.err;
Packit 67cb25
    return status;  
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    double sgn;
Packit 67cb25
    double ratio;
Packit 67cb25
    /* The CF1 call will hit 10000 iterations for x > 10000 + l */
Packit 67cb25
    int stat_CF1 = gsl_sf_bessel_J_CF1(l+0.5, x, &ratio, &sgn);
Packit 67cb25
    const double BESSEL_J_SMALL = GSL_DBL_MIN / GSL_DBL_EPSILON;
Packit 67cb25
    double jellp1 = BESSEL_J_SMALL * ratio;
Packit 67cb25
    double jell   = BESSEL_J_SMALL;
Packit 67cb25
    double jellm1;
Packit 67cb25
    int ell;
Packit 67cb25
    for(ell = l; ell > 0; ell--) {
Packit 67cb25
      jellm1 = -jellp1 + (2*ell + 1)/x * jell;
Packit 67cb25
      jellp1 = jell;
Packit 67cb25
      jell   = jellm1;
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    if(fabs(jell) > fabs(jellp1)) {
Packit 67cb25
      gsl_sf_result j0_result;
Packit 67cb25
      int stat_j0  = gsl_sf_bessel_j0_e(x, &j0_result);
Packit 67cb25
      double pre   = BESSEL_J_SMALL / jell;
Packit 67cb25
      result->val  = j0_result.val * pre;
Packit 67cb25
      result->err  = j0_result.err * fabs(pre);
Packit 67cb25
      result->err += 4.0 * GSL_DBL_EPSILON * (0.5*l + 1.0) * fabs(result->val);
Packit 67cb25
      return GSL_ERROR_SELECT_2(stat_j0, stat_CF1);
Packit 67cb25
    }
Packit 67cb25
    else {
Packit 67cb25
      gsl_sf_result j1_result;
Packit 67cb25
      int stat_j1  = gsl_sf_bessel_j1_e(x, &j1_result);
Packit 67cb25
      double pre   = BESSEL_J_SMALL / jellp1;
Packit 67cb25
      result->val  = j1_result.val * pre;
Packit 67cb25
      result->err  = j1_result.err * fabs(pre);
Packit 67cb25
      result->err += 4.0 * GSL_DBL_EPSILON * (0.5*l + 1.0) * fabs(result->val);
Packit 67cb25
      return GSL_ERROR_SELECT_2(stat_j1, stat_CF1);
Packit 67cb25
    }
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_sf_bessel_jl_array(const int lmax, const double x, double * result_array)
Packit 67cb25
{
Packit 67cb25
  /* CHECK_POINTER(result_array) */
Packit 67cb25
Packit 67cb25
  if(lmax < 0 || x < 0.0) {
Packit 67cb25
    int j;
Packit 67cb25
    for(j=0; j<=lmax; j++) result_array[j] = 0.0;
Packit 67cb25
    GSL_ERROR ("error", GSL_EDOM);
Packit 67cb25
  }
Packit 67cb25
  else if(x == 0.0) {
Packit 67cb25
    int j;
Packit 67cb25
    for(j=1; j<=lmax; j++) result_array[j] = 0.0;
Packit 67cb25
    result_array[0] = 1.0;
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    gsl_sf_result r_jellp1;
Packit 67cb25
    gsl_sf_result r_jell;
Packit 67cb25
    int stat_0 = gsl_sf_bessel_jl_e(lmax+1, x, &r_jellp1);
Packit 67cb25
    int stat_1 = gsl_sf_bessel_jl_e(lmax,   x, &r_jell);
Packit 67cb25
    double jellp1 = r_jellp1.val;
Packit 67cb25
    double jell   = r_jell.val;
Packit 67cb25
    double jellm1;
Packit 67cb25
    int ell;
Packit 67cb25
Packit 67cb25
    result_array[lmax] = jell;
Packit 67cb25
    for(ell = lmax; ell >= 1; ell--) {
Packit 67cb25
      jellm1 = -jellp1 + (2*ell + 1)/x * jell;
Packit 67cb25
      jellp1 = jell;
Packit 67cb25
      jell   = jellm1;
Packit 67cb25
      result_array[ell-1] = jellm1;
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    return GSL_ERROR_SELECT_2(stat_0, stat_1);
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int gsl_sf_bessel_jl_steed_array(const int lmax, const double x, double * jl_x)
Packit 67cb25
{
Packit 67cb25
  /* CHECK_POINTER(jl_x) */
Packit 67cb25
Packit 67cb25
  if(lmax < 0 || x < 0.0) {
Packit 67cb25
    int j;
Packit 67cb25
    for(j=0; j<=lmax; j++) jl_x[j] = 0.0;
Packit 67cb25
    GSL_ERROR ("error", GSL_EDOM);
Packit 67cb25
  }
Packit 67cb25
  else if(x == 0.0) {
Packit 67cb25
    int j;
Packit 67cb25
    for(j=1; j<=lmax; j++) jl_x[j] = 0.0;
Packit 67cb25
    jl_x[0] = 1.0;
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
  else if(x < 2.0*GSL_ROOT4_DBL_EPSILON) {
Packit 67cb25
    /* first two terms of Taylor series */
Packit 67cb25
    double inv_fact = 1.0;  /* 1/(1 3 5 ... (2l+1)) */
Packit 67cb25
    double x_l      = 1.0;  /* x^l */
Packit 67cb25
    int l;
Packit 67cb25
    for(l=0; l<=lmax; l++) {
Packit 67cb25
      jl_x[l]  = x_l * inv_fact;
Packit 67cb25
      jl_x[l] *= 1.0 - 0.5*x*x/(2.0*l+3.0);
Packit 67cb25
      inv_fact /= 2.0*l+3.0;
Packit 67cb25
      x_l      *= x;
Packit 67cb25
    }
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    /* Steed/Barnett algorithm [Comp. Phys. Comm. 21, 297 (1981)] */
Packit 67cb25
    double x_inv = 1.0/x;
Packit 67cb25
    double W = 2.0*x_inv;
Packit 67cb25
    double F = 1.0;
Packit 67cb25
    double FP = (lmax+1.0) * x_inv;
Packit 67cb25
    double B = 2.0*FP + x_inv;
Packit 67cb25
    double end = B + 20000.0*W;
Packit 67cb25
    double D = 1.0/B;
Packit 67cb25
    double del = -D;
Packit 67cb25
    
Packit 67cb25
    FP += del;
Packit 67cb25
    
Packit 67cb25
    /* continued fraction */
Packit 67cb25
    do {
Packit 67cb25
      B += W;
Packit 67cb25
      D = 1.0/(B-D);
Packit 67cb25
      del *= (B*D - 1.);
Packit 67cb25
      FP += del;
Packit 67cb25
      if(D < 0.0) F = -F;
Packit 67cb25
      if(B > end) {
Packit 67cb25
        GSL_ERROR ("error", GSL_EMAXITER);
Packit 67cb25
      }
Packit 67cb25
    }
Packit 67cb25
    while(fabs(del) >= fabs(FP) * GSL_DBL_EPSILON);
Packit 67cb25
    
Packit 67cb25
    FP *= F;
Packit 67cb25
    
Packit 67cb25
    if(lmax > 0) {
Packit 67cb25
      /* downward recursion */
Packit 67cb25
      double XP2 = FP;
Packit 67cb25
      double PL = lmax * x_inv;
Packit 67cb25
      int L  = lmax;
Packit 67cb25
      int LP;
Packit 67cb25
      jl_x[lmax] = F;
Packit 67cb25
      for(LP = 1; LP<=lmax; LP++) {
Packit 67cb25
        jl_x[L-1] = PL * jl_x[L] + XP2;
Packit 67cb25
        FP = PL*jl_x[L-1] - jl_x[L];
Packit 67cb25
        XP2 = FP;
Packit 67cb25
        PL -= x_inv;
Packit 67cb25
        --L;
Packit 67cb25
      }
Packit 67cb25
      F = jl_x[0];
Packit 67cb25
    }
Packit 67cb25
    
Packit 67cb25
    /* normalization */
Packit 67cb25
    W = x_inv / hypot(FP, F);
Packit 67cb25
    jl_x[0] = W*F;
Packit 67cb25
    if(lmax > 0) {
Packit 67cb25
      int L;
Packit 67cb25
      for(L=1; L<=lmax; L++) {
Packit 67cb25
        jl_x[L] *= W;
Packit 67cb25
      }
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
/*-*-*-*-*-*-*-*-*-* Functions w/ Natural Prototypes *-*-*-*-*-*-*-*-*-*-*/
Packit 67cb25
Packit 67cb25
#include "eval.h"
Packit 67cb25
Packit 67cb25
double gsl_sf_bessel_j0(const double x)
Packit 67cb25
{
Packit 67cb25
  EVAL_RESULT(gsl_sf_bessel_j0_e(x, &result));
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double gsl_sf_bessel_j1(const double x)
Packit 67cb25
{
Packit 67cb25
  EVAL_RESULT(gsl_sf_bessel_j1_e(x, &result));
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double gsl_sf_bessel_j2(const double x)
Packit 67cb25
{
Packit 67cb25
  EVAL_RESULT(gsl_sf_bessel_j2_e(x, &result));
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double gsl_sf_bessel_jl(const int l, const double x)
Packit 67cb25
{
Packit 67cb25
  EVAL_RESULT(gsl_sf_bessel_jl_e(l, x, &result));
Packit 67cb25
}
Packit 67cb25