gsl_randist__dirichlet.c 4.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164
  1. /* randist/dirichlet.c
  2. *
  3. * Copyright (C) 2007 Brian Gough
  4. * Copyright (C) 2002 Gavin E. Crooks <gec@compbio.berkeley.edu>
  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_rng.h"
  24. #include "gsl_randist.h"
  25. #include "gsl_sf_gamma.h"
  26. /* The Dirichlet probability distribution of order K-1 is
  27. p(\theta_1,...,\theta_K) d\theta_1 ... d\theta_K =
  28. (1/Z) \prod_i=1,K \theta_i^{alpha_i - 1} \delta(1 -\sum_i=1,K \theta_i)
  29. The normalization factor Z can be expressed in terms of gamma functions:
  30. Z = {\prod_i=1,K \Gamma(\alpha_i)} / {\Gamma( \sum_i=1,K \alpha_i)}
  31. The K constants, \alpha_1,...,\alpha_K, must be positive. The K parameters,
  32. \theta_1,...,\theta_K are nonnegative and sum to 1.
  33. The random variates are generated by sampling K values from gamma
  34. distributions with parameters a=\alpha_i, b=1, and renormalizing.
  35. See A.M. Law, W.D. Kelton, Simulation Modeling and Analysis (1991).
  36. Gavin E. Crooks <gec@compbio.berkeley.edu> (2002)
  37. */
  38. static void ran_dirichlet_small (const gsl_rng * r, const size_t K, const double alpha[], double theta[]);
  39. void
  40. gsl_ran_dirichlet (const gsl_rng * r, const size_t K,
  41. const double alpha[], double theta[])
  42. {
  43. size_t i;
  44. double norm = 0.0;
  45. for (i = 0; i < K; i++)
  46. {
  47. theta[i] = gsl_ran_gamma (r, alpha[i], 1.0);
  48. }
  49. for (i = 0; i < K; i++)
  50. {
  51. norm += theta[i];
  52. }
  53. if (norm < GSL_SQRT_DBL_MIN) /* Handle underflow */
  54. {
  55. ran_dirichlet_small (r, K, alpha, theta);
  56. return;
  57. }
  58. for (i = 0; i < K; i++)
  59. {
  60. theta[i] /= norm;
  61. }
  62. }
  63. /* When the values of alpha[] are small, scale the variates to avoid
  64. underflow so that the result is not 0/0. Note that the Dirichlet
  65. distribution is defined by a ratio of gamma functions so we can
  66. take out an arbitrary factor to keep the values in the range of
  67. double precision. */
  68. static void
  69. ran_dirichlet_small (const gsl_rng * r, const size_t K,
  70. const double alpha[], double theta[])
  71. {
  72. size_t i;
  73. double norm = 0.0, umax = 0;
  74. for (i = 0; i < K; i++)
  75. {
  76. double u = log(gsl_rng_uniform_pos (r)) / alpha[i];
  77. theta[i] = u;
  78. if (u > umax || i == 0) {
  79. umax = u;
  80. }
  81. }
  82. for (i = 0; i < K; i++)
  83. {
  84. theta[i] = exp(theta[i] - umax);
  85. }
  86. for (i = 0; i < K; i++)
  87. {
  88. theta[i] = theta[i] * gsl_ran_gamma (r, alpha[i] + 1.0, 1.0);
  89. }
  90. for (i = 0; i < K; i++)
  91. {
  92. norm += theta[i];
  93. }
  94. for (i = 0; i < K; i++)
  95. {
  96. theta[i] /= norm;
  97. }
  98. }
  99. double
  100. gsl_ran_dirichlet_pdf (const size_t K,
  101. const double alpha[], const double theta[])
  102. {
  103. return exp (gsl_ran_dirichlet_lnpdf (K, alpha, theta));
  104. }
  105. double
  106. gsl_ran_dirichlet_lnpdf (const size_t K,
  107. const double alpha[], const double theta[])
  108. {
  109. /*We calculate the log of the pdf to minimize the possibility of overflow */
  110. size_t i;
  111. double log_p = 0.0;
  112. double sum_alpha = 0.0;
  113. for (i = 0; i < K; i++)
  114. {
  115. log_p += (alpha[i] - 1.0) * log (theta[i]);
  116. }
  117. for (i = 0; i < K; i++)
  118. {
  119. sum_alpha += alpha[i];
  120. }
  121. log_p += gsl_sf_lngamma (sum_alpha);
  122. for (i = 0; i < K; i++)
  123. {
  124. log_p -= gsl_sf_lngamma (alpha[i]);
  125. }
  126. return log_p;
  127. }