Blame randist/gausstail.c

Packit 67cb25
/* randist/gausstail.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 James Theiler, 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
#include <config.h>
Packit 67cb25
#include <math.h>
Packit 67cb25
#include <gsl/gsl_math.h>
Packit 67cb25
#include <gsl/gsl_rng.h>
Packit 67cb25
#include <gsl/gsl_randist.h>
Packit 67cb25
#include <gsl/gsl_sf_erf.h>
Packit 67cb25
Packit 67cb25
double
Packit 67cb25
gsl_ran_gaussian_tail (const gsl_rng * r, const double a, const double sigma)
Packit 67cb25
{
Packit 67cb25
  /* Returns a gaussian random variable larger than a
Packit 67cb25
   * This implementation does one-sided upper-tailed deviates.
Packit 67cb25
   */
Packit 67cb25
Packit 67cb25
  double s = a / sigma;
Packit 67cb25
Packit 67cb25
  if (s < 1)
Packit 67cb25
    {
Packit 67cb25
      /* For small s, use a direct rejection method. The limit s < 1
Packit 67cb25
         can be adjusted to optimise the overall efficiency */
Packit 67cb25
Packit 67cb25
      double x;
Packit 67cb25
Packit 67cb25
      do
Packit 67cb25
        {
Packit 67cb25
          x = gsl_ran_gaussian (r, 1.0);
Packit 67cb25
        }
Packit 67cb25
      while (x < s);
Packit 67cb25
      return x * sigma;
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      /* Use the "supertail" deviates from the last two steps
Packit 67cb25
       * of Marsaglia's rectangle-wedge-tail method, as described
Packit 67cb25
       * in Knuth, v2, 3rd ed, pp 123-128.  (See also exercise 11, p139,
Packit 67cb25
       * and the solution, p586.)
Packit 67cb25
       */
Packit 67cb25
Packit 67cb25
      double u, v, x;
Packit 67cb25
Packit 67cb25
      do
Packit 67cb25
        {
Packit 67cb25
          u = gsl_rng_uniform (r);
Packit 67cb25
          do
Packit 67cb25
            {
Packit 67cb25
              v = gsl_rng_uniform (r);
Packit 67cb25
            }
Packit 67cb25
          while (v == 0.0);
Packit 67cb25
          x = sqrt (s * s - 2 * log (v));
Packit 67cb25
        }
Packit 67cb25
      while (x * u > s);
Packit 67cb25
      return x * sigma;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double
Packit 67cb25
gsl_ran_gaussian_tail_pdf (const double x, const double a, const double sigma)
Packit 67cb25
{
Packit 67cb25
  if (x < a)
Packit 67cb25
    {
Packit 67cb25
      return 0;
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      double N, p;
Packit 67cb25
      double u = x / sigma ;
Packit 67cb25
Packit 67cb25
      double f = gsl_sf_erfc (a / (sqrt (2.0) * sigma));
Packit 67cb25
Packit 67cb25
      N = 0.5 * f;
Packit 67cb25
Packit 67cb25
      p = (1 / (N * sqrt (2 * M_PI) * sigma)) * exp (-u * u / 2);
Packit 67cb25
Packit 67cb25
      return p;
Packit 67cb25
    }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double
Packit 67cb25
gsl_ran_ugaussian_tail (const gsl_rng * r, const double a)
Packit 67cb25
{
Packit 67cb25
  return gsl_ran_gaussian_tail (r, a, 1.0) ;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
double
Packit 67cb25
gsl_ran_ugaussian_tail_pdf (const double x, const double a)
Packit 67cb25
{
Packit 67cb25
  return gsl_ran_gaussian_tail_pdf (x, a, 1.0) ;
Packit 67cb25
}