gsl_randist__exppow.c 3.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123
  1. /* randist/exppow.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2006, 2007 James Theiler, Brian Gough
  4. * Copyright (C) 2006 Giulio Bottazzi
  5. *
  6. * This program is free software; you can redistribute it and/or modify
  7. * it under the terms of the GNU General Public License as published by
  8. * the Free Software Foundation; either version 3 of the License, or (at
  9. * your option) any later version.
  10. *
  11. * This program is distributed in the hope that it will be useful, but
  12. * WITHOUT ANY WARRANTY; without even the implied warranty of
  13. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  14. * General Public License for more details.
  15. *
  16. * You should have received a copy of the GNU General Public License
  17. * along with this program; if not, write to the Free Software
  18. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
  19. */
  20. #include "gsl__config.h"
  21. #include <math.h>
  22. #include "gsl_math.h"
  23. #include "gsl_sf_gamma.h"
  24. #include "gsl_rng.h"
  25. #include "gsl_randist.h"
  26. /* The exponential power probability distribution is
  27. p(x) dx = (1/(2 a Gamma(1+1/b))) * exp(-|x/a|^b) dx
  28. for -infty < x < infty. For b = 1 it reduces to the Laplace
  29. distribution.
  30. The exponential power distribution is related to the gamma
  31. distribution by E = a * pow(G(1/b),1/b), where E is an exponential
  32. power variate and G is a gamma variate.
  33. We use this relation for b < 1. For b >=1 we use rejection methods
  34. based on the laplace and gaussian distributions which should be
  35. faster. For b>4 we revert to the gamma method.
  36. See P. R. Tadikamalla, "Random Sampling from the Exponential Power
  37. Distribution", Journal of the American Statistical Association,
  38. September 1980, Volume 75, Number 371, pages 683-686.
  39. */
  40. double
  41. gsl_ran_exppow (const gsl_rng * r, const double a, const double b)
  42. {
  43. if (b < 1 || b > 4)
  44. {
  45. double u = gsl_rng_uniform (r);
  46. double v = gsl_ran_gamma (r, 1 / b, 1.0);
  47. double z = a * pow (v, 1 / b);
  48. if (u > 0.5)
  49. {
  50. return z;
  51. }
  52. else
  53. {
  54. return -z;
  55. }
  56. }
  57. else if (b == 1)
  58. {
  59. /* Laplace distribution */
  60. return gsl_ran_laplace (r, a);
  61. }
  62. else if (b < 2)
  63. {
  64. /* Use laplace distribution for rejection method, from Tadikamalla */
  65. double x, h, u;
  66. double B = pow (1 / b, 1 / b);
  67. do
  68. {
  69. x = gsl_ran_laplace (r, B);
  70. u = gsl_rng_uniform_pos (r);
  71. h = -pow (fabs (x), b) + fabs (x) / B - 1 + (1 / b);
  72. }
  73. while (log (u) > h);
  74. return a * x;
  75. }
  76. else if (b == 2)
  77. {
  78. /* Gaussian distribution */
  79. return gsl_ran_gaussian (r, a / sqrt (2.0));
  80. }
  81. else
  82. {
  83. /* Use gaussian for rejection method, from Tadikamalla */
  84. double x, h, u;
  85. double B = pow (1 / b, 1 / b);
  86. do
  87. {
  88. x = gsl_ran_gaussian (r, B);
  89. u = gsl_rng_uniform_pos (r);
  90. h = -pow (fabs (x), b) + (x * x) / (2 * B * B) + (1 / b) - 0.5;
  91. }
  92. while (log (u) > h);
  93. return a * x;
  94. }
  95. }
  96. double
  97. gsl_ran_exppow_pdf (const double x, const double a, const double b)
  98. {
  99. double p;
  100. double lngamma = gsl_sf_lngamma (1 + 1 / b);
  101. p = (1 / (2 * a)) * exp (-pow (fabs (x / a), b) - lngamma);
  102. return p;
  103. }