Blame specfunc/hyperg_2F1.c

Packit 67cb25
/* specfunc/hyperg_2F1.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2004 Gerard Jungman
Packit 67cb25
 * Copyright (C) 2009 Brian Gough
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_exp.h>
Packit 67cb25
#include <gsl/gsl_sf_pow_int.h>
Packit 67cb25
#include <gsl/gsl_sf_gamma.h>
Packit 67cb25
#include <gsl/gsl_sf_psi.h>
Packit 67cb25
#include <gsl/gsl_sf_hyperg.h>
Packit 67cb25
Packit 67cb25
#include "error.h"
Packit 67cb25
Packit 67cb25
#define locEPS (1000.0*GSL_DBL_EPSILON)
Packit 67cb25
Packit 67cb25
Packit 67cb25
/* Assumes c != negative integer.
Packit 67cb25
 */
Packit 67cb25
static int
Packit 67cb25
hyperg_2F1_series(const double a, const double b, const double c,
Packit 67cb25
                  const double x, 
Packit 67cb25
                  gsl_sf_result * result
Packit 67cb25
                  )
Packit 67cb25
{
Packit 67cb25
  double sum_pos = 1.0;
Packit 67cb25
  double sum_neg = 0.0;
Packit 67cb25
  double del_pos = 1.0;
Packit 67cb25
  double del_neg = 0.0;
Packit 67cb25
  double del = 1.0;
Packit 67cb25
  double del_prev;
Packit 67cb25
  double k = 0.0;
Packit 67cb25
  int i = 0;
Packit 67cb25
Packit 67cb25
  if(fabs(c) < GSL_DBL_EPSILON) {
Packit 67cb25
    result->val = 0.0; /* FIXME: ?? */
Packit 67cb25
    result->err = 1.0;
Packit 67cb25
    GSL_ERROR ("error", GSL_EDOM);
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  do {
Packit 67cb25
    if(++i > 30000) {
Packit 67cb25
      result->val  = sum_pos - sum_neg;
Packit 67cb25
      result->err  = del_pos + del_neg;
Packit 67cb25
      result->err += 2.0 * GSL_DBL_EPSILON * (sum_pos + sum_neg);
Packit 67cb25
      result->err += 2.0 * GSL_DBL_EPSILON * (2.0*sqrt(k)+1.0) * fabs(result->val);
Packit 67cb25
      GSL_ERROR ("error", GSL_EMAXITER);
Packit 67cb25
    }
Packit 67cb25
    del_prev = del;
Packit 67cb25
    del *= (a+k)*(b+k) * x / ((c+k) * (k+1.0));  /* Gauss series */
Packit 67cb25
Packit 67cb25
    if(del > 0.0) {
Packit 67cb25
      del_pos  =  del;
Packit 67cb25
      sum_pos +=  del;
Packit 67cb25
    }
Packit 67cb25
    else if(del == 0.0) {
Packit 67cb25
      /* Exact termination (a or b was a negative integer).
Packit 67cb25
       */
Packit 67cb25
      del_pos = 0.0;
Packit 67cb25
      del_neg = 0.0;
Packit 67cb25
      break;
Packit 67cb25
    }
Packit 67cb25
    else {
Packit 67cb25
      del_neg  = -del;
Packit 67cb25
      sum_neg -=  del;
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    /*
Packit 67cb25
     * This stopping criteria is taken from the thesis
Packit 67cb25
     * "Computation of Hypergeometic Functions" by J. Pearson, pg. 31
Packit 67cb25
     * (http://people.maths.ox.ac.uk/porterm/research/pearson_final.pdf)
Packit 67cb25
     * and fixes bug #45926
Packit 67cb25
     */
Packit 67cb25
    if (fabs(del_prev / (sum_pos - sum_neg)) < GSL_DBL_EPSILON &&
Packit 67cb25
        fabs(del / (sum_pos - sum_neg)) < GSL_DBL_EPSILON)
Packit 67cb25
      break;
Packit 67cb25
Packit 67cb25
    k += 1.0;
Packit 67cb25
  } while(fabs((del_pos + del_neg)/(sum_pos-sum_neg)) > GSL_DBL_EPSILON);
Packit 67cb25
Packit 67cb25
  result->val  = sum_pos - sum_neg;
Packit 67cb25
  result->err  = del_pos + del_neg;
Packit 67cb25
  result->err += 2.0 * GSL_DBL_EPSILON * (sum_pos + sum_neg);
Packit 67cb25
  result->err += 2.0 * GSL_DBL_EPSILON * (2.0*sqrt(k) + 1.0) * fabs(result->val);
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
/* a = aR + i aI, b = aR - i aI */
Packit 67cb25
static
Packit 67cb25
int
Packit 67cb25
hyperg_2F1_conj_series(const double aR, const double aI, const double c,
Packit 67cb25
                       double x,
Packit 67cb25
                       gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  if(c == 0.0) {
Packit 67cb25
    result->val = 0.0; /* FIXME: should be Inf */
Packit 67cb25
    result->err = 0.0;
Packit 67cb25
    GSL_ERROR ("error", GSL_EDOM);
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    double sum_pos = 1.0;
Packit 67cb25
    double sum_neg = 0.0;
Packit 67cb25
    double del_pos = 1.0;
Packit 67cb25
    double del_neg = 0.0;
Packit 67cb25
    double del = 1.0;
Packit 67cb25
    double k = 0.0;
Packit 67cb25
    do {
Packit 67cb25
      del *= ((aR+k)*(aR+k) + aI*aI)/((k+1.0)*(c+k)) * x;
Packit 67cb25
Packit 67cb25
      if(del >= 0.0) {
Packit 67cb25
        del_pos  =  del;
Packit 67cb25
        sum_pos +=  del;
Packit 67cb25
      }
Packit 67cb25
      else {
Packit 67cb25
        del_neg  = -del;
Packit 67cb25
        sum_neg -=  del;
Packit 67cb25
      }
Packit 67cb25
Packit 67cb25
      if(k > 30000) {
Packit 67cb25
        result->val  = sum_pos - sum_neg;
Packit 67cb25
        result->err  = del_pos + del_neg;
Packit 67cb25
        result->err += 2.0 * GSL_DBL_EPSILON * (sum_pos + sum_neg);
Packit 67cb25
        result->err += 2.0 * GSL_DBL_EPSILON * (2.0*sqrt(k)+1.0) * fabs(result->val);
Packit 67cb25
        GSL_ERROR ("error", GSL_EMAXITER);
Packit 67cb25
      }
Packit 67cb25
Packit 67cb25
      k += 1.0;
Packit 67cb25
    } while(fabs((del_pos + del_neg)/(sum_pos - sum_neg)) > GSL_DBL_EPSILON);
Packit 67cb25
Packit 67cb25
    result->val  = sum_pos - sum_neg;
Packit 67cb25
    result->err  = del_pos + del_neg;
Packit 67cb25
    result->err += 2.0 * GSL_DBL_EPSILON * (sum_pos + sum_neg);
Packit 67cb25
    result->err += 2.0 * GSL_DBL_EPSILON * (2.0*sqrt(k) + 1.0) * fabs(result->val);
Packit 67cb25
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
/* Luke's rational approximation. The most accesible
Packit 67cb25
 * discussion is in [Kolbig, CPC 23, 51 (1981)].
Packit 67cb25
 * The convergence is supposedly guaranteed for x < 0.
Packit 67cb25
 * You have to read Luke's books to see this and other
Packit 67cb25
 * results. Unfortunately, the stability is not so
Packit 67cb25
 * clear to me, although it seems very efficient when
Packit 67cb25
 * it works.
Packit 67cb25
 */
Packit 67cb25
static
Packit 67cb25
int
Packit 67cb25
hyperg_2F1_luke(const double a, const double b, const double c,
Packit 67cb25
                const double xin, 
Packit 67cb25
                gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  int stat_iter;
Packit 67cb25
  const double RECUR_BIG = 1.0e+50;
Packit 67cb25
  const int nmax = 20000;
Packit 67cb25
  int n = 3;
Packit 67cb25
  const double x  = -xin;
Packit 67cb25
  const double x3 = x*x*x;
Packit 67cb25
  const double t0 = a*b/c;
Packit 67cb25
  const double t1 = (a+1.0)*(b+1.0)/(2.0*c);
Packit 67cb25
  const double t2 = (a+2.0)*(b+2.0)/(2.0*(c+1.0));
Packit 67cb25
  double F = 1.0;
Packit 67cb25
  double prec;
Packit 67cb25
Packit 67cb25
  double Bnm3 = 1.0;                                  /* B0 */
Packit 67cb25
  double Bnm2 = 1.0 + t1 * x;                         /* B1 */
Packit 67cb25
  double Bnm1 = 1.0 + t2 * x * (1.0 + t1/3.0 * x);    /* B2 */
Packit 67cb25
 
Packit 67cb25
  double Anm3 = 1.0;                                                      /* A0 */
Packit 67cb25
  double Anm2 = Bnm2 - t0 * x;                                            /* A1 */
Packit 67cb25
  double Anm1 = Bnm1 - t0*(1.0 + t2*x)*x + t0 * t1 * (c/(c+1.0)) * x*x;   /* A2 */
Packit 67cb25
Packit 67cb25
  while(1) {
Packit 67cb25
    double npam1 = n + a - 1;
Packit 67cb25
    double npbm1 = n + b - 1;
Packit 67cb25
    double npcm1 = n + c - 1;
Packit 67cb25
    double npam2 = n + a - 2;
Packit 67cb25
    double npbm2 = n + b - 2;
Packit 67cb25
    double npcm2 = n + c - 2;
Packit 67cb25
    double tnm1  = 2*n - 1;
Packit 67cb25
    double tnm3  = 2*n - 3;
Packit 67cb25
    double tnm5  = 2*n - 5;
Packit 67cb25
    double n2 = n*n;
Packit 67cb25
    double F1 =  (3.0*n2 + (a+b-6)*n + 2 - a*b - 2*(a+b)) / (2*tnm3*npcm1);
Packit 67cb25
    double F2 = -(3.0*n2 - (a+b+6)*n + 2 - a*b)*npam1*npbm1/(4*tnm1*tnm3*npcm2*npcm1);
Packit 67cb25
    double F3 = (npam2*npam1*npbm2*npbm1*(n-a-2)*(n-b-2)) / (8*tnm3*tnm3*tnm5*(n+c-3)*npcm2*npcm1);
Packit 67cb25
    double E  = -npam1*npbm1*(n-c-1) / (2*tnm3*npcm2*npcm1);
Packit 67cb25
Packit 67cb25
    double An = (1.0+F1*x)*Anm1 + (E + F2*x)*x*Anm2 + F3*x3*Anm3;
Packit 67cb25
    double Bn = (1.0+F1*x)*Bnm1 + (E + F2*x)*x*Bnm2 + F3*x3*Bnm3;
Packit 67cb25
    double r = An/Bn;
Packit 67cb25
Packit 67cb25
    prec = fabs((F - r)/F);
Packit 67cb25
    F = r;
Packit 67cb25
Packit 67cb25
    if(prec < GSL_DBL_EPSILON || n > nmax) break;
Packit 67cb25
Packit 67cb25
    if(fabs(An) > RECUR_BIG || fabs(Bn) > RECUR_BIG) {
Packit 67cb25
      An   /= RECUR_BIG;
Packit 67cb25
      Bn   /= RECUR_BIG;
Packit 67cb25
      Anm1 /= RECUR_BIG;
Packit 67cb25
      Bnm1 /= RECUR_BIG;
Packit 67cb25
      Anm2 /= RECUR_BIG;
Packit 67cb25
      Bnm2 /= RECUR_BIG;
Packit 67cb25
      Anm3 /= RECUR_BIG;
Packit 67cb25
      Bnm3 /= RECUR_BIG;
Packit 67cb25
    }
Packit 67cb25
    else if(fabs(An) < 1.0/RECUR_BIG || fabs(Bn) < 1.0/RECUR_BIG) {
Packit 67cb25
      An   *= RECUR_BIG;
Packit 67cb25
      Bn   *= RECUR_BIG;
Packit 67cb25
      Anm1 *= RECUR_BIG;
Packit 67cb25
      Bnm1 *= RECUR_BIG;
Packit 67cb25
      Anm2 *= RECUR_BIG;
Packit 67cb25
      Bnm2 *= RECUR_BIG;
Packit 67cb25
      Anm3 *= RECUR_BIG;
Packit 67cb25
      Bnm3 *= RECUR_BIG;
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    n++;
Packit 67cb25
    Bnm3 = Bnm2;
Packit 67cb25
    Bnm2 = Bnm1;
Packit 67cb25
    Bnm1 = Bn;
Packit 67cb25
    Anm3 = Anm2;
Packit 67cb25
    Anm2 = Anm1;
Packit 67cb25
    Anm1 = An;
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  result->val  = F;
Packit 67cb25
  result->err  = 2.0 * fabs(prec * F);
Packit 67cb25
  result->err += 2.0 * GSL_DBL_EPSILON * (n+1.0) * fabs(F);
Packit 67cb25
Packit 67cb25
  /* FIXME: just a hack: there's a lot of shit going on here */
Packit 67cb25
  result->err *= 8.0 * (fabs(a) + fabs(b) + 1.0);
Packit 67cb25
Packit 67cb25
  stat_iter = (n >= nmax ? GSL_EMAXITER : GSL_SUCCESS );
Packit 67cb25
Packit 67cb25
  return stat_iter;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
/* Luke's rational approximation for the
Packit 67cb25
 * case a = aR + i aI, b = aR - i aI.
Packit 67cb25
 */
Packit 67cb25
static
Packit 67cb25
int
Packit 67cb25
hyperg_2F1_conj_luke(const double aR, const double aI, const double c,
Packit 67cb25
                     const double xin, 
Packit 67cb25
                     gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  int stat_iter;
Packit 67cb25
  const double RECUR_BIG = 1.0e+50;
Packit 67cb25
  const int nmax = 10000;
Packit 67cb25
  int n = 3;
Packit 67cb25
  const double x = -xin;
Packit 67cb25
  const double x3 = x*x*x;
Packit 67cb25
  const double atimesb = aR*aR + aI*aI;
Packit 67cb25
  const double apb     = 2.0*aR;
Packit 67cb25
  const double t0 = atimesb/c;
Packit 67cb25
  const double t1 = (atimesb +     apb + 1.0)/(2.0*c);
Packit 67cb25
  const double t2 = (atimesb + 2.0*apb + 4.0)/(2.0*(c+1.0));
Packit 67cb25
  double F = 1.0;
Packit 67cb25
  double prec;
Packit 67cb25
Packit 67cb25
  double Bnm3 = 1.0;                                  /* B0 */
Packit 67cb25
  double Bnm2 = 1.0 + t1 * x;                         /* B1 */
Packit 67cb25
  double Bnm1 = 1.0 + t2 * x * (1.0 + t1/3.0 * x);    /* B2 */
Packit 67cb25
 
Packit 67cb25
  double Anm3 = 1.0;                                                      /* A0 */
Packit 67cb25
  double Anm2 = Bnm2 - t0 * x;                                            /* A1 */
Packit 67cb25
  double Anm1 = Bnm1 - t0*(1.0 + t2*x)*x + t0 * t1 * (c/(c+1.0)) * x*x;   /* A2 */
Packit 67cb25
Packit 67cb25
  while(1) {
Packit 67cb25
    double nm1 = n - 1;
Packit 67cb25
    double nm2 = n - 2;
Packit 67cb25
    double npam1_npbm1 = atimesb + nm1*apb + nm1*nm1;
Packit 67cb25
    double npam2_npbm2 = atimesb + nm2*apb + nm2*nm2;
Packit 67cb25
    double npcm1 = nm1 + c;
Packit 67cb25
    double npcm2 = nm2 + c;
Packit 67cb25
    double tnm1  = 2*n - 1;
Packit 67cb25
    double tnm3  = 2*n - 3;
Packit 67cb25
    double tnm5  = 2*n - 5;
Packit 67cb25
    double n2 = n*n;
Packit 67cb25
    double F1 =  (3.0*n2 + (apb-6)*n + 2 - atimesb - 2*apb) / (2*tnm3*npcm1);
Packit 67cb25
    double F2 = -(3.0*n2 - (apb+6)*n + 2 - atimesb)*npam1_npbm1/(4*tnm1*tnm3*npcm2*npcm1);
Packit 67cb25
    double F3 = (npam2_npbm2*npam1_npbm1*(nm2*nm2 - nm2*apb + atimesb)) / (8*tnm3*tnm3*tnm5*(n+c-3)*npcm2*npcm1);
Packit 67cb25
    double E  = -npam1_npbm1*(n-c-1) / (2*tnm3*npcm2*npcm1);
Packit 67cb25
Packit 67cb25
    double An = (1.0+F1*x)*Anm1 + (E + F2*x)*x*Anm2 + F3*x3*Anm3;
Packit 67cb25
    double Bn = (1.0+F1*x)*Bnm1 + (E + F2*x)*x*Bnm2 + F3*x3*Bnm3;
Packit 67cb25
    double r = An/Bn;
Packit 67cb25
Packit 67cb25
    prec = fabs(F - r)/fabs(F);
Packit 67cb25
    F = r;
Packit 67cb25
Packit 67cb25
    if(prec < GSL_DBL_EPSILON || n > nmax) break;
Packit 67cb25
Packit 67cb25
    if(fabs(An) > RECUR_BIG || fabs(Bn) > RECUR_BIG) {
Packit 67cb25
      An   /= RECUR_BIG;
Packit 67cb25
      Bn   /= RECUR_BIG;
Packit 67cb25
      Anm1 /= RECUR_BIG;
Packit 67cb25
      Bnm1 /= RECUR_BIG;
Packit 67cb25
      Anm2 /= RECUR_BIG;
Packit 67cb25
      Bnm2 /= RECUR_BIG;
Packit 67cb25
      Anm3 /= RECUR_BIG;
Packit 67cb25
      Bnm3 /= RECUR_BIG;
Packit 67cb25
    }
Packit 67cb25
    else if(fabs(An) < 1.0/RECUR_BIG || fabs(Bn) < 1.0/RECUR_BIG) {
Packit 67cb25
      An   *= RECUR_BIG;
Packit 67cb25
      Bn   *= RECUR_BIG;
Packit 67cb25
      Anm1 *= RECUR_BIG;
Packit 67cb25
      Bnm1 *= RECUR_BIG;
Packit 67cb25
      Anm2 *= RECUR_BIG;
Packit 67cb25
      Bnm2 *= RECUR_BIG;
Packit 67cb25
      Anm3 *= RECUR_BIG;
Packit 67cb25
      Bnm3 *= RECUR_BIG;
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    n++;
Packit 67cb25
    Bnm3 = Bnm2;
Packit 67cb25
    Bnm2 = Bnm1;
Packit 67cb25
    Bnm1 = Bn;
Packit 67cb25
    Anm3 = Anm2;
Packit 67cb25
    Anm2 = Anm1;
Packit 67cb25
    Anm1 = An;
Packit 67cb25
  }
Packit 67cb25
  
Packit 67cb25
  result->val  = F;
Packit 67cb25
  result->err  = 2.0 * fabs(prec * F);
Packit 67cb25
  result->err += 2.0 * GSL_DBL_EPSILON * (n+1.0) * fabs(F);
Packit 67cb25
Packit 67cb25
  /* FIXME: see above */
Packit 67cb25
  result->err *= 8.0 * (fabs(aR) + fabs(aI) + 1.0);
Packit 67cb25
Packit 67cb25
  stat_iter = (n >= nmax ? GSL_EMAXITER : GSL_SUCCESS );
Packit 67cb25
Packit 67cb25
  return stat_iter;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
/* Do the reflection described in [Moshier, p. 334].
Packit 67cb25
 * Assumes a,b,c != neg integer.
Packit 67cb25
 */
Packit 67cb25
static
Packit 67cb25
int
Packit 67cb25
hyperg_2F1_reflect(const double a, const double b, const double c,
Packit 67cb25
                   const double x, gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  const double d = c - a - b;
Packit 67cb25
  const int intd  = floor(d+0.5);
Packit 67cb25
  const int d_integer = ( fabs(d - intd) < locEPS );
Packit 67cb25
Packit 67cb25
  if(d_integer) {
Packit 67cb25
    const double ln_omx = log(1.0 - x);
Packit 67cb25
    const double ad = fabs(d);
Packit 67cb25
    int stat_F2 = GSL_SUCCESS;
Packit 67cb25
    double sgn_2;
Packit 67cb25
    gsl_sf_result F1;
Packit 67cb25
    gsl_sf_result F2;
Packit 67cb25
    double d1, d2;
Packit 67cb25
    gsl_sf_result lng_c;
Packit 67cb25
    gsl_sf_result lng_ad2;
Packit 67cb25
    gsl_sf_result lng_bd2;
Packit 67cb25
    int stat_c;
Packit 67cb25
    int stat_ad2;
Packit 67cb25
    int stat_bd2;
Packit 67cb25
Packit 67cb25
    if(d >= 0.0) {
Packit 67cb25
      d1 = d;
Packit 67cb25
      d2 = 0.0;
Packit 67cb25
    }
Packit 67cb25
    else {
Packit 67cb25
      d1 = 0.0;
Packit 67cb25
      d2 = d;
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    stat_ad2 = gsl_sf_lngamma_e(a+d2, &lng_ad2);
Packit 67cb25
    stat_bd2 = gsl_sf_lngamma_e(b+d2, &lng_bd2);
Packit 67cb25
    stat_c   = gsl_sf_lngamma_e(c,    &lng_c);
Packit 67cb25
Packit 67cb25
    /* Evaluate F1.
Packit 67cb25
     */
Packit 67cb25
    if(ad < GSL_DBL_EPSILON) {
Packit 67cb25
      /* d = 0 */
Packit 67cb25
      F1.val = 0.0;
Packit 67cb25
      F1.err = 0.0;
Packit 67cb25
    }
Packit 67cb25
    else {
Packit 67cb25
      gsl_sf_result lng_ad;
Packit 67cb25
      gsl_sf_result lng_ad1;
Packit 67cb25
      gsl_sf_result lng_bd1;
Packit 67cb25
      int stat_ad  = gsl_sf_lngamma_e(ad,   &lng_ad);
Packit 67cb25
      int stat_ad1 = gsl_sf_lngamma_e(a+d1, &lng_ad1);
Packit 67cb25
      int stat_bd1 = gsl_sf_lngamma_e(b+d1, &lng_bd1);
Packit 67cb25
Packit 67cb25
      if(stat_ad1 == GSL_SUCCESS && stat_bd1 == GSL_SUCCESS && stat_ad == GSL_SUCCESS) {
Packit 67cb25
        /* Gamma functions in the denominator are ok.
Packit 67cb25
         * Proceed with evaluation.
Packit 67cb25
         */
Packit 67cb25
        int i;
Packit 67cb25
        double sum1 = 1.0;
Packit 67cb25
        double term = 1.0;
Packit 67cb25
        double ln_pre1_val = lng_ad.val + lng_c.val + d2*ln_omx - lng_ad1.val - lng_bd1.val;
Packit 67cb25
        double ln_pre1_err = lng_ad.err + lng_c.err + lng_ad1.err + lng_bd1.err + GSL_DBL_EPSILON * fabs(ln_pre1_val);
Packit 67cb25
        int stat_e;
Packit 67cb25
Packit 67cb25
        /* Do F1 sum.
Packit 67cb25
         */
Packit 67cb25
        for(i=1; i
Packit 67cb25
          int j = i-1;
Packit 67cb25
          term *= (a + d2 + j) * (b + d2 + j) / (1.0 + d2 + j) / i * (1.0-x);
Packit 67cb25
          sum1 += term;
Packit 67cb25
        }
Packit 67cb25
        
Packit 67cb25
        stat_e = gsl_sf_exp_mult_err_e(ln_pre1_val, ln_pre1_err,
Packit 67cb25
                                       sum1, GSL_DBL_EPSILON*fabs(sum1),
Packit 67cb25
                                       &F1;;
Packit 67cb25
        if(stat_e == GSL_EOVRFLW) {
Packit 67cb25
          OVERFLOW_ERROR(result);
Packit 67cb25
        }
Packit 67cb25
      }
Packit 67cb25
      else {
Packit 67cb25
        /* Gamma functions in the denominator were not ok.
Packit 67cb25
         * So the F1 term is zero.
Packit 67cb25
         */
Packit 67cb25
        F1.val = 0.0;
Packit 67cb25
        F1.err = 0.0;
Packit 67cb25
      }
Packit 67cb25
    } /* end F1 evaluation */
Packit 67cb25
Packit 67cb25
Packit 67cb25
    /* Evaluate F2.
Packit 67cb25
     */
Packit 67cb25
    if(stat_ad2 == GSL_SUCCESS && stat_bd2 == GSL_SUCCESS) {
Packit 67cb25
      /* Gamma functions in the denominator are ok.
Packit 67cb25
       * Proceed with evaluation.
Packit 67cb25
       */
Packit 67cb25
      const int maxiter = 2000;
Packit 67cb25
      double psi_1 = -M_EULER;
Packit 67cb25
      gsl_sf_result psi_1pd; 
Packit 67cb25
      gsl_sf_result psi_apd1;
Packit 67cb25
      gsl_sf_result psi_bpd1;
Packit 67cb25
      int stat_1pd  = gsl_sf_psi_e(1.0 + ad, &psi_1pd);
Packit 67cb25
      int stat_apd1 = gsl_sf_psi_e(a + d1,   &psi_apd1);
Packit 67cb25
      int stat_bpd1 = gsl_sf_psi_e(b + d1,   &psi_bpd1);
Packit 67cb25
      int stat_dall = GSL_ERROR_SELECT_3(stat_1pd, stat_apd1, stat_bpd1);
Packit 67cb25
Packit 67cb25
      double psi_val = psi_1 + psi_1pd.val - psi_apd1.val - psi_bpd1.val - ln_omx;
Packit 67cb25
      double psi_err = psi_1pd.err + psi_apd1.err + psi_bpd1.err + GSL_DBL_EPSILON*fabs(psi_val);
Packit 67cb25
      double fact = 1.0;
Packit 67cb25
      double sum2_val = psi_val;
Packit 67cb25
      double sum2_err = psi_err;
Packit 67cb25
      double ln_pre2_val = lng_c.val + d1*ln_omx - lng_ad2.val - lng_bd2.val;
Packit 67cb25
      double ln_pre2_err = lng_c.err + lng_ad2.err + lng_bd2.err + GSL_DBL_EPSILON * fabs(ln_pre2_val);
Packit 67cb25
      int stat_e;
Packit 67cb25
Packit 67cb25
      int j;
Packit 67cb25
Packit 67cb25
      /* Do F2 sum.
Packit 67cb25
       */
Packit 67cb25
      for(j=1; j
Packit 67cb25
        /* values for psi functions use recurrence; Abramowitz+Stegun 6.3.5 */
Packit 67cb25
        double term1 = 1.0/(double)j  + 1.0/(ad+j);
Packit 67cb25
        double term2 = 1.0/(a+d1+j-1.0) + 1.0/(b+d1+j-1.0);
Packit 67cb25
        double delta = 0.0;
Packit 67cb25
        psi_val += term1 - term2;
Packit 67cb25
        psi_err += GSL_DBL_EPSILON * (fabs(term1) + fabs(term2));
Packit 67cb25
        fact *= (a+d1+j-1.0)*(b+d1+j-1.0)/((ad+j)*j) * (1.0-x);
Packit 67cb25
        delta = fact * psi_val;
Packit 67cb25
        sum2_val += delta;
Packit 67cb25
        sum2_err += fabs(fact * psi_err) + GSL_DBL_EPSILON*fabs(delta);
Packit 67cb25
        if(fabs(delta) < GSL_DBL_EPSILON * fabs(sum2_val)) break;
Packit 67cb25
      }
Packit 67cb25
Packit 67cb25
      if(j == maxiter) stat_F2 = GSL_EMAXITER;
Packit 67cb25
Packit 67cb25
      if(sum2_val == 0.0) {
Packit 67cb25
        F2.val = 0.0;
Packit 67cb25
        F2.err = 0.0;
Packit 67cb25
      }
Packit 67cb25
      else {
Packit 67cb25
        stat_e = gsl_sf_exp_mult_err_e(ln_pre2_val, ln_pre2_err,
Packit 67cb25
                                       sum2_val, sum2_err,
Packit 67cb25
                                       &F2;;
Packit 67cb25
        if(stat_e == GSL_EOVRFLW) {
Packit 67cb25
          result->val = 0.0;
Packit 67cb25
          result->err = 0.0;
Packit 67cb25
          GSL_ERROR ("error", GSL_EOVRFLW);
Packit 67cb25
        }
Packit 67cb25
      }
Packit 67cb25
      stat_F2 = GSL_ERROR_SELECT_2(stat_F2, stat_dall);
Packit 67cb25
    }
Packit 67cb25
    else {
Packit 67cb25
      /* Gamma functions in the denominator not ok.
Packit 67cb25
       * So the F2 term is zero.
Packit 67cb25
       */
Packit 67cb25
      F2.val = 0.0;
Packit 67cb25
      F2.err = 0.0;
Packit 67cb25
    } /* end F2 evaluation */
Packit 67cb25
Packit 67cb25
    sgn_2 = ( GSL_IS_ODD(intd) ? -1.0 : 1.0 );
Packit 67cb25
    result->val  = F1.val + sgn_2 * F2.val;
Packit 67cb25
    result->err  = F1.err + F2. err;
Packit 67cb25
    result->err += 2.0 * GSL_DBL_EPSILON * (fabs(F1.val) + fabs(F2.val));
Packit 67cb25
    result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
    return stat_F2;
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    /* d not an integer */
Packit 67cb25
Packit 67cb25
    gsl_sf_result pre1, pre2;
Packit 67cb25
    double sgn1, sgn2;
Packit 67cb25
    gsl_sf_result F1, F2;
Packit 67cb25
    int status_F1, status_F2;
Packit 67cb25
Packit 67cb25
    /* These gamma functions appear in the denominator, so we
Packit 67cb25
     * catch their harmless domain errors and set the terms to zero.
Packit 67cb25
     */
Packit 67cb25
    gsl_sf_result ln_g1ca,  ln_g1cb,  ln_g2a,  ln_g2b;
Packit 67cb25
    double sgn_g1ca, sgn_g1cb, sgn_g2a, sgn_g2b;
Packit 67cb25
    int stat_1ca = gsl_sf_lngamma_sgn_e(c-a, &ln_g1ca, &sgn_g1ca);
Packit 67cb25
    int stat_1cb = gsl_sf_lngamma_sgn_e(c-b, &ln_g1cb, &sgn_g1cb);
Packit 67cb25
    int stat_2a  = gsl_sf_lngamma_sgn_e(a, &ln_g2a, &sgn_g2a);
Packit 67cb25
    int stat_2b  = gsl_sf_lngamma_sgn_e(b, &ln_g2b, &sgn_g2b);
Packit 67cb25
    int ok1 = (stat_1ca == GSL_SUCCESS && stat_1cb == GSL_SUCCESS);
Packit 67cb25
    int ok2 = (stat_2a  == GSL_SUCCESS && stat_2b  == GSL_SUCCESS);
Packit 67cb25
    
Packit 67cb25
    gsl_sf_result ln_gc,  ln_gd,  ln_gmd;
Packit 67cb25
    double sgn_gc, sgn_gd, sgn_gmd;
Packit 67cb25
    gsl_sf_lngamma_sgn_e( c, &ln_gc,  &sgn_gc);
Packit 67cb25
    gsl_sf_lngamma_sgn_e( d, &ln_gd,  &sgn_gd);
Packit 67cb25
    gsl_sf_lngamma_sgn_e(-d, &ln_gmd, &sgn_gmd);
Packit 67cb25
    
Packit 67cb25
    sgn1 = sgn_gc * sgn_gd  * sgn_g1ca * sgn_g1cb;
Packit 67cb25
    sgn2 = sgn_gc * sgn_gmd * sgn_g2a  * sgn_g2b;
Packit 67cb25
Packit 67cb25
    if(ok1 && ok2) {
Packit 67cb25
      double ln_pre1_val = ln_gc.val + ln_gd.val  - ln_g1ca.val - ln_g1cb.val;
Packit 67cb25
      double ln_pre2_val = ln_gc.val + ln_gmd.val - ln_g2a.val  - ln_g2b.val + d*log(1.0-x);
Packit 67cb25
      double ln_pre1_err = ln_gc.err + ln_gd.err + ln_g1ca.err + ln_g1cb.err;
Packit 67cb25
      double ln_pre2_err = ln_gc.err + ln_gmd.err + ln_g2a.err  + ln_g2b.err;
Packit 67cb25
      if(ln_pre1_val < GSL_LOG_DBL_MAX && ln_pre2_val < GSL_LOG_DBL_MAX) {
Packit 67cb25
        gsl_sf_exp_err_e(ln_pre1_val, ln_pre1_err, &pre1;;
Packit 67cb25
        gsl_sf_exp_err_e(ln_pre2_val, ln_pre2_err, &pre2;;
Packit 67cb25
        pre1.val *= sgn1;
Packit 67cb25
        pre2.val *= sgn2;
Packit 67cb25
      }
Packit 67cb25
      else {
Packit 67cb25
        OVERFLOW_ERROR(result);
Packit 67cb25
      }
Packit 67cb25
    }
Packit 67cb25
    else if(ok1 && !ok2) {
Packit 67cb25
      double ln_pre1_val = ln_gc.val + ln_gd.val - ln_g1ca.val - ln_g1cb.val;
Packit 67cb25
      double ln_pre1_err = ln_gc.err + ln_gd.err + ln_g1ca.err + ln_g1cb.err;
Packit 67cb25
      if(ln_pre1_val < GSL_LOG_DBL_MAX) {
Packit 67cb25
        gsl_sf_exp_err_e(ln_pre1_val, ln_pre1_err, &pre1;;
Packit 67cb25
        pre1.val *= sgn1;
Packit 67cb25
        pre2.val = 0.0;
Packit 67cb25
        pre2.err = 0.0;
Packit 67cb25
      }
Packit 67cb25
      else {
Packit 67cb25
        OVERFLOW_ERROR(result);
Packit 67cb25
      }
Packit 67cb25
    }
Packit 67cb25
    else if(!ok1 && ok2) {
Packit 67cb25
      double ln_pre2_val = ln_gc.val + ln_gmd.val - ln_g2a.val - ln_g2b.val + d*log(1.0-x);
Packit 67cb25
      double ln_pre2_err = ln_gc.err + ln_gmd.err + ln_g2a.err + ln_g2b.err;
Packit 67cb25
      if(ln_pre2_val < GSL_LOG_DBL_MAX) {
Packit 67cb25
        pre1.val = 0.0;
Packit 67cb25
        pre1.err = 0.0;
Packit 67cb25
        gsl_sf_exp_err_e(ln_pre2_val, ln_pre2_err, &pre2;;
Packit 67cb25
        pre2.val *= sgn2;
Packit 67cb25
      }
Packit 67cb25
      else {
Packit 67cb25
        OVERFLOW_ERROR(result);
Packit 67cb25
      }
Packit 67cb25
    }
Packit 67cb25
    else {
Packit 67cb25
      pre1.val = 0.0;
Packit 67cb25
      pre2.val = 0.0;
Packit 67cb25
      UNDERFLOW_ERROR(result);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    status_F1 = hyperg_2F1_series(  a,   b, 1.0-d, 1.0-x, &F1;;
Packit 67cb25
    status_F2 = hyperg_2F1_series(c-a, c-b, 1.0+d, 1.0-x, &F2;;
Packit 67cb25
Packit 67cb25
    result->val  = pre1.val*F1.val + pre2.val*F2.val;
Packit 67cb25
    result->err  = fabs(pre1.val*F1.err) + fabs(pre2.val*F2.err);
Packit 67cb25
    result->err += fabs(pre1.err*F1.val) + fabs(pre2.err*F2.val);
Packit 67cb25
    result->err += 2.0 * GSL_DBL_EPSILON * (fabs(pre1.val*F1.val) + fabs(pre2.val*F2.val));
Packit 67cb25
    result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
Packit 67cb25
Packit 67cb25
    if (status_F1)
Packit 67cb25
      return status_F1;
Packit 67cb25
Packit 67cb25
    if (status_F2)
Packit 67cb25
      return status_F2;
Packit 67cb25
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
static int pow_omx(const double x, const double p, gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  double ln_omx;
Packit 67cb25
  double ln_result;
Packit 67cb25
  if(fabs(x) < GSL_ROOT5_DBL_EPSILON) {
Packit 67cb25
    ln_omx = -x*(1.0 + x*(1.0/2.0 + x*(1.0/3.0 + x/4.0 + x*x/5.0)));
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    ln_omx = log(1.0-x);
Packit 67cb25
  }
Packit 67cb25
  ln_result = p * ln_omx;
Packit 67cb25
  return gsl_sf_exp_err_e(ln_result, GSL_DBL_EPSILON * fabs(ln_result), result);
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
/*-*-*-*-*-*-*-*-*-*-*-* Functions with Error Codes *-*-*-*-*-*-*-*-*-*-*-*/
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_sf_hyperg_2F1_e(double a, double b, const double c,
Packit 67cb25
                       const double x,
Packit 67cb25
                       gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  const double d = c - a - b;
Packit 67cb25
  const double rinta = floor(a + 0.5);
Packit 67cb25
  const double rintb = floor(b + 0.5);
Packit 67cb25
  const double rintc = floor(c + 0.5);
Packit 67cb25
  const int a_neg_integer = ( a < 0.0  &&  fabs(a - rinta) < locEPS );
Packit 67cb25
  const int b_neg_integer = ( b < 0.0  &&  fabs(b - rintb) < locEPS );
Packit 67cb25
  const int c_neg_integer = ( c < 0.0  &&  fabs(c - rintc) < locEPS );
Packit 67cb25
Packit 67cb25
  result->val = 0.0;
Packit 67cb25
  result->err = 0.0;
Packit 67cb25
Packit 67cb25
   /* Handle x == 1.0 RJM */
Packit 67cb25
Packit 67cb25
  if (fabs (x - 1.0) < locEPS && (c - a - b) > 0 && c != 0 && !c_neg_integer) {
Packit 67cb25
    gsl_sf_result lngamc, lngamcab, lngamca, lngamcb;
Packit 67cb25
    double lngamc_sgn, lngamca_sgn, lngamcb_sgn;
Packit 67cb25
    int status;
Packit 67cb25
    int stat1 = gsl_sf_lngamma_sgn_e (c, &lngamc, &lngamc_sgn);
Packit 67cb25
    int stat2 = gsl_sf_lngamma_e (c - a - b, &lngamcab);
Packit 67cb25
    int stat3 = gsl_sf_lngamma_sgn_e (c - a, &lngamca, &lngamca_sgn);
Packit 67cb25
    int stat4 = gsl_sf_lngamma_sgn_e (c - b, &lngamcb, &lngamcb_sgn);
Packit 67cb25
    
Packit 67cb25
    if (stat1 != GSL_SUCCESS || stat2 != GSL_SUCCESS
Packit 67cb25
        || stat3 != GSL_SUCCESS || stat4 != GSL_SUCCESS)
Packit 67cb25
      {
Packit 67cb25
        DOMAIN_ERROR (result);
Packit 67cb25
      }
Packit 67cb25
    
Packit 67cb25
    status =
Packit 67cb25
      gsl_sf_exp_err_e (lngamc.val + lngamcab.val - lngamca.val - lngamcb.val,
Packit 67cb25
                        lngamc.err + lngamcab.err + lngamca.err + lngamcb.err,
Packit 67cb25
                        result);
Packit 67cb25
    
Packit 67cb25
    result->val *= lngamc_sgn / (lngamca_sgn * lngamcb_sgn);
Packit 67cb25
      return status;
Packit 67cb25
  }
Packit 67cb25
  
Packit 67cb25
  if(x < -1.0 || 1.0 <= x) {
Packit 67cb25
    DOMAIN_ERROR(result);
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  if(c_neg_integer) {
Packit 67cb25
    /* If c is a negative integer, then either a or b must be a
Packit 67cb25
       negative integer of smaller magnitude than c to ensure
Packit 67cb25
       cancellation of the series. */
Packit 67cb25
    if(! (a_neg_integer && a > c + 0.1) && ! (b_neg_integer && b > c + 0.1)) {
Packit 67cb25
      DOMAIN_ERROR(result);
Packit 67cb25
    }
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  if(fabs(c-b) < locEPS || fabs(c-a) < locEPS) {
Packit 67cb25
    return pow_omx(x, d, result);  /* (1-x)^(c-a-b) */
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  if(a >= 0.0 && b >= 0.0 && c >=0.0 && x >= 0.0 && x < 0.995) {
Packit 67cb25
    /* Series has all positive definite
Packit 67cb25
     * terms and x is not close to 1.
Packit 67cb25
     */
Packit 67cb25
    return hyperg_2F1_series(a, b, c, x, result);
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  if(fabs(a) < 10.0 && fabs(b) < 10.0) {
Packit 67cb25
    /* a and b are not too large, so we attempt
Packit 67cb25
     * variations on the series summation.
Packit 67cb25
     */
Packit 67cb25
    if(a_neg_integer) {
Packit 67cb25
      return hyperg_2F1_series(rinta, b, c, x, result);
Packit 67cb25
    }
Packit 67cb25
    if(b_neg_integer) {
Packit 67cb25
      return hyperg_2F1_series(a, rintb, c, x, result);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    if(x < -0.25) {
Packit 67cb25
      return hyperg_2F1_luke(a, b, c, x, result);
Packit 67cb25
    }
Packit 67cb25
    else if(x < 0.5) {
Packit 67cb25
      return hyperg_2F1_series(a, b, c, x, result);
Packit 67cb25
    }
Packit 67cb25
    else {
Packit 67cb25
      if(fabs(c) > 10.0) {
Packit 67cb25
        return hyperg_2F1_series(a, b, c, x, result);
Packit 67cb25
      }
Packit 67cb25
      else {
Packit 67cb25
        return hyperg_2F1_reflect(a, b, c, x, result);
Packit 67cb25
      }
Packit 67cb25
    }
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    /* Either a or b or both large.
Packit 67cb25
     * Introduce some new variables ap,bp so that bp is
Packit 67cb25
     * the larger in magnitude.
Packit 67cb25
     */
Packit 67cb25
    double ap, bp; 
Packit 67cb25
    if(fabs(a) > fabs(b)) {
Packit 67cb25
      bp = a;
Packit 67cb25
      ap = b;
Packit 67cb25
    }
Packit 67cb25
    else {
Packit 67cb25
      bp = b;
Packit 67cb25
      ap = a;
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    if(x < 0.0) {
Packit 67cb25
      /* What the hell, maybe Luke will converge.
Packit 67cb25
       */
Packit 67cb25
      return hyperg_2F1_luke(a, b, c, x, result);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    if(GSL_MAX_DBL(fabs(ap),1.0)*fabs(bp)*fabs(x) < 2.0*fabs(c)) {
Packit 67cb25
      /* If c is large enough or x is small enough,
Packit 67cb25
       * we can attempt the series anyway.
Packit 67cb25
       */
Packit 67cb25
      return hyperg_2F1_series(a, b, c, x, result);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    if(fabs(bp*bp*x*x) < 0.001*fabs(bp) && fabs(ap) < 10.0) {
Packit 67cb25
      /* The famous but nearly worthless "large b" asymptotic.
Packit 67cb25
       */
Packit 67cb25
      int stat = gsl_sf_hyperg_1F1_e(ap, c, bp*x, result);
Packit 67cb25
      result->err = 0.001 * fabs(result->val);
Packit 67cb25
      return stat;
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    /* We give up. */
Packit 67cb25
    result->val = 0.0;
Packit 67cb25
    result->err = 0.0;
Packit 67cb25
    GSL_ERROR ("error", GSL_EUNIMPL);
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_sf_hyperg_2F1_conj_e(const double aR, const double aI, const double c,
Packit 67cb25
                            const double x,
Packit 67cb25
                            gsl_sf_result * result)
Packit 67cb25
{
Packit 67cb25
  const double ax = fabs(x);
Packit 67cb25
  const double rintc = floor(c + 0.5);
Packit 67cb25
  const int c_neg_integer = ( c < 0.0  &&  fabs(c - rintc) < locEPS );
Packit 67cb25
Packit 67cb25
  result->val = 0.0;
Packit 67cb25
  result->err = 0.0;
Packit 67cb25
Packit 67cb25
  if(ax >= 1.0 || c_neg_integer || c == 0.0) {
Packit 67cb25
    DOMAIN_ERROR(result);
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  if(   (ax < 0.25 && fabs(aR) < 20.0 && fabs(aI) < 20.0)
Packit 67cb25
     || (c > 0.0 && x > 0.0)
Packit 67cb25
    ) {
Packit 67cb25
    return hyperg_2F1_conj_series(aR, aI, c, x, result);
Packit 67cb25
  }
Packit 67cb25
  else if(fabs(aR) < 10.0 && fabs(aI) < 10.0) {
Packit 67cb25
    if(x < -0.25) {
Packit 67cb25
      return hyperg_2F1_conj_luke(aR, aI, c, x, result);
Packit 67cb25
    }
Packit 67cb25
    else {
Packit 67cb25
      return hyperg_2F1_conj_series(aR, aI, c, x, result);
Packit 67cb25
    }
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    if(x < 0.0) {
Packit 67cb25
      /* What the hell, maybe Luke will converge.
Packit 67cb25
       */
Packit 67cb25
      return hyperg_2F1_conj_luke(aR, aI, c, x, result); 
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
    /* Give up. */
Packit 67cb25
    result->val = 0.0;
Packit 67cb25
    result->err = 0.0;
Packit 67cb25
    GSL_ERROR ("error", GSL_EUNIMPL);
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_sf_hyperg_2F1_renorm_e(const double a, const double b, const double c,
Packit 67cb25
                              const double x,
Packit 67cb25
                              gsl_sf_result * result
Packit 67cb25
                              )
Packit 67cb25
{
Packit 67cb25
  const double rinta = floor(a + 0.5);
Packit 67cb25
  const double rintb = floor(b + 0.5);
Packit 67cb25
  const double rintc = floor(c + 0.5);
Packit 67cb25
  const int a_neg_integer = ( a < 0.0  &&  fabs(a - rinta) < locEPS );
Packit 67cb25
  const int b_neg_integer = ( b < 0.0  &&  fabs(b - rintb) < locEPS );
Packit 67cb25
  const int c_neg_integer = ( c < 0.0  &&  fabs(c - rintc) < locEPS );
Packit 67cb25
  
Packit 67cb25
  if(c_neg_integer) {
Packit 67cb25
    if((a_neg_integer && a > c+0.1) || (b_neg_integer && b > c+0.1)) {
Packit 67cb25
      /* 2F1 terminates early */
Packit 67cb25
      result->val = 0.0;
Packit 67cb25
      result->err = 0.0;
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
    else {
Packit 67cb25
      /* 2F1 does not terminate early enough, so something survives */
Packit 67cb25
      /* [Abramowitz+Stegun, 15.1.2] */
Packit 67cb25
      gsl_sf_result g1, g2, g3, g4, g5;
Packit 67cb25
      double s1, s2, s3, s4, s5;
Packit 67cb25
      int stat = 0;
Packit 67cb25
      stat += gsl_sf_lngamma_sgn_e(a-c+1, &g1, &s1;;
Packit 67cb25
      stat += gsl_sf_lngamma_sgn_e(b-c+1, &g2, &s2;;
Packit 67cb25
      stat += gsl_sf_lngamma_sgn_e(a, &g3, &s3;;
Packit 67cb25
      stat += gsl_sf_lngamma_sgn_e(b, &g4, &s4;;
Packit 67cb25
      stat += gsl_sf_lngamma_sgn_e(-c+2, &g5, &s5;;
Packit 67cb25
      if(stat != 0) {
Packit 67cb25
        DOMAIN_ERROR(result);
Packit 67cb25
      }
Packit 67cb25
      else {
Packit 67cb25
        gsl_sf_result F;
Packit 67cb25
        int stat_F = gsl_sf_hyperg_2F1_e(a-c+1, b-c+1, -c+2, x, &F);
Packit 67cb25
        double ln_pre_val = g1.val + g2.val - g3.val - g4.val - g5.val;
Packit 67cb25
        double ln_pre_err = g1.err + g2.err + g3.err + g4.err + g5.err;
Packit 67cb25
        double sg  = s1 * s2 * s3 * s4 * s5;
Packit 67cb25
        int stat_e = gsl_sf_exp_mult_err_e(ln_pre_val, ln_pre_err,
Packit 67cb25
                                              sg * F.val, F.err,
Packit 67cb25
                                              result);
Packit 67cb25
        return GSL_ERROR_SELECT_2(stat_e, stat_F);
Packit 67cb25
      }
Packit 67cb25
    }
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    /* generic c */
Packit 67cb25
    gsl_sf_result F;
Packit 67cb25
    gsl_sf_result lng;
Packit 67cb25
    double sgn;
Packit 67cb25
    int stat_g = gsl_sf_lngamma_sgn_e(c, &lng, &sgn);
Packit 67cb25
    int stat_F = gsl_sf_hyperg_2F1_e(a, b, c, x, &F);
Packit 67cb25
    int stat_e = gsl_sf_exp_mult_err_e(-lng.val, lng.err,
Packit 67cb25
                                          sgn*F.val, F.err,
Packit 67cb25
                                          result);
Packit 67cb25
    return GSL_ERROR_SELECT_3(stat_e, stat_F, stat_g);
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_sf_hyperg_2F1_conj_renorm_e(const double aR, const double aI, const double c,
Packit 67cb25
                                   const double x,
Packit 67cb25
                                   gsl_sf_result * result
Packit 67cb25
                                   )
Packit 67cb25
{
Packit 67cb25
  const double rintc = floor(c  + 0.5);
Packit 67cb25
  const double rinta = floor(aR + 0.5);
Packit 67cb25
  const int a_neg_integer = ( aR < 0.0 && fabs(aR-rinta) < locEPS && aI == 0.0);
Packit 67cb25
  const int c_neg_integer = (  c < 0.0 && fabs(c - rintc) < locEPS );
Packit 67cb25
Packit 67cb25
  if(c_neg_integer) {
Packit 67cb25
    if(a_neg_integer && aR > c+0.1) {
Packit 67cb25
      /* 2F1 terminates early */
Packit 67cb25
      result->val = 0.0;
Packit 67cb25
      result->err = 0.0;
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
    else {
Packit 67cb25
      /* 2F1 does not terminate early enough, so something survives */
Packit 67cb25
      /* [Abramowitz+Stegun, 15.1.2] */
Packit 67cb25
      gsl_sf_result g1, g2;
Packit 67cb25
      gsl_sf_result g3;
Packit 67cb25
      gsl_sf_result a1, a2;
Packit 67cb25
      int stat = 0;
Packit 67cb25
      stat += gsl_sf_lngamma_complex_e(aR-c+1, aI, &g1, &a1;;
Packit 67cb25
      stat += gsl_sf_lngamma_complex_e(aR, aI, &g2, &a2;;
Packit 67cb25
      stat += gsl_sf_lngamma_e(-c+2.0, &g3;;
Packit 67cb25
      if(stat != 0) {
Packit 67cb25
        DOMAIN_ERROR(result);
Packit 67cb25
      }
Packit 67cb25
      else {
Packit 67cb25
        gsl_sf_result F;
Packit 67cb25
        int stat_F = gsl_sf_hyperg_2F1_conj_e(aR-c+1, aI, -c+2, x, &F);
Packit 67cb25
        double ln_pre_val = 2.0*(g1.val - g2.val) - g3.val;
Packit 67cb25
        double ln_pre_err = 2.0 * (g1.err + g2.err) + g3.err;
Packit 67cb25
        int stat_e = gsl_sf_exp_mult_err_e(ln_pre_val, ln_pre_err,
Packit 67cb25
                                              F.val, F.err,
Packit 67cb25
                                              result);
Packit 67cb25
        return GSL_ERROR_SELECT_2(stat_e, stat_F);
Packit 67cb25
      }
Packit 67cb25
    }
Packit 67cb25
  }
Packit 67cb25
  else {
Packit 67cb25
    /* generic c */
Packit 67cb25
    gsl_sf_result F;
Packit 67cb25
    gsl_sf_result lng;
Packit 67cb25
    double sgn;
Packit 67cb25
    int stat_g = gsl_sf_lngamma_sgn_e(c, &lng, &sgn);
Packit 67cb25
    int stat_F = gsl_sf_hyperg_2F1_conj_e(aR, aI, c, x, &F);
Packit 67cb25
    int stat_e = gsl_sf_exp_mult_err_e(-lng.val, lng.err,
Packit 67cb25
                                          sgn*F.val, F.err,
Packit 67cb25
                                          result);
Packit 67cb25
    return GSL_ERROR_SELECT_3(stat_e, stat_F, stat_g);
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_hyperg_2F1(double a, double b, double c, double x)
Packit 67cb25
{
Packit 67cb25
  EVAL_RESULT(gsl_sf_hyperg_2F1_e(a, b, c, x, &result));
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double gsl_sf_hyperg_2F1_conj(double aR, double aI, double c, double x)
Packit 67cb25
{
Packit 67cb25
  EVAL_RESULT(gsl_sf_hyperg_2F1_conj_e(aR, aI, c, x, &result));
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double gsl_sf_hyperg_2F1_renorm(double a, double b, double c, double x)
Packit 67cb25
{
Packit 67cb25
  EVAL_RESULT(gsl_sf_hyperg_2F1_renorm_e(a, b, c, x, &result));
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double gsl_sf_hyperg_2F1_conj_renorm(double aR, double aI, double c, double x)
Packit 67cb25
{
Packit 67cb25
  EVAL_RESULT(gsl_sf_hyperg_2F1_conj_renorm_e(aR, aI, c, x, &result));
Packit 67cb25
}