Blame specfunc/bessel_k.c

Packit 67cb25
/* specfunc/bessel_k.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 1996, 1997, 1998, 1999, 2000 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_gamma.h>
Packit 67cb25
#include <gsl/gsl_sf_bessel.h>
Packit 67cb25
Packit 67cb25
#include "error.h"
Packit 67cb25
#include "check.h"
Packit 67cb25
Packit 67cb25
#include "bessel.h"
Packit 67cb25
Packit 67cb25
/*-*-*-*-*-*-*-*-*-*-*-* Private Section *-*-*-*-*-*-*-*-*-*-*-*/
Packit 67cb25
Packit 67cb25
/* [Abramowitz+Stegun, 10.2.4 + 10.2.6]
Packit 67cb25
 * with lmax=15, precision ~ 15D for x < 3
Packit 67cb25
 *
Packit 67cb25
 * assumes l >= 1
Packit 67cb25
 */
Packit 67cb25
static int bessel_kl_scaled_small_x(int l, const double x, gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  gsl_sf_result num_fact;
Packit 67cb25
  double den  = gsl_sf_pow_int(x, l+1);
Packit 67cb25
  int stat_df = gsl_sf_doublefact_e((unsigned int) (2*l-1), &num_fact);
Packit 67cb25
Packit 67cb25
  if(stat_df != GSL_SUCCESS || den == 0.0) {
Packit 67cb25
    OVERFLOW_ERROR(result);
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    const int lmax = 50;
Packit 67cb25
    gsl_sf_result ipos_term;
Packit 67cb25
    double ineg_term;
Packit 67cb25
    double sgn = (GSL_IS_ODD(l) ? -1.0 : 1.0);
Packit 67cb25
    double ex  = exp(x);
Packit 67cb25
    double t = 0.5*x*x;
Packit 67cb25
    double sum = 1.0;
Packit 67cb25
    double t_coeff = 1.0;
Packit 67cb25
    double t_power = 1.0;
Packit 67cb25
    double delta;
Packit 67cb25
    int stat_il;
Packit 67cb25
    int i;
Packit 67cb25
Packit 67cb25
    for(i=1; i
Packit 67cb25
      t_coeff /= i*(2*(i-l) - 1);
Packit 67cb25
      t_power *= t;
Packit 67cb25
      delta = t_power*t_coeff;
Packit 67cb25
      sum += delta;
Packit 67cb25
      if(fabs(delta/sum) < GSL_DBL_EPSILON) break;
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    stat_il = gsl_sf_bessel_il_scaled_e(l, x, &ipos_term);
Packit 67cb25
    ineg_term =  sgn * num_fact.val/den * sum;
Packit 67cb25
    result->val = -sgn * 0.5*M_PI * (ex*ipos_term.val - ineg_term);
Packit 67cb25
    result->val *= ex;
Packit 67cb25
    result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    return stat_il;
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
/*-*-*-*-*-*-*-*-*-*-*-* Functions with Error Codes *-*-*-*-*-*-*-*-*-*-*-*/
Packit 67cb25
Packit 67cb25
int gsl_sf_bessel_k0_scaled_e(const double x, gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  /* CHECK_POINTER(result) */
Packit 67cb25
Packit 67cb25
  if(x <= 0.0) {
Packit 67cb25
    DOMAIN_ERROR(result);
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    result->val = M_PI/(2.0*x);
Packit 67cb25
    result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    CHECK_UNDERFLOW(result);
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int gsl_sf_bessel_k1_scaled_e(const double x, gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  /* CHECK_POINTER(result) */
Packit 67cb25
Packit 67cb25
  if(x <= 0.0) {
Packit 67cb25
    DOMAIN_ERROR(result);
Packit 67cb25
  }
Packit 67cb25
  else if(x < (M_SQRTPI+1.0)/(M_SQRT2*GSL_SQRT_DBL_MAX)) {
Packit 67cb25
    OVERFLOW_ERROR(result);
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    result->val = M_PI/(2.0*x) * (1.0 + 1.0/x);
Packit 67cb25
    result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    CHECK_UNDERFLOW(result);
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int gsl_sf_bessel_k2_scaled_e(const double x, gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  /* CHECK_POINTER(result) */
Packit 67cb25
Packit 67cb25
  if(x <= 0.0) {
Packit 67cb25
    DOMAIN_ERROR(result);
Packit 67cb25
  }
Packit 67cb25
  else if(x < 2.0/GSL_ROOT3_DBL_MAX) {
Packit 67cb25
    OVERFLOW_ERROR(result);
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    result->val = M_PI/(2.0*x) * (1.0 + 3.0/x * (1.0 + 1.0/x));
Packit 67cb25
    result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    CHECK_UNDERFLOW(result);
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int gsl_sf_bessel_kl_scaled_e(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(l == 0) {
Packit 67cb25
    return gsl_sf_bessel_k0_scaled_e(x, result);
Packit 67cb25
  }
Packit 67cb25
  else if(l == 1) {
Packit 67cb25
    return gsl_sf_bessel_k1_scaled_e(x, result);
Packit 67cb25
  }
Packit 67cb25
  else if(l == 2) {
Packit 67cb25
    return gsl_sf_bessel_k2_scaled_e(x, result);
Packit 67cb25
  }
Packit 67cb25
  else if(x < 3.0) {
Packit 67cb25
    return bessel_kl_scaled_small_x(l, x, result);
Packit 67cb25
  }
Packit 67cb25
  else if(GSL_ROOT3_DBL_EPSILON * x > (l*l + l + 1)) {
Packit 67cb25
    int status = gsl_sf_bessel_Knu_scaled_asympx_e(l + 0.5, x, result);
Packit 67cb25
    double pre = sqrt((0.5*M_PI)/x);
Packit 67cb25
    result->val *= pre;
Packit 67cb25
    result->err *= pre;
Packit 67cb25
    return status;
Packit 67cb25
  }
Packit 67cb25
  else if(GSL_MIN(0.29/(l*l+1.0), 0.5/(l*l+1.0+x*x)) < GSL_ROOT3_DBL_EPSILON) {
Packit 67cb25
    int status = gsl_sf_bessel_Knu_scaled_asymp_unif_e(l + 0.5, x, result);
Packit 67cb25
    double pre = sqrt((0.5*M_PI)/x);
Packit 67cb25
    result->val *= pre;
Packit 67cb25
    result->err *= pre;
Packit 67cb25
    return status;
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    /* recurse upward */
Packit 67cb25
    gsl_sf_result r_bk;
Packit 67cb25
    gsl_sf_result r_bkm;
Packit 67cb25
    int stat_1 = gsl_sf_bessel_k1_scaled_e(x, &r_bk);
Packit 67cb25
    int stat_0 = gsl_sf_bessel_k0_scaled_e(x, &r_bkm);
Packit 67cb25
    double bkp;
Packit 67cb25
    double bk  = r_bk.val;
Packit 67cb25
    double bkm = r_bkm.val;
Packit 67cb25
    int j;
Packit 67cb25
    for(j=1; j
Packit 67cb25
      bkp = (2*j+1)/x*bk + bkm;
Packit 67cb25
      bkm = bk;
Packit 67cb25
      bk  = bkp;
Packit 67cb25
    }
Packit 67cb25
    result->val  = bk;
Packit 67cb25
    result->err  = fabs(bk) * (fabs(r_bk.err/r_bk.val) + fabs(r_bkm.err/r_bkm.val));
Packit 67cb25
    result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
Packit 67cb25
    return GSL_ERROR_SELECT_2(stat_1, stat_0);
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
int 
Packit 67cb25
gsl_sf_bessel_kl_scaled_array(const int lmax, const double x, double * result_array)
Packit 67cb25
{
Packit 67cb25
  if(lmax < 0 || x <= 0.0) {
Packit 67cb25
    GSL_ERROR("domain error", GSL_EDOM);
Packit 67cb25
  } else if (lmax == 0) {
Packit 67cb25
    gsl_sf_result result;
Packit 67cb25
    int stat = gsl_sf_bessel_k0_scaled_e(x, &result);
Packit 67cb25
    result_array[0] = result.val;
Packit 67cb25
    return stat;
Packit 67cb25
  } else {
Packit 67cb25
    int ell;
Packit 67cb25
    double kellp1, kell, kellm1;
Packit 67cb25
    gsl_sf_result r_kell;
Packit 67cb25
    gsl_sf_result r_kellm1;
Packit 67cb25
    gsl_sf_bessel_k1_scaled_e(x, &r_kell);
Packit 67cb25
    gsl_sf_bessel_k0_scaled_e(x, &r_kellm1);
Packit 67cb25
    kell   = r_kell.val;
Packit 67cb25
    kellm1 = r_kellm1.val;
Packit 67cb25
    result_array[0] = kellm1;
Packit 67cb25
    result_array[1] = kell;
Packit 67cb25
    for(ell = 1; ell < lmax; ell++) {
Packit 67cb25
      kellp1 = (2*ell+1)/x * kell + kellm1;
Packit 67cb25
      result_array[ell+1] = kellp1;
Packit 67cb25
      kellm1 = kell;
Packit 67cb25
      kell   = kellp1;
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_k0_scaled(const double x)
Packit 67cb25
{
Packit 67cb25
  EVAL_RESULT(gsl_sf_bessel_k0_scaled_e(x, &result));
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double gsl_sf_bessel_k1_scaled(const double x)
Packit 67cb25
{
Packit 67cb25
  EVAL_RESULT(gsl_sf_bessel_k1_scaled_e(x, &result));
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double gsl_sf_bessel_k2_scaled(const double x)
Packit 67cb25
{
Packit 67cb25
  EVAL_RESULT(gsl_sf_bessel_k2_scaled_e(x, &result));
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double gsl_sf_bessel_kl_scaled(const int l, const double x)
Packit 67cb25
{
Packit 67cb25
  EVAL_RESULT(gsl_sf_bessel_kl_scaled_e(l, x, &result));
Packit 67cb25
}
Packit 67cb25
Packit 67cb25