gsl_randist__hyperg.c 2.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125
  1. /* randist/hyperg.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 James Theiler, Brian Gough
  4. *
  5. * This program is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU General Public License as published by
  7. * the Free Software Foundation; either version 3 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * This program is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  13. * General Public License for more details.
  14. *
  15. * You should have received a copy of the GNU General Public License
  16. * along with this program; if not, write to the Free Software
  17. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
  18. */
  19. #include "gsl__config.h"
  20. #include <math.h>
  21. #include "gsl_rng.h"
  22. #include "gsl_randist.h"
  23. #include "gsl_sf_gamma.h"
  24. /* The hypergeometric distribution has the form,
  25. prob(k) = choose(n1,t) choose(n2, t-k) / choose(n1+n2,t)
  26. where choose(a,b) = a!/(b!(a-b)!)
  27. n1 + n2 is the total population (tagged plus untagged)
  28. n1 is the tagged population
  29. t is the number of samples taken (without replacement)
  30. k is the number of tagged samples found
  31. */
  32. unsigned int
  33. gsl_ran_hypergeometric (const gsl_rng * r, unsigned int n1, unsigned int n2,
  34. unsigned int t)
  35. {
  36. const unsigned int n = n1 + n2;
  37. unsigned int i = 0;
  38. unsigned int a = n1;
  39. unsigned int b = n1 + n2;
  40. unsigned int k = 0;
  41. if (t > n)
  42. {
  43. t = n ;
  44. }
  45. if (t < n / 2)
  46. {
  47. for (i = 0 ; i < t ; i++)
  48. {
  49. double u = gsl_rng_uniform(r) ;
  50. if (b * u < a)
  51. {
  52. k++ ;
  53. if (k == n1)
  54. return k ;
  55. a-- ;
  56. }
  57. b-- ;
  58. }
  59. return k;
  60. }
  61. else
  62. {
  63. for (i = 0 ; i < n - t ; i++)
  64. {
  65. double u = gsl_rng_uniform(r) ;
  66. if (b * u < a)
  67. {
  68. k++ ;
  69. if (k == n1)
  70. return n1 - k ;
  71. a-- ;
  72. }
  73. b-- ;
  74. }
  75. return n1 - k;
  76. }
  77. }
  78. double
  79. gsl_ran_hypergeometric_pdf (const unsigned int k,
  80. const unsigned int n1,
  81. const unsigned int n2,
  82. unsigned int t)
  83. {
  84. if (t > n1 + n2)
  85. {
  86. t = n1 + n2 ;
  87. }
  88. if (k > n1 || k > t)
  89. {
  90. return 0 ;
  91. }
  92. else if (t > n2 && k + n2 < t )
  93. {
  94. return 0 ;
  95. }
  96. else
  97. {
  98. double p;
  99. double c1 = gsl_sf_lnchoose(n1,k);
  100. double c2 = gsl_sf_lnchoose(n2,t-k);
  101. double c3 = gsl_sf_lnchoose(n1+n2,t);
  102. p = exp(c1 + c2 - c3) ;
  103. return p;
  104. }
  105. }