gsl_randist.h 8.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186
  1. /* randist/gsl_randist.h
  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. #ifndef __GSL_RANDIST_H__
  20. #define __GSL_RANDIST_H__
  21. #include "gsl_rng.h"
  22. #undef __BEGIN_DECLS
  23. #undef __END_DECLS
  24. #ifdef __cplusplus
  25. # define __BEGIN_DECLS extern "C" {
  26. # define __END_DECLS }
  27. #else
  28. # define __BEGIN_DECLS /* empty */
  29. # define __END_DECLS /* empty */
  30. #endif
  31. __BEGIN_DECLS
  32. unsigned int gsl_ran_bernoulli (const gsl_rng * r, double p);
  33. double gsl_ran_bernoulli_pdf (const unsigned int k, double p);
  34. double gsl_ran_beta (const gsl_rng * r, const double a, const double b);
  35. double gsl_ran_beta_pdf (const double x, const double a, const double b);
  36. unsigned int gsl_ran_binomial (const gsl_rng * r, double p, unsigned int n);
  37. unsigned int gsl_ran_binomial_knuth (const gsl_rng * r, double p, unsigned int n);
  38. unsigned int gsl_ran_binomial_tpe (const gsl_rng * r, double p, unsigned int n);
  39. double gsl_ran_binomial_pdf (const unsigned int k, const double p, const unsigned int n);
  40. double gsl_ran_exponential (const gsl_rng * r, const double mu);
  41. double gsl_ran_exponential_pdf (const double x, const double mu);
  42. double gsl_ran_exppow (const gsl_rng * r, const double a, const double b);
  43. double gsl_ran_exppow_pdf (const double x, const double a, const double b);
  44. double gsl_ran_cauchy (const gsl_rng * r, const double a);
  45. double gsl_ran_cauchy_pdf (const double x, const double a);
  46. double gsl_ran_chisq (const gsl_rng * r, const double nu);
  47. double gsl_ran_chisq_pdf (const double x, const double nu);
  48. void gsl_ran_dirichlet (const gsl_rng * r, const size_t K, const double alpha[], double theta[]);
  49. double gsl_ran_dirichlet_pdf (const size_t K, const double alpha[], const double theta[]);
  50. double gsl_ran_dirichlet_lnpdf (const size_t K, const double alpha[], const double theta[]);
  51. double gsl_ran_erlang (const gsl_rng * r, const double a, const double n);
  52. double gsl_ran_erlang_pdf (const double x, const double a, const double n);
  53. double gsl_ran_fdist (const gsl_rng * r, const double nu1, const double nu2);
  54. double gsl_ran_fdist_pdf (const double x, const double nu1, const double nu2);
  55. double gsl_ran_flat (const gsl_rng * r, const double a, const double b);
  56. double gsl_ran_flat_pdf (double x, const double a, const double b);
  57. double gsl_ran_gamma (const gsl_rng * r, const double a, const double b);
  58. double gsl_ran_gamma_int (const gsl_rng * r, const unsigned int a);
  59. double gsl_ran_gamma_pdf (const double x, const double a, const double b);
  60. double gsl_ran_gamma_mt (const gsl_rng * r, const double a, const double b);
  61. double gsl_ran_gamma_knuth (const gsl_rng * r, const double a, const double b);
  62. double gsl_ran_gaussian (const gsl_rng * r, const double sigma);
  63. double gsl_ran_gaussian_ratio_method (const gsl_rng * r, const double sigma);
  64. double gsl_ran_gaussian_ziggurat (const gsl_rng * r, const double sigma);
  65. double gsl_ran_gaussian_pdf (const double x, const double sigma);
  66. double gsl_ran_ugaussian (const gsl_rng * r);
  67. double gsl_ran_ugaussian_ratio_method (const gsl_rng * r);
  68. double gsl_ran_ugaussian_pdf (const double x);
  69. double gsl_ran_gaussian_tail (const gsl_rng * r, const double a, const double sigma);
  70. double gsl_ran_gaussian_tail_pdf (const double x, const double a, const double sigma);
  71. double gsl_ran_ugaussian_tail (const gsl_rng * r, const double a);
  72. double gsl_ran_ugaussian_tail_pdf (const double x, const double a);
  73. void gsl_ran_bivariate_gaussian (const gsl_rng * r, double sigma_x, double sigma_y, double rho, double *x, double *y);
  74. double gsl_ran_bivariate_gaussian_pdf (const double x, const double y, const double sigma_x, const double sigma_y, const double rho);
  75. double gsl_ran_landau (const gsl_rng * r);
  76. double gsl_ran_landau_pdf (const double x);
  77. unsigned int gsl_ran_geometric (const gsl_rng * r, const double p);
  78. double gsl_ran_geometric_pdf (const unsigned int k, const double p);
  79. unsigned int gsl_ran_hypergeometric (const gsl_rng * r, unsigned int n1, unsigned int n2, unsigned int t);
  80. double gsl_ran_hypergeometric_pdf (const unsigned int k, const unsigned int n1, const unsigned int n2, unsigned int t);
  81. double gsl_ran_gumbel1 (const gsl_rng * r, const double a, const double b);
  82. double gsl_ran_gumbel1_pdf (const double x, const double a, const double b);
  83. double gsl_ran_gumbel2 (const gsl_rng * r, const double a, const double b);
  84. double gsl_ran_gumbel2_pdf (const double x, const double a, const double b);
  85. double gsl_ran_logistic (const gsl_rng * r, const double a);
  86. double gsl_ran_logistic_pdf (const double x, const double a);
  87. double gsl_ran_lognormal (const gsl_rng * r, const double zeta, const double sigma);
  88. double gsl_ran_lognormal_pdf (const double x, const double zeta, const double sigma);
  89. unsigned int gsl_ran_logarithmic (const gsl_rng * r, const double p);
  90. double gsl_ran_logarithmic_pdf (const unsigned int k, const double p);
  91. void gsl_ran_multinomial (const gsl_rng * r, const size_t K,
  92. const unsigned int N, const double p[],
  93. unsigned int n[] );
  94. double gsl_ran_multinomial_pdf (const size_t K,
  95. const double p[], const unsigned int n[] );
  96. double gsl_ran_multinomial_lnpdf (const size_t K,
  97. const double p[], const unsigned int n[] );
  98. unsigned int gsl_ran_negative_binomial (const gsl_rng * r, double p, double n);
  99. double gsl_ran_negative_binomial_pdf (const unsigned int k, const double p, double n);
  100. unsigned int gsl_ran_pascal (const gsl_rng * r, double p, unsigned int n);
  101. double gsl_ran_pascal_pdf (const unsigned int k, const double p, unsigned int n);
  102. double gsl_ran_pareto (const gsl_rng * r, double a, const double b);
  103. double gsl_ran_pareto_pdf (const double x, const double a, const double b);
  104. unsigned int gsl_ran_poisson (const gsl_rng * r, double mu);
  105. void gsl_ran_poisson_array (const gsl_rng * r, size_t n, unsigned int array[],
  106. double mu);
  107. double gsl_ran_poisson_pdf (const unsigned int k, const double mu);
  108. double gsl_ran_rayleigh (const gsl_rng * r, const double sigma);
  109. double gsl_ran_rayleigh_pdf (const double x, const double sigma);
  110. double gsl_ran_rayleigh_tail (const gsl_rng * r, const double a, const double sigma);
  111. double gsl_ran_rayleigh_tail_pdf (const double x, const double a, const double sigma);
  112. double gsl_ran_tdist (const gsl_rng * r, const double nu);
  113. double gsl_ran_tdist_pdf (const double x, const double nu);
  114. double gsl_ran_laplace (const gsl_rng * r, const double a);
  115. double gsl_ran_laplace_pdf (const double x, const double a);
  116. double gsl_ran_levy (const gsl_rng * r, const double c, const double alpha);
  117. double gsl_ran_levy_skew (const gsl_rng * r, const double c, const double alpha, const double beta);
  118. double gsl_ran_weibull (const gsl_rng * r, const double a, const double b);
  119. double gsl_ran_weibull_pdf (const double x, const double a, const double b);
  120. void gsl_ran_dir_2d (const gsl_rng * r, double * x, double * y);
  121. void gsl_ran_dir_2d_trig_method (const gsl_rng * r, double * x, double * y);
  122. void gsl_ran_dir_3d (const gsl_rng * r, double * x, double * y, double * z);
  123. void gsl_ran_dir_nd (const gsl_rng * r, size_t n, double * x);
  124. void gsl_ran_shuffle (const gsl_rng * r, void * base, size_t nmembm, size_t size);
  125. int gsl_ran_choose (const gsl_rng * r, void * dest, size_t k, void * src, size_t n, size_t size) ;
  126. void gsl_ran_sample (const gsl_rng * r, void * dest, size_t k, void * src, size_t n, size_t size) ;
  127. typedef struct { /* struct for Walker algorithm */
  128. size_t K;
  129. size_t *A;
  130. double *F;
  131. } gsl_ran_discrete_t;
  132. gsl_ran_discrete_t * gsl_ran_discrete_preproc (size_t K, const double *P);
  133. void gsl_ran_discrete_free(gsl_ran_discrete_t *g);
  134. size_t gsl_ran_discrete (const gsl_rng *r, const gsl_ran_discrete_t *g);
  135. double gsl_ran_discrete_pdf (size_t k, const gsl_ran_discrete_t *g);
  136. __END_DECLS
  137. #endif /* __GSL_RANDIST_H__ */