/*******************************************************************************
*
* This file is part of the General Hidden Markov Model Library,
* GHMM version __VERSION__, see http://ghmm.org
*
* Filename: ghmm/ghmm/randvar.c
* Authors: Bernhard Knab, Benjamin Rich, Janne Grunau
*
* Copyright (C) 1998-2004 Alexander Schliep
* Copyright (C) 1998-2001 ZAIK/ZPR, Universitaet zu Koeln
* Copyright (C) 2002-2004 Max-Planck-Institut fuer Molekulare Genetik,
* Berlin
*
* Contact: schliep@ghmm.org
*
* This library is free software; you can redistribute it and/or
* modify it under the terms of the GNU Library General Public
* License as published by the Free Software Foundation; either
* version 2 of the License, or (at your option) any later version.
*
* This library is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* Library General Public License for more details.
*
* You should have received a copy of the GNU Library General Public
* License along with this library; if not, write to the Free
* Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
*
*
* This file is version $Revision$
* from $Date$
* last change by $Author$.
*
*******************************************************************************/
#ifdef HAVE_CONFIG_H
# include "../config.h"
#endif
#include <math.h>
#include <float.h>
#include <string.h>
#ifdef HAVE_LIBPTHREAD
# include <pthread.h>
#endif /* HAVE_LIBPTHREAD */
#include "ghmm.h"
#include "mes.h"
#include "mprintf.h"
#include "randvar.h"
#include "rng.h"
#ifdef DO_WITH_GSL
# include <gsl/gsl_math.h>
# include <gsl/gsl_sf_erf.h>
# include <gsl/gsl_randist.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_blas.h>
#endif /* DO_WITH_GSL */
#if defined(__STDC_VERSION__) && (__STDC_VERSION__ >= 199901L)
#else
static double ighmm_erf (double x);
static double ighmm_erfc (double x);
#endif /* check for ISO C99 */
/* A list of already calculated values of the density function of a
N(0,1)-distribution, with x in [0.00, 19.99] */
#define PDFLEN 2000
#define X_STEP_PDF 0.01 /* step size */
#define X_FAKT_PDF 100 /* equivalent to step size */
static double pdf_stdnormal[PDFLEN];
static int pdf_stdnormal_exists = 0;
/* A list of already calulated values PHI of the Gauss distribution is
read in, x in [-9.999, 0] */
#define X_STEP_PHI 0.001 /* step size */
#define X_FAKT_PHI 1000 /* equivalent to step size */
static double x_PHI_1 = -1.0;
#ifndef M_SQRT1_2
#define M_SQRT1_2 0.70710678118654752440084436210
#endif
/*============================================================================*/
/* needed by ighmm_gtail_pmue_interpol */
double ighmm_rand_get_xfaktphi ()
{
return X_FAKT_PHI;
}
double ighmm_rand_get_xstepphi ()
{
return X_STEP_PHI;
}
double ighmm_rand_get_philen ()
{
#ifdef DO_WITH_GSL
return 0 /*PHI_len*/;
#else
return ighmm_rand_get_xPHIless1 () / X_STEP_PHI;
#endif
}
/*============================================================================*/
double ighmm_rand_get_PHI (double x)
{
#if defined(__STDC_VERSION__) && (__STDC_VERSION__ >= 199901L)
return (erf (x * M_SQRT1_2) + 1.0) / 2.0;
#else
return (ighmm_erf (x * M_SQRT1_2) + 1.0) / 2.0;
#endif
} /* randvar_get_PHI */
/*============================================================================*/
/* When is PHI[x,0,1] == 1? */
double ighmm_rand_get_xPHIless1 ()
{
# define CUR_PROC "ighmm_rand_get_xPHIless1"
if (x_PHI_1 == -1) {
double low, up, half;
low = 0;
up = 100;
while (up - low > 0.001) {
half = (low + up) / 2.0;
if (ighmm_rand_get_PHI (half) < 1.0)
low = half;
else
up = half;
}
x_PHI_1 = low;
}
return (x_PHI_1);
# undef CUR_PROC
}
/*============================================================================*/
double ighmm_rand_get_1overa (double x, double mean, double u)
{
/* Calulates 1/a(x, mean, u), with a = the integral from x til \infty over
the Gauss density function */
# define CUR_PROC "ighmm_rand_get_1overa"
double erfc_value;
if (u <= 0.0) {
GHMM_LOG(LCONVERTED, "u <= 0.0 not allowed\n");
goto STOP;
}
#if defined(__STDC_VERSION__) && (__STDC_VERSION__ >= 199901L)
erfc_value = erfc ((x - mean) / sqrt (u * 2));
#else
erfc_value = ighmm_erfc ((x - mean) / sqrt (u * 2));
#endif
if (erfc_value <= DBL_MIN) {
ighmm_mes (MES_WIN, "a ~= 0.0 critical! (mue = %.2f, u =%.2f)\n", mean, u);
return (erfc_value);
}
else
return (2.0 / erfc_value);
STOP:
return (-1.0);
# undef CUR_PROC
} /* ighmm_rand_get_1overa */
/*============================================================================*/
/* REMARK:
The calulation of this density function was testet, by calculating the
following integral sum for arbitrary mue and u:
for (x = 0, x < ..., x += step(=0.01/0.001/0.0001))
isum += step * ighmm_rand_normal_density_pos(x, mue, u);
In each case, the sum "converged" evidently towards 1!
(BK, 14.6.99)
CHANGE:
Truncate at -EPS_NDT (const.h), so that x = 0 doesn't lead to a problem.
(BK, 15.3.2000)
*/
double ighmm_rand_normal_density_pos (double x, double mean, double u)
{
# define CUR_PROC "ighmm_rand_normal_density_pos"
return ighmm_rand_normal_density_trunc (x, mean, u, -GHMM_EPS_NDT);
# undef CUR_PROC
} /* double ighmm_rand_normal_density_pos */
/*============================================================================*/
double ighmm_rand_normal_density_trunc(double x, double mean, double u,
double a)
{
# define CUR_PROC "ighmm_rand_normal_density_trunc"
#ifndef DO_WITH_GSL
double c;
#endif /* DO_WITH_GSL */
if (u <= 0.0) {
GHMM_LOG(LERROR, "u <= 0.0 not allowed");
goto STOP;
}
if (x < a)
return 0.0;
#ifdef DO_WITH_GSL
/* move mean to the right position */
return gsl_ran_gaussian_tail_pdf(x - mean, a - mean, sqrt(u));
#else
if ((c = ighmm_rand_get_1overa(a, mean, u)) == -1) {
GHMM_LOG_QUEUED(LERROR);
goto STOP;
};
return c * ighmm_rand_normal_density(x, mean, u);
#endif /* DO_WITH_GSL */
STOP:
return -1.0;
# undef CUR_PROC
} /* double ighmm_rand_normal_density_trunc */
/*============================================================================*/
double ighmm_rand_normal_density (double x, double mean, double u)
{
# define CUR_PROC "ighmm_rand_normal_density"
#ifndef DO_WITH_GSL
double expo;
#endif
if (u <= 0.0) {
GHMM_LOG(LCONVERTED, "u <= 0.0 not allowed\n");
goto STOP;
}
/* The denominator is possibly < EPS??? Check that ? */
#ifdef DO_WITH_GSL
/* double gsl_ran_gaussian_pdf (double x, double sigma) */
return gsl_ran_gaussian_pdf (x - mean, sqrt (u));
#else
expo = exp (-1 * m_sqr (mean - x) / (2 * u));
return (1 / (sqrt (2 * PI * u)) * expo);
#endif
STOP:
return (-1.0);
# undef CUR_PROC
} /* double ighmm_rand_normal_density */
/*============================================================================*/
/* covariance matrix is linearized */
double ighmm_rand_binormal_density(const double *x, double *mean, double *cov)
{
# define CUR_PROC "ighmm_rand_binormal_density"
double rho;
#ifndef DO_WITH_GSL
double numerator,part1,part2,part3;
#endif
if (cov[0] <= 0.0 || cov[2 + 1] <= 0.0) {
GHMM_LOG(LCONVERTED, "variance <= 0.0 not allowed\n");
goto STOP;
}
rho = cov[1] / ( sqrt (cov[0]) * sqrt (cov[2 + 1]) );
/* The denominator is possibly < EPS??? Check that ? */
#ifdef DO_WITH_GSL
/* double gsl_ran_bivariate_gaussian_pdf (double x, double y, double sigma_x,
double sigma_y, double rho) */
return gsl_ran_bivariate_gaussian_pdf (x[0], x[1], sqrt (cov[0]),
sqrt (cov[2 + 1]), rho);
#else
part1 = (x[0] - mean[0]) / sqrt (cov[0]);
part2 = (x[1] - mean[1]) / sqrt (cov[2 + 1]);
part3 = m_sqr (part1) - 2 * part1 * part2 + m_sqr (part2);
numerator = exp ( -1 * (part3) / ( 2 * (1 - m_sqr(rho)) ) );
return (numerator / ( 2 * PI * sqrt(1 - m_sqr(rho)) ));
#endif
STOP:
return (-1.0);
# undef CUR_PROC
} /* double ighmm_rand_binormal_density */
/*============================================================================*/
/* matrices are linearized */
double ighmm_rand_multivariate_normal_density(int length, const double *x, double *mean, double *sigmainv, double det)
{
# define CUR_PROC "ighmm_rand_multivariate_normal_density"
/* multivariate normal density function */
/*
* length dimension of the random vetor
* x point at which to evaluate the pdf
* mean vector of means of size n
* sigmainv inverse variance matrix of dimension n x n
* det determinant of covariance matrix
*/
#ifdef DO_WITH_GSL
int i, j;
double ax,ay;
gsl_vector *ym, *xm, *gmean;
gsl_matrix *inv = gsl_matrix_alloc(length, length);
for (i=0; i<length; ++i) {
for (j=0; j<length; ++j) {
gsl_matrix_set(inv, i, j, sigmainv[i*length+j]);
}
}
xm = gsl_vector_alloc(length);
gmean = gsl_vector_alloc(length);
/*gsl_vector_memcpy(xm, x);*/
for (i=0; i<length; ++i) {
gsl_vector_set(xm, i, x[i]);
gsl_vector_set(gmean, i, mean[i]);
}
gsl_vector_sub(xm, gmean);
ym = gsl_vector_alloc(length);
gsl_blas_dsymv(CblasUpper, 1.0, inv, xm, 0.0, ym);
gsl_matrix_free(inv);
gsl_blas_ddot(xm, ym, &ay);
gsl_vector_free(xm);
gsl_vector_free(ym);
ay = exp(-0.5*ay) / sqrt(pow((2*M_PI), length) * det);
return ay;
#else
/* do without GSL */
int i, j;
double ay, tempv;
ay = 0;
for (i=0; i<length; ++i) {
tempv = 0;
for (j=0; j<length; ++j) {
tempv += (x[j]-mean[j])*sigmainv[j*length+i];
}
ay += tempv*(x[i]-mean[i]);
}
ay = exp(-0.5*ay) / sqrt(pow((2*PI), length) * det);
return ay;
#endif
# undef CUR_PROC
} /* double ighmm_rand_multivariate_normal_density */
/*============================================================================*/
double ighmm_rand_uniform_density (double x, double max, double min)
{
# define CUR_PROC "ighmm_rand_uniform_density"
double prob;
if (max <= min) {
GHMM_LOG(LCONVERTED, "max <= min not allowed \n");
goto STOP;
}
prob = 1.0/(max-min);
if ( (x <= max) && (x>=min) ){
return prob;
}else{
return 0.0;
}
STOP:
return (-1.0);
# undef CUR_PROC
} /* double ighmm_rand_uniform_density */
/*============================================================================*/
/* special ghmm_cmodel pdf need it: smo->density==normal_approx: */
/* generates a table of of aequidistant samples of gaussian pdf */
static int randvar_init_pdf_stdnormal ()
{
# define CUR_PROC "randvar_init_pdf_stdnormal"
int i;
double x = 0.00;
for (i = 0; i < PDFLEN; i++) {
pdf_stdnormal[i] = 1 / (sqrt (2 * PI)) * exp (-1 * x * x / 2);
x += (double) X_STEP_PDF;
}
pdf_stdnormal_exists = 1;
/* printf("pdf_stdnormal_exists = %d\n", pdf_stdnormal_exists); */
return (0);
# undef CUR_PROC
} /* randvar_init_pdf_stdnormal */
double ighmm_rand_normal_density_approx (double x, double mean, double u)
{
# define CUR_PROC "ighmm_rand_normal_density_approx"
#ifdef HAVE_LIBPTHREAD
static pthread_mutex_t lock;
#endif /* HAVE_LIBPTHREAD */
int i;
double y, z, pdf_x;
if (u <= 0.0) {
GHMM_LOG(LCONVERTED, "u <= 0.0 not allowed\n");
goto STOP;
}
if (!pdf_stdnormal_exists) {
#ifdef HAVE_LIBPTHREAD
pthread_mutex_lock (&lock); /* Put on a lock, because the clustering is parallel */
#endif /* HAVE_LIBPTHREAD */
randvar_init_pdf_stdnormal ();
#ifdef HAVE_LIBPTHREAD
pthread_mutex_unlock (&lock); /* Take the lock off */
#endif /* HAVE_LIBPTHREAD */
}
y = 1 / sqrt (u);
z = fabs ((x - mean) * y);
i = (int) (z * X_FAKT_PDF);
/* linear interpolation: */
if (i >= PDFLEN - 1) {
i = PDFLEN - 1;
pdf_x = y * pdf_stdnormal[i];
}
else
pdf_x = y * (pdf_stdnormal[i] +
(z - i * X_STEP_PDF) *
(pdf_stdnormal[i + 1] - pdf_stdnormal[i]) / X_STEP_PDF);
return (pdf_x);
STOP:
return (-1.0);
# undef CUR_PROC
} /* double ighmm_rand_normal_density_approx */
double ighmm_rand_dirichlet(int seed, int len, double *alpha, double *theta){
if (seed != 0) {
GHMM_RNG_SET(RNG, seed);
}
#ifdef DO_WITH_GSL
gsl_ran_dirichlet(RNG, len, alpha, theta);
#else
printf("not implemted without gsl. Compile with gsl to use dirichlet");
return 0;
#endif
}
double ighmm_rand_gamma(double a, double b, int seed){
if (seed != 0) {
GHMM_RNG_SET(RNG, seed);
}
#ifdef DO_WITH_GSL
return (gsl_ran_gamma_knuth(RNG, a, b));
#else
printf("not implemted without gsl. Compile with gsl to use gamma");
return 0;
#endif
}
/*============================================================================*/
double ighmm_rand_std_normal (int seed)
{
# define CUR_PROC "ighmm_rand_std_normal"
if (seed != 0) {
GHMM_RNG_SET (RNG, seed);
}
#ifdef DO_WITH_GSL
return (gsl_ran_gaussian (RNG, 1.0));
#else
/* Use the polar Box-Mueller transform */
/*
double x, y, r2;
do {
x = 2.0 * GHMM_RNG_UNIFORM(RNG) - 1.0;
y = 2.0 * GHMM_RNG_UNIFORM(RNG) - 1.0;
r2 = (x * x) + (y * y);
} while (r2 >= 1.0);
return x * sqrt((-2.0 * log(r2)) / r2);
*/
double r2, theta;
r2 = -2.0 * log (GHMM_RNG_UNIFORM (RNG)); /* r2 ~ chi-square(2) */
theta = 2.0 * PI * GHMM_RNG_UNIFORM (RNG); /* theta ~ uniform(0, 2 \pi) */
return sqrt (r2) * cos (theta);
#endif
# undef CUR_PROC
} /* ighmm_rand_std_normal */
/*============================================================================*/
double ighmm_rand_normal(double mue, double u, int seed)
{
# define CUR_POC "ighmm_rand_normal"
if (seed != 0) {
GHMM_RNG_SET(RNG, seed);
}
#ifdef DO_WITH_GSL
double x = gsl_ran_gaussian(RNG, sqrt (u)) + mue ;
return x;
#else
double x;
x = sqrt(u) * ighmm_rand_std_normal(seed) + mue;
return x;
#endif
# undef CUR_PROC
} /* ighmm_rand_normal */
/*============================================================================*/
int ighmm_rand_multivariate_normal (int dim, double *x, double *mue, double *sigmacd, int seed)
{
# define CUR_PROC "ighmm_rand_multivariate_normal"
/* generate random vector of multivariate normal
*
* dim number of dimensions
* x space to store resulting vector in
* mue vector of means
* sigmacd linearized cholesky decomposition of cov matrix
* seed RNG seed
*
* see Barr & Slezak, A Comparison of Multivariate Normal Generators */
int i, j;
#ifdef DO_WITH_GSL
gsl_vector *y = gsl_vector_alloc(dim);
gsl_vector *xgsl = gsl_vector_alloc(dim);
gsl_matrix *cd = gsl_matrix_alloc(dim, dim);
#endif
if (seed != 0) {
GHMM_RNG_SET (RNG, seed);
/* do something here */
return 0;
}
else {
#ifdef DO_WITH_GSL
/* cholesky decomposition matrix */
for (i=0;i<dim;i++) {
for (j=0;j<dim;j++) {
gsl_matrix_set(cd, i, j, sigmacd[i*dim+j]);
}
}
/* generate a random vector N(O,I) */
for (i=0;i<dim;i++) {
gsl_vector_set(y, i, ighmm_rand_std_normal(seed));
}
/* multiply cd with y */
gsl_blas_dgemv(CblasNoTrans, 1.0, cd, y, 0.0, xgsl);
for (i=0;i<dim;i++) {
x[i] = gsl_vector_get(xgsl, i) + mue[i];
}
gsl_vector_free(y);
gsl_vector_free(xgsl);
gsl_matrix_free(cd);
#else
/* multivariate random numbers without gsl */
double randuni;
for (i=0;i<dim;i++) {
randuni = ighmm_rand_std_normal(seed);
for (j=0;j<dim;j++) {
if (i==0)
x[j] = mue[j];
x[j] += randuni * sigmacd[j*dim+i];
}
}
#endif
return 0;
}
# undef CUR_PROC
} /* ighmm_rand_multivariate_normal */
/*============================================================================*/
#ifndef DO_WITH_GSL
# define C0 2.515517
# define C1 0.802853
# define C2 0.010328
# define D1 1.432788
# define D2 0.189269
# define D3 0.001308
#endif
double ighmm_rand_normal_right (double a, double mue, double u, int seed)
{
# define CUR_PROC "ighmm_rand_normal_right"
double x = -1;
double sigma;
#ifdef DO_WITH_GSL
double s;
#else
double U, Us, Us1, Feps, t, T;
#endif
if (u <= 0.0) {
GHMM_LOG(LCONVERTED, "u <= 0.0 not allowed\n");
goto STOP;
}
sigma = sqrt(u);
if (seed != 0) {
GHMM_RNG_SET (RNG, seed);
}
#ifdef DO_WITH_GSL
/* move boundary to lower values in order to achieve maximum at mue
gsl_ran_gaussian_tail(generator, lower_boundary, sigma)
*/
return mue + gsl_ran_gaussian_tail(RNG, a - mue, sqrt (u));
#else /* DO_WITH_GSL */
/* Inverse transformation with restricted sampling by Fishman */
U = GHMM_RNG_UNIFORM(RNG);
Feps = ighmm_rand_get_PHI((a-mue) / sigma);
Us = Feps + (1-Feps) * U;
Us1 = 1-Us;
t = m_min (Us, Us1);
t = sqrt (-log (t * t));
T =
sigma * (t - (C0 + t * (C1 + t * C2))
/ (1 + t * (D1 + t * (D2 + t * D3))));
if (Us < Us1)
x = mue - T;
else
x = mue + T;
#endif /* DO_WITH_GSL */
STOP:
return x;
# undef CUR_PROC
} /* randvar_normal_pos */
/*============================================================================*/
double ighmm_rand_uniform_int (int seed, int K)
{
# define CUR_PROC "ighmm_rand_uniform_int"
if (seed != 0) {
GHMM_RNG_SET (RNG, seed);
}
#ifdef DO_WITH_GSL
/* more direct solution than old version ! */
return (double) gsl_rng_uniform_int (RNG, K);
#else
return (double) ((int) (((double) K) * GHMM_RNG_UNIFORM (RNG)));
#endif
# undef CUR_PROC
} /* ighmm_rand_uniform_int */
/*===========================================================================*/
double ighmm_rand_uniform_cont (int seed, double max, double min)
{
# define CUR_PROC "ighmm_rand_uniform_cont"
if (max <= min) {
GHMM_LOG(LCONVERTED, "max <= min not allowed\n");
goto STOP;
}
if (seed != 0) {
GHMM_RNG_SET (RNG, seed);
}
#ifdef DO_WITH_GSL
return (double)(((double)gsl_rng_uniform (RNG)*(max-min)) + min);
#else
return (double)((GHMM_RNG_UNIFORM (RNG))*(max-min) + min );
#endif
STOP:
return (-1.0);
# undef CUR_PROC
} /* ighmm_rand_uniform_cont */
/*============================================================================*/
/* cumalative distribution function of N(mean, u) */
double ighmm_rand_normal_cdf (double x, double mean, double u)
{
# define CUR_PROC "ighmm_rand_normal_cdf"
if (u <= 0.0) {
GHMM_LOG(LCONVERTED, "u <= 0.0 not allowed\n");
goto STOP;
}
#if defined(__STDC_VERSION__) && (__STDC_VERSION__ >= 199901L)
/* PHI(x)=erf(x/sqrt(2))/2+0.5 */
return (erf ((x - mean) / sqrt (u * 2.0)) + 1.0) / 2.0;
#else
return (ighmm_erf ((x - mean) / sqrt (u * 2.0)) + 1.0) / 2.0;
#endif /* Check for ISO C99 */
STOP:
return (-1.0);
# undef CUR_PROC
} /* double ighmm_rand_normal_cdf */
/*============================================================================*/
/* cumalative distribution function of a-truncated N(mean, u) */
double ighmm_rand_normal_right_cdf (double x, double mean, double u, double a)
{
# define CUR_PROC "ighmm_rand_normal_right_cdf"
if (x <= a)
return (0.0);
if (u <= a) {
GHMM_LOG(LCONVERTED, "u <= a not allowed\n");
goto STOP;
}
#if defined(__STDC_VERSION__) && (__STDC_VERSION__ >= 199901L)
/*
Function: int erfc (double x, gsl_sf_result * result)
These routines compute the complementary error function
erfc(x) = 1 - erf(x) = 2/\sqrt(\pi) \int_x^\infty \exp(-t^2).
*/
return 1.0 + (erf ((x - mean) / sqrt (u * 2)) -
1.0) / erfc ((a - mean) / sqrt (u * 2));
#else
return 1.0 + (ighmm_erf ((x - mean) / sqrt (u * 2)) -
1.0) / ighmm_erfc ((a - mean) / sqrt (u * 2));
#endif /* Check for ISO C99 */
STOP:
return (-1.0);
# undef CUR_PROC
} /* double ighmm_rand_normal_cdf */
/*============================================================================*/
/* cumalative distribution function of a uniform distribution in the range [min,max] */
double ighmm_rand_uniform_cdf (double x, double max, double min)
{
# define CUR_PROC "ighmm_rand_uniform_cdf"
if (max <= min) {
GHMM_LOG(LCONVERTED, "max <= min not allowed\n");
goto STOP;
}
if (x < min) {
return 0.0;
}
if (x >= max) {
return 1.0;
}
return (x-min)/(max-min);
STOP:
return (-1.0);
# undef CUR_PROC
} /* ighmm_rand_uniform_cdf */
#if defined(__STDC_VERSION__) && (__STDC_VERSION__ >= 199901L)
#else
/* THIS WILL BE OBSOLETE WHEN WE USE ISO C99 */
/*===========================================================================
*
* The following functions for the error function and the complementory
* error function are taken from
* http://www.mathematik.uni-bielefeld.de/~sillke/ALGORITHMS/special-functions/erf.c
* and have following different copyright.
* ====================================================
* Copyright (C) 1993 by Sun Microsystems, Inc. All rights reserved.
*
* Developed at SunPro, a Sun Microsystems, Inc. business.
* Permission to use, copy, modify, and distribute this
* software is freely granted, provided that this notice
* is preserved.
* ====================================================
*/
/*
* ====================================================
*
* Reference:
*
* W. J. Cody,
* Rational Chebychev approximations for the
* error function.
* Mathematics of Computations 23 (1969) 631-637
*
* W. J. Cody,
* Performance evaluations of programs for the
* error and complementary error function.
* Transactions of the ACM on Mathematical Software,
* 16:1 (March 1990) 38-46
*
* W. J. Cody,
* SPECFUN - A portable special function package,
* In: New Computing environments;
* Microcomputers in Large-Scale Scientific Computing,
* A. Wouk, SIAM, 1987, 1-12
*
* W. J. Cody,
* http://www.netlib.org/specfun/erf.
*
* For calculations with complex arguments see:
*
* Walter Gautschi
* "Efficient computation of the complex error function"
* SIAM J. Numer. Anal.
* 7:1 (1970), 187-198
*
* J.A.C. Weideman
* "Computation of the complex error function"
* SIAM J. Numer. Anal.
* 31:5 (1994), 1497-1518
*
* ====================================================
*/
/* double erf(double x)
* double erfc(double x)
* x
* 2 |\
* erf(x) = --------- | exp(-t*t)dt
* sqrt(pi) \|
* 0
*
* erfc(x) = 1-erf(x)
* Note that
* erf(-x) = -erf(x)
* erfc(-x) = 2 - erfc(x)
*
* Method:
* 1. For |x| in [0, 0.84375]
* erf(x) = x + x*R(x^2)
* erfc(x) = 1 - erf(x) if x in [-.84375,0.25]
* = 0.5 + ((0.5-x)-x*R) if x in [0.25,0.84375]
* where R = P/Q where P is an odd poly of degree 8 and
* Q is an odd poly of degree 10.
* -57.90
* | R - (erf(x)-x)/x | <= 2
*
*
* Remark. The formula is derived by noting
* erf(x) = (2/sqrt(pi))*(x - x^3/3 + x^5/10 - x^7/42 + ....)
* and that
* 2/sqrt(pi) = 1.128379167095512573896158903121545171688
* is close to one. The interval is chosen because the fix
* point of erf(x) is near 0.6174 (i.e., erf(x)=x when x is
* near 0.6174), and by some experiment, 0.84375 is chosen to
* guarantee the error is less than one ulp for erf.
*
* 2. For |x| in [0.84375,1.25], let s = |x| - 1, and
* c = 0.84506291151 rounded to single (24 bits)
* erf(x) = sign(x) * (c + P1(s)/Q1(s))
* erfc(x) = (1-c) - P1(s)/Q1(s) if x > 0
* 1+(c+P1(s)/Q1(s)) if x < 0
* |P1/Q1 - (erf(|x|)-c)| <= 2**-59.06
* Remark: here we use the taylor series expansion at x=1.
* erf(1+s) = erf(1) + s*Poly(s)
* = 0.845.. + P1(s)/Q1(s)
* That is, we use rational approximation to approximate
* erf(1+s) - (c = (single)0.84506291151)
* Note that |P1/Q1|< 0.078 for x in [0.84375,1.25]
* where
* P1(s) = degree 6 poly in s
* Q1(s) = degree 6 poly in s
*
* 3. For x in [1.25,1/0.35(~2.857143)],
* erfc(x) = (1/x)*exp(-x*x-0.5625+R1/S1)
* erf(x) = 1 - erfc(x)
* where
* R1(z) = degree 7 poly in z, (z=1/x^2)
* S1(z) = degree 8 poly in z
*
* 4. For x in [1/0.35,28]
* erfc(x) = (1/x)*exp(-x*x-0.5625+R2/S2) if x > 0
* = 2.0 - (1/x)*exp(-x*x-0.5625+R2/S2) if -6<x<0
* = 2.0 - tiny (if x <= -6)
* erf(x) = sign(x)*(1.0 - erfc(x)) if x < 6, else
* erf(x) = sign(x)*(1.0 - tiny)
* where
* R2(z) = degree 6 poly in z, (z=1/x^2)
* S2(z) = degree 7 poly in z
*
* Note1:
* To compute exp(-x*x-0.5625+R/S), let s be a single
* precision number and s := x; then
* -x*x = -s*s + (s-x)*(s+x)
* exp(-x*x-0.5626+R/S) =
* exp(-s*s-0.5625)*exp((s-x)*(s+x)+R/S);
* Note2:
* Here 4 and 5 make use of the asymptotic series
* exp(-x*x)
* erfc(x) ~ ---------- * ( 1 + Poly(1/x^2) )
* x*sqrt(pi)
* We use rational approximation to approximate
* g(s)=f(1/x^2) = log(erfc(x)*x) - x*x + 0.5625
* Here is the error bound for R1/S1 and R2/S2
* |R1/S1 - f(x)| < 2**(-62.57)
* |R2/S2 - f(x)| < 2**(-61.52)
*
* 5. For inf > x >= 28
* erf(x) = sign(x) *(1 - tiny) (raise inexact)
* erfc(x) = tiny*tiny (raise underflow) if x > 0
* = 2 - tiny if x<0
*
* 7. Special case:
* erf(0) = 0, erf(inf) = 1, erf(-inf) = -1,
* erfc(0) = 1, erfc(inf) = 0, erfc(-inf) = 2,
* erfc/erf(NaN) is NaN
*/
static const double
tiny = 1e-300,
half = 5.00000000000000000000e-01, /* 0x3FE00000, 0x00000000 */
one = 1.00000000000000000000e+00, /* 0x3FF00000, 0x00000000 */
two = 2.00000000000000000000e+00, /* 0x40000000, 0x00000000 */
erx = 8.45062911510467529297e-01, /* 0x3FEB0AC1, 0x60000000 */
/*
* Coefficients for approximation to erf on [0,0.84375]
*/
efx = 1.28379167095512586316e-01, /* 0x3FC06EBA, 0x8214DB69 */
efx8 = 1.02703333676410069053e+00, /* 0x3FF06EBA, 0x8214DB69 */
pp0 = 1.28379167095512558561e-01, /* 0x3FC06EBA, 0x8214DB68 */
pp1 = -3.25042107247001499370e-01, /* 0xBFD4CD7D, 0x691CB913 */
pp2 = -2.84817495755985104766e-02, /* 0xBF9D2A51, 0xDBD7194F */
pp3 = -5.77027029648944159157e-03, /* 0xBF77A291, 0x236668E4 */
pp4 = -2.37630166566501626084e-05, /* 0xBEF8EAD6, 0x120016AC */
qq1 = 3.97917223959155352819e-01, /* 0x3FD97779, 0xCDDADC09 */
qq2 = 6.50222499887672944485e-02, /* 0x3FB0A54C, 0x5536CEBA */
qq3 = 5.08130628187576562776e-03, /* 0x3F74D022, 0xC4D36B0F */
qq4 = 1.32494738004321644526e-04, /* 0x3F215DC9, 0x221C1A10 */
qq5 = -3.96022827877536812320e-06, /* 0xBED09C43, 0x42A26120 */
/*
* Coefficients for approximation to erf in [0.84375,1.25]
*/
pa0 = -2.36211856075265944077e-03, /* 0xBF6359B8, 0xBEF77538 */
pa1 = 4.14856118683748331666e-01, /* 0x3FDA8D00, 0xAD92B34D */
pa2 = -3.72207876035701323847e-01, /* 0xBFD7D240, 0xFBB8C3F1 */
pa3 = 3.18346619901161753674e-01, /* 0x3FD45FCA, 0x805120E4 */
pa4 = -1.10894694282396677476e-01, /* 0xBFBC6398, 0x3D3E28EC */
pa5 = 3.54783043256182359371e-02, /* 0x3FA22A36, 0x599795EB */
pa6 = -2.16637559486879084300e-03, /* 0xBF61BF38, 0x0A96073F */
qa1 = 1.06420880400844228286e-01, /* 0x3FBB3E66, 0x18EEE323 */
qa2 = 5.40397917702171048937e-01, /* 0x3FE14AF0, 0x92EB6F33 */
qa3 = 7.18286544141962662868e-02, /* 0x3FB2635C, 0xD99FE9A7 */
qa4 = 1.26171219808761642112e-01, /* 0x3FC02660, 0xE763351F */
qa5 = 1.36370839120290507362e-02, /* 0x3F8BEDC2, 0x6B51DD1C */
qa6 = 1.19844998467991074170e-02, /* 0x3F888B54, 0x5735151D */
/*
* Coefficients for approximation to erfc in [1.25,1/0.35]
*/
ra0 = -9.86494403484714822705e-03, /* 0xBF843412, 0x600D6435 */
ra1 = -6.93858572707181764372e-01, /* 0xBFE63416, 0xE4BA7360 */
ra2 = -1.05586262253232909814e+01, /* 0xC0251E04, 0x41B0E726 */
ra3 = -6.23753324503260060396e+01, /* 0xC04F300A, 0xE4CBA38D */
ra4 = -1.62396669462573470355e+02, /* 0xC0644CB1, 0x84282266 */
ra5 = -1.84605092906711035994e+02, /* 0xC067135C, 0xEBCCABB2 */
ra6 = -8.12874355063065934246e+01, /* 0xC0545265, 0x57E4D2F2 */
ra7 = -9.81432934416914548592e+00, /* 0xC023A0EF, 0xC69AC25C */
sa1 = 1.96512716674392571292e+01, /* 0x4033A6B9, 0xBD707687 */
sa2 = 1.37657754143519042600e+02, /* 0x4061350C, 0x526AE721 */
sa3 = 4.34565877475229228821e+02, /* 0x407B290D, 0xD58A1A71 */
sa4 = 6.45387271733267880336e+02, /* 0x40842B19, 0x21EC2868 */
sa5 = 4.29008140027567833386e+02, /* 0x407AD021, 0x57700314 */
sa6 = 1.08635005541779435134e+02, /* 0x405B28A3, 0xEE48AE2C */
sa7 = 6.57024977031928170135e+00, /* 0x401A47EF, 0x8E484A93 */
sa8 = -6.04244152148580987438e-02, /* 0xBFAEEFF2, 0xEE749A62 */
/*
* Coefficients for approximation to erfc in [1/.35,28]
*/
rb0 = -9.86494292470009928597e-03, /* 0xBF843412, 0x39E86F4A */
rb1 = -7.99283237680523006574e-01, /* 0xBFE993BA, 0x70C285DE */
rb2 = -1.77579549177547519889e+01, /* 0xC031C209, 0x555F995A */
rb3 = -1.60636384855821916062e+02, /* 0xC064145D, 0x43C5ED98 */
rb4 = -6.37566443368389627722e+02, /* 0xC083EC88, 0x1375F228 */
rb5 = -1.02509513161107724954e+03, /* 0xC0900461, 0x6A2E5992 */
rb6 = -4.83519191608651397019e+02, /* 0xC07E384E, 0x9BDC383F */
sb1 = 3.03380607434824582924e+01, /* 0x403E568B, 0x261D5190 */
sb2 = 3.25792512996573918826e+02, /* 0x40745CAE, 0x221B9F0A */
sb3 = 1.53672958608443695994e+03, /* 0x409802EB, 0x189D5118 */
sb4 = 3.19985821950859553908e+03, /* 0x40A8FFB7, 0x688C246A */
sb5 = 2.55305040643316442583e+03, /* 0x40A3F219, 0xCEDF3BE6 */
sb6 = 4.74528541206955367215e+02, /* 0x407DA874, 0xE79FE763 */
sb7 = -2.24409524465858183362e+01; /* 0xC03670E2, 0x42712D62 */
static double ighmm_erf (double x) {
double R,S,P,Q,s,y,z,r;
double ax = fabs(x);
#ifdef HAVE_IEEE754
if (!isfinite(x)) {
if (isnan(x)) return x; /* erf(nan)=nan */
return (x==ax) ? 1 : -1; /* erf(+-inf)=+-1 */
}
#endif
if (ax < 0.84375) { /* |x|<0.84375 */
if (ax < 3.7252903e-9) { /* |x|<2**-28 */
if (ax < tiny)
return 0.125*(8.0*x+efx8*x); /*avoid underflow */
return x + efx*x;
}
z = x*x;
r = pp0+z*(pp1+z*(pp2+z*(pp3+z*pp4)));
s = one+z*(qq1+z*(qq2+z*(qq3+z*(qq4+z*qq5))));
y = r/s;
return x + x*y;
}
if (ax < 1.25) { /* 0.84375 <= |x| < 1.25 */
s = ax-one;
P = pa0+s*(pa1+s*(pa2+s*(pa3+s*(pa4+s*(pa5+s*pa6)))));
Q = one+s*(qa1+s*(qa2+s*(qa3+s*(qa4+s*(qa5+s*qa6)))));
if (x>=0) return erx + P/Q;
else return -erx - P/Q;
}
if (ax >= 6.0) { /* inf>|x|>=6 */
if (x>=0) return one-tiny;
else return tiny-one;
}
s = one/(x*x);
if (ax < 2.857142857) { /* |x| < 1/0.35 */
R=ra0+s*(ra1+s*(ra2+s*(ra3+s*(ra4+s*(ra5+s*(ra6+s*ra7))))));
S=one+s*(sa1+s*(sa2+s*(sa3+s*(sa4+s*(sa5+s*(sa6+s*(sa7+s*sa8)))))));
}
else { /* |x| >= 1/0.35 */
R=rb0+s*(rb1+s*(rb2+s*(rb3+s*(rb4+s*(rb5+s*rb6)))));
S=one+s*(sb1+s*(sb2+s*(sb3+s*(sb4+s*(sb5+s*(sb6+s*sb7))))));
}
z = (double)(float) ax;
r = exp(-z*z-0.5625)*exp((z-ax)*(z+ax)+R/S);
if (x>=0) return one-r/ax;
else return r/ax-one;
}
static double ighmm_erfc (double x) {
double R,S,P,Q,s,y,z,r;
double ax = fabs(x);
#ifdef HAVE_IEEE754
if (!isfinite(x)) {
if (isnan(x)) return x; /* erfc(nan)=nan */
return (x==ax) ? 0 : 2; /* erfc(+-inf)=0,2 */
}
#endif
if (ax < 0.84375) { /* |x|<0.84375 */
if (ax < 13.8777878e-18) /* |x|<2**-56 */
return one-x;
z = x*x;
r = pp0+z*(pp1+z*(pp2+z*(pp3+z*pp4)));
s = one+z*(qq1+z*(qq2+z*(qq3+z*(qq4+z*qq5))));
y = r/s;
if (ax < 0.25) { /* x<1/4 */
return one-(x+x*y);
}
else {
r = x*y;
r += (x-half);
return half - r ;
}
}
if (ax < 1.25) { /* 0.84375 <= |x| < 1.25 */
s = fabs(x)-one;
P = pa0+s*(pa1+s*(pa2+s*(pa3+s*(pa4+s*(pa5+s*pa6)))));
Q = one+s*(qa1+s*(qa2+s*(qa3+s*(qa4+s*(qa5+s*qa6)))));
if (x>=0) {
z = one-erx; return z - P/Q;
}
else {
z = erx+P/Q; return one+z;
}
}
if (ax < 28.0) { /* |x|<28 */
s = one/(x*x);
if (ax < 2.857142857) { /* |x| < 1/.35 ~ 2.857143*/
R=ra0+s*(ra1+s*(ra2+s*(ra3+s*(ra4+s*(ra5+s*(ra6+s*ra7))))));
S=one+s*(sa1+s*(sa2+s*(sa3+s*(sa4+s*(sa5+s*(sa6+s*(sa7+s*sa8)))))));
}
else { /* |x| >= 1/.35 ~ 2.857143 */
if (x < -6.0) return two-tiny;/* x < -6 */
R=rb0+s*(rb1+s*(rb2+s*(rb3+s*(rb4+s*(rb5+s*rb6)))));
S=one+s*(sb1+s*(sb2+s*(sb3+s*(sb4+s*(sb5+s*(sb6+s*sb7))))));
}
z = (double)(float) ax;
r = exp(-z*z-0.5625)*exp((z-ax)*(z+ax)+R/S);
if (x>0) return r/ax;
else return two-r/ax;
}
else {
if(x>0) return tiny*tiny;
else return two-tiny;
}
}
#endif /* check for ISO C99 */