Ergebnis für URL: http://www.gnu.org/software/gsl/doc/html/randist.html#references-and-further-reading #[1]Index [2]Search [3]Statistics [4]Quasi-Random Sequences
[5]GSL
2.8
____________________
* [6]Introduction
* [7]Using the Library
* [8]Error Handling
* [9]Mathematical Functions
* [10]Complex Numbers
* [11]Polynomials
* [12]Special Functions
* [13]Vectors and Matrices
* [14]Permutations
* [15]Combinations
* [16]Multisets
* [17]Sorting
* [18]BLAS Support
* [19]Linear Algebra
* [20]Eigensystems
* [21]Fast Fourier Transforms (FFTs)
* [22]Numerical Integration
* [23]Random Number Generation
* [24]Quasi-Random Sequences
* [25]Random Number Distributions
+ [26]Introduction
+ [27]The Gaussian Distribution
+ [28]The Gaussian Tail Distribution
+ [29]The Bivariate Gaussian Distribution
+ [30]The Multivariate Gaussian Distribution
+ [31]The Exponential Distribution
+ [32]The Laplace Distribution
+ [33]The Exponential Power Distribution
+ [34]The Cauchy Distribution
+ [35]The Rayleigh Distribution
+ [36]The Rayleigh Tail Distribution
+ [37]The Landau Distribution
+ [38]The Levy alpha-Stable Distributions
+ [39]The Levy skew alpha-Stable Distribution
+ [40]The Gamma Distribution
+ [41]The Flat (Uniform) Distribution
+ [42]The Lognormal Distribution
+ [43]The Chi-squared Distribution
+ [44]The F-distribution
+ [45]The t-distribution
+ [46]The Beta Distribution
+ [47]The Logistic Distribution
+ [48]The Pareto Distribution
+ [49]Spherical Vector Distributions
+ [50]The Weibull Distribution
+ [51]The Type-1 Gumbel Distribution
+ [52]The Type-2 Gumbel Distribution
+ [53]The Dirichlet Distribution
+ [54]General Discrete Distributions
+ [55]The Poisson Distribution
+ [56]The Bernoulli Distribution
+ [57]The Binomial Distribution
+ [58]The Multinomial Distribution
+ [59]The Negative Binomial Distribution
+ [60]The Pascal Distribution
+ [61]The Geometric Distribution
+ [62]The Hypergeometric Distribution
+ [63]The Logarithmic Distribution
+ [64]The Wishart Distribution
+ [65]Shuffling and Sampling
+ [66]Examples
+ [67]References and Further Reading
* [68]Statistics
* [69]Running Statistics
* [70]Moving Window Statistics
* [71]Digital Filtering
* [72]Histograms
* [73]N-tuples
* [74]Monte Carlo Integration
* [75]Simulated Annealing
* [76]Ordinary Differential Equations
* [77]Interpolation
* [78]Numerical Differentiation
* [79]Chebyshev Approximations
* [80]Series Acceleration
* [81]Wavelet Transforms
* [82]Discrete Hankel Transforms
* [83]One Dimensional Root-Finding
* [84]One Dimensional Minimization
* [85]Multidimensional Root-Finding
* [86]Multidimensional Minimization
* [87]Linear Least-Squares Fitting
* [88]Nonlinear Least-Squares Fitting
* [89]Basis Splines
* [90]Sparse Matrices
* [91]Sparse BLAS Support
* [92]Sparse Linear Algebra
* [93]Physical Constants
* [94]IEEE floating-point arithmetic
* [95]Debugging Numerical Programs
* [96]Contributors to GSL
* [97]Autoconf Macros
* [98]GSL CBLAS Library
* [99]GNU General Public License
* [100]GNU Free Documentation License
[101]GSL
* »
* Random Number Distributions
* [102]View page source
[103]Next [104]Previous
____________________________________________________________________________
Random Number Distributions[105]¶
This chapter describes functions for generating random variates and computing
their probability distributions. Samples from the distributions described in this
chapter can be obtained using any of the random number generators in the library
as an underlying source of randomness.
In the simplest cases a non-uniform distribution can be obtained analytically
from the uniform distribution of a random number generator by applying an
appropriate transformation. This method uses one call to the random number
generator. More complicated distributions are created by the acceptance-rejection
method, which compares the desired distribution against a distribution which is
similar and known analytically. This usually requires several samples from the
generator.
The library also provides cumulative distribution functions and inverse
cumulative distribution functions, sometimes referred to as quantile functions.
The cumulative distribution functions and their inverses are computed separately
for the upper and lower tails of the distribution, allowing full accuracy to be
retained for small results.
The functions for random variates and probability density functions described in
this section are declared in gsl_randist.h. The corresponding cumulative
distribution functions are declared in gsl_cdf.h.
Note that the discrete random variate functions always return a value of type
unsigned int, and on most platforms this has a maximum value of
2^{32}-1 \approx 4.29 \times 10^9
They should only be called with a safe range of parameters (where there is a
negligible probability of a variate exceeding this limit) to prevent incorrect
results due to overflow.
Introduction[106]¶
Continuous random number distributions are defined by a probability density
function, p(x) , such that the probability of x occurring in the infinitesimal
range x to x + dx is p(x) dx .
The cumulative distribution function for the lower tail P(x) is defined by the
integral,
P(x) = \int_{-\infty}^{x} dx' p(x')
and gives the probability of a variate taking a value less than x .
The cumulative distribution function for the upper tail Q(x) is defined by the
integral,
Q(x) = \int_{x}^{+\infty} dx' p(x')
and gives the probability of a variate taking a value greater than x .
The upper and lower cumulative distribution functions are related by P(x) + Q(x)
= 1 and satisfy 0 \le P(x) \le 1 , 0 \le Q(x) \le 1 .
The inverse cumulative distributions, x = P^{-1}(P) and x = Q^{-1}(Q) give the
values of x which correspond to a specific value of P or Q . They can be used to
find confidence limits from probability values.
For discrete distributions the probability of sampling the integer value k is
given by p(k) , where \sum_k p(k) = 1 . The cumulative distribution for the lower
tail P(k) of a discrete distribution is defined as,
P(k) = \sum_{i \le k} p(i)
where the sum is over the allowed range of the distribution less than or equal to
k .
The cumulative distribution for the upper tail of a discrete distribution Q(k) is
defined as
Q(k) = \sum_{i > k} p(i)
giving the sum of probabilities for all values greater than k . These two
definitions satisfy the identity P(k)+Q(k)=1 .
If the range of the distribution is 1 to n inclusive then P(n) = 1 , Q(n) = 0
while P(1) = p(1) , Q(1) = 1 - p(1) .
The Gaussian Distribution[107]¶
double gsl_ran_gaussian(const [108]gsl_rng *r, double sigma)[109]¶
This function returns a Gaussian random variate, with mean zero and
standard deviation [110]sigma. The probability distribution for Gaussian
random variates is,
p(x) dx = {1 \over \sqrt{2 \pi \sigma^2}} \exp (-x^2 / 2\sigma^2) dx
for x in the range -\infty to +\infty . Use the transformation z = \mu + x
on the numbers returned by [111]gsl_ran_gaussian() to obtain a Gaussian
distribution with mean \mu . This function uses the Box-Muller algorithm
which requires two calls to the random number generator [112]r.
double gsl_ran_gaussian_pdf(double x, double sigma)[113]¶
This function computes the probability density p(x) at [114]x for a
Gaussian distribution with standard deviation [115]sigma, using the
formula given above.
_images/rand-gaussian.png
double gsl_ran_gaussian_ziggurat(const [116]gsl_rng *r, double sigma)[117]¶
double gsl_ran_gaussian_ratio_method(const [118]gsl_rng *r, double sigma)[119]¶
This function computes a Gaussian random variate using the alternative
Marsaglia-Tsang ziggurat and Kinderman-Monahan-Leva ratio methods. The
Ziggurat algorithm is the fastest available algorithm in most cases.
double gsl_ran_ugaussian(const [120]gsl_rng *r)[121]¶
double gsl_ran_ugaussian_pdf(double x)[122]¶
double gsl_ran_ugaussian_ratio_method(const [123]gsl_rng *r)[124]¶
These functions compute results for the unit Gaussian distribution. They
are equivalent to the functions above with a standard deviation of one,
sigma = 1.
double gsl_cdf_gaussian_P(double x, double sigma)[125]¶
double gsl_cdf_gaussian_Q(double x, double sigma)[126]¶
double gsl_cdf_gaussian_Pinv(double P, double sigma)[127]¶
double gsl_cdf_gaussian_Qinv(double Q, double sigma)[128]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the Gaussian distribution with standard deviation
[129]sigma.
double gsl_cdf_ugaussian_P(double x)[130]¶
double gsl_cdf_ugaussian_Q(double x)[131]¶
double gsl_cdf_ugaussian_Pinv(double P)[132]¶
double gsl_cdf_ugaussian_Qinv(double Q)[133]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the unit Gaussian distribution.
The Gaussian Tail Distribution[134]¶
double gsl_ran_gaussian_tail(const [135]gsl_rng *r, double a, double sigma)[136]¶
This function provides random variates from the upper tail of a Gaussian
distribution with standard deviation [137]sigma. The values returned are
larger than the lower limit [138]a, which must be positive. The method is
based on Marsaglia's famous rectangle-wedge-tail algorithm (Ann. Math.
Stat. 32, 894-899 (1961)), with this aspect explained in Knuth, v2, 3rd
ed, p139,586 (exercise 11).
The probability distribution for Gaussian tail random variates is,
p(x) dx = {1 \over N(a;\sigma) \sqrt{2 \pi \sigma^2}} \exp (- x^2 /
2\sigma^2) dx
for x > a where N(a;\sigma) is the normalization constant,
N(a;\sigma) = {1 \over 2} \hbox{erfc}\left({a \over \sqrt{2
\sigma^2}}\right).
double gsl_ran_gaussian_tail_pdf(double x, double a, double sigma)[139]¶
This function computes the probability density p(x) at [140]x for a
Gaussian tail distribution with standard deviation [141]sigma and lower
limit [142]a, using the formula given above.
_images/rand-gaussian-tail.png
double gsl_ran_ugaussian_tail(const [143]gsl_rng *r, double a)[144]¶
double gsl_ran_ugaussian_tail_pdf(double x, double a)[145]¶
These functions compute results for the tail of a unit Gaussian
distribution. They are equivalent to the functions above with a standard
deviation of one, sigma = 1.
The Bivariate Gaussian Distribution[146]¶
void gsl_ran_bivariate_gaussian(const [147]gsl_rng *r, double sigma_x, double
sigma_y, double rho, double *x, double *y)[148]¶
This function generates a pair of correlated Gaussian variates, with mean
zero, correlation coefficient [149]rho and standard deviations
[150]sigma_x and [151]sigma_y in the x and y directions. The probability
distribution for bivariate Gaussian random variates is,
p(x,y) dx dy = {1 \over 2 \pi \sigma_x \sigma_y \sqrt{1-\rho^2}} \exp
\left(-{(x^2/\sigma_x^2 + y^2/\sigma_y^2 - 2 \rho x y/(\sigma_x\sigma_y))
\over 2(1-\rho^2)}\right) dx dy
for x,y in the range -\infty to +\infty . The correlation coefficient
[152]rho should lie between 1 and -1 .
double gsl_ran_bivariate_gaussian_pdf(double x, double y, double sigma_x, double
sigma_y, double rho)[153]¶
This function computes the probability density p(x,y) at ([154]x, [155]y)
for a bivariate Gaussian distribution with standard deviations
[156]sigma_x, [157]sigma_y and correlation coefficient [158]rho, using the
formula given above.
_images/rand-bivariate-gaussian.png
The Multivariate Gaussian Distribution[159]¶
int gsl_ran_multivariate_gaussian(const [160]gsl_rng *r, const [161]gsl_vector
*mu, const [162]gsl_matrix *L, [163]gsl_vector *result)[164]¶
This function generates a random vector satisfying the k -dimensional
multivariate Gaussian distribution with mean \mu and variance-covariance
matrix \Sigma . On input, the k -vector \mu is given in [165]mu, and the
Cholesky factor of the k -by- k matrix \Sigma = L L^T is given in the
lower triangle of [166]L, as output from
[167]gsl_linalg_cholesky_decomp(). The random vector is stored in
[168]result on output. The probability distribution for multivariate
Gaussian random variates is
p(x_1,\dots,x_k) dx_1 \dots dx_k = {1 \over \sqrt{(2 \pi)^k |\Sigma|}}
\exp \left(-{1 \over 2} (x - \mu)^T \Sigma^{-1} (x - \mu)\right) dx_1
\dots dx_k
int gsl_ran_multivariate_gaussian_pdf(const [169]gsl_vector *x, const
[170]gsl_vector *mu, const [171]gsl_matrix *L, double *result,
[172]gsl_vector *work)[173]¶
int gsl_ran_multivariate_gaussian_log_pdf(const [174]gsl_vector *x, const
[175]gsl_vector *mu, const [176]gsl_matrix *L, double *result,
[177]gsl_vector *work)[178]¶
These functions compute p(x) or \log{p(x)} at the point [179]x, using mean
vector [180]mu and variance-covariance matrix specified by its Cholesky
factor [181]L using the formula above. Additional workspace of length k is
required in [182]work.
int gsl_ran_multivariate_gaussian_mean(const [183]gsl_matrix *X, [184]gsl_vector
*mu_hat)[185]¶
Given a set of n samples X_j from a k -dimensional multivariate Gaussian
distribution, this function computes the maximum likelihood estimate of
the mean of the distribution, given by
\Hat{\mu} = {1 \over n} \sum_{j=1}^n X_j
The samples X_1,X_2,\dots,X_n are given in the n -by- k matrix [186]X, and
the maximum likelihood estimate of the mean is stored in [187]mu_hat on
output.
int gsl_ran_multivariate_gaussian_vcov(const [188]gsl_matrix *X, [189]gsl_matrix
*sigma_hat)[190]¶
Given a set of n samples X_j from a k -dimensional multivariate Gaussian
distribution, this function computes the maximum likelihood estimate of
the variance-covariance matrix of the distribution, given by
\Hat{\Sigma} = {1 \over n} \sum_{j=1}^n \left( X_j - \Hat{\mu} \right)
\left( X_j - \Hat{\mu} \right)^T
The samples X_1,X_2,\dots,X_n are given in the n -by- k matrix [191]X and
the maximum likelihood estimate of the variance-covariance matrix is
stored in [192]sigma_hat on output.
The Exponential Distribution[193]¶
double gsl_ran_exponential(const [194]gsl_rng *r, double mu)[195]¶
This function returns a random variate from the exponential distribution
with mean [196]mu. The distribution is,
p(x) dx = {1 \over \mu} \exp(-x/\mu) dx
for x \ge 0 .
double gsl_ran_exponential_pdf(double x, double mu)[197]¶
This function computes the probability density p(x) at [198]x for an
exponential distribution with mean [199]mu, using the formula given above.
_images/rand-exponential.png
double gsl_cdf_exponential_P(double x, double mu)[200]¶
double gsl_cdf_exponential_Q(double x, double mu)[201]¶
double gsl_cdf_exponential_Pinv(double P, double mu)[202]¶
double gsl_cdf_exponential_Qinv(double Q, double mu)[203]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the exponential distribution with mean [204]mu.
The Laplace Distribution[205]¶
double gsl_ran_laplace(const [206]gsl_rng *r, double a)[207]¶
This function returns a random variate from the Laplace distribution with
width [208]a. The distribution is,
p(x) dx = {1 \over 2 a} \exp(-|x/a|) dx
for -\infty < x < \infty .
double gsl_ran_laplace_pdf(double x, double a)[209]¶
This function computes the probability density p(x) at [210]x for a
Laplace distribution with width [211]a, using the formula given above.
_images/rand-laplace.png
double gsl_cdf_laplace_P(double x, double a)[212]¶
double gsl_cdf_laplace_Q(double x, double a)[213]¶
double gsl_cdf_laplace_Pinv(double P, double a)[214]¶
double gsl_cdf_laplace_Qinv(double Q, double a)[215]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the Laplace distribution with width [216]a.
The Exponential Power Distribution[217]¶
double gsl_ran_exppow(const [218]gsl_rng *r, double a, double b)[219]¶
This function returns a random variate from the exponential power
distribution with scale parameter [220]a and exponent [221]b. The
distribution is,
p(x) dx = {1 \over 2 a \Gamma(1+1/b)} \exp(-|x/a|^b) dx
for x \ge 0 . For b = 1 this reduces to the Laplace distribution. For b =
2 it has the same form as a Gaussian distribution, but with a = \sqrt{2}
\sigma .
double gsl_ran_exppow_pdf(double x, double a, double b)[222]¶
This function computes the probability density p(x) at [223]x for an
exponential power distribution with scale parameter [224]a and exponent
[225]b, using the formula given above.
_images/rand-exppow.png
double gsl_cdf_exppow_P(double x, double a, double b)[226]¶
double gsl_cdf_exppow_Q(double x, double a, double b)[227]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
for the exponential power distribution with parameters [228]a and [229]b.
The Cauchy Distribution[230]¶
double gsl_ran_cauchy(const [231]gsl_rng *r, double a)[232]¶
This function returns a random variate from the Cauchy distribution with
scale parameter [233]a. The probability distribution for Cauchy random
variates is,
p(x) dx = {1 \over a\pi (1 + (x/a)^2) } dx
for x in the range -\infty to +\infty . The Cauchy distribution is also
known as the Lorentz distribution.
double gsl_ran_cauchy_pdf(double x, double a)[234]¶
This function computes the probability density p(x) at [235]x for a Cauchy
distribution with scale parameter [236]a, using the formula given above.
_images/rand-cauchy.png
double gsl_cdf_cauchy_P(double x, double a)[237]¶
double gsl_cdf_cauchy_Q(double x, double a)[238]¶
double gsl_cdf_cauchy_Pinv(double P, double a)[239]¶
double gsl_cdf_cauchy_Qinv(double Q, double a)[240]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the Cauchy distribution with scale parameter
[241]a.
The Rayleigh Distribution[242]¶
double gsl_ran_rayleigh(const [243]gsl_rng *r, double sigma)[244]¶
This function returns a random variate from the Rayleigh distribution with
scale parameter [245]sigma. The distribution is,
p(x) dx = {x \over \sigma^2} \exp(- x^2/(2 \sigma^2)) dx
for x > 0 .
double gsl_ran_rayleigh_pdf(double x, double sigma)[246]¶
This function computes the probability density p(x) at [247]x for a
Rayleigh distribution with scale parameter [248]sigma, using the formula
given above.
_images/rand-rayleigh.png
double gsl_cdf_rayleigh_P(double x, double sigma)[249]¶
double gsl_cdf_rayleigh_Q(double x, double sigma)[250]¶
double gsl_cdf_rayleigh_Pinv(double P, double sigma)[251]¶
double gsl_cdf_rayleigh_Qinv(double Q, double sigma)[252]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the Rayleigh distribution with scale parameter
[253]sigma.
The Rayleigh Tail Distribution[254]¶
double gsl_ran_rayleigh_tail(const [255]gsl_rng *r, double a, double sigma)[256]¶
This function returns a random variate from the tail of the Rayleigh
distribution with scale parameter [257]sigma and a lower limit of [258]a.
The distribution is,
p(x) dx = {x \over \sigma^2} \exp ((a^2 - x^2) /(2 \sigma^2)) dx
for x > a .
double gsl_ran_rayleigh_tail_pdf(double x, double a, double sigma)[259]¶
This function computes the probability density p(x) at [260]x for a
Rayleigh tail distribution with scale parameter [261]sigma and lower limit
[262]a, using the formula given above.
_images/rand-rayleigh-tail.png
The Landau Distribution[263]¶
double gsl_ran_landau(const [264]gsl_rng *r)[265]¶
This function returns a random variate from the Landau distribution. The
probability distribution for Landau random variates is defined
analytically by the complex integral,
p(x) = {1 \over {2 \pi i}} \int_{c-i\infty}^{c+i\infty} ds\, \exp(s
\log(s) + x s)
For numerical purposes it is more convenient to use the following
equivalent form of the integral,
p(x) = (1/\pi) \int_0^\infty dt \exp(-t \log(t) - x t) \sin(\pi t).
double gsl_ran_landau_pdf(double x)[266]¶
This function computes the probability density p(x) at [267]x for the
Landau distribution using an approximation to the formula given above.
_images/rand-landau.png
The Levy alpha-Stable Distributions[268]¶
double gsl_ran_levy(const [269]gsl_rng *r, double c, double alpha)[270]¶
This function returns a random variate from the Levy symmetric stable
distribution with scale [271]c and exponent [272]alpha. The symmetric
stable probability distribution is defined by a Fourier transform,
p(x) = {1 \over 2 \pi} \int_{-\infty}^{+\infty} dt \exp(-it x - |c
t|^\alpha)
There is no explicit solution for the form of p(x) and the library does
not define a corresponding pdf function. For \alpha = 1 the distribution
reduces to the Cauchy distribution. For \alpha = 2 it is a Gaussian
distribution with \sigma = \sqrt{2} c . For \alpha < 1 the tails of the
distribution become extremely wide.
The algorithm only works for 0 < \alpha \le 2 .
_images/rand-levy.png
The Levy skew alpha-Stable Distribution[273]¶
double gsl_ran_levy_skew(const [274]gsl_rng *r, double c, double alpha, double
beta)[275]¶
This function returns a random variate from the Levy skew stable
distribution with scale [276]c, exponent [277]alpha and skewness parameter
[278]beta. The skewness parameter must lie in the range [-1,1] . The Levy
skew stable probability distribution is defined by a Fourier transform,
p(x) = {1 \over 2 \pi} \int_{-\infty}^{+\infty} dt \exp(-it x - |c
t|^\alpha (1-i \beta \sgn(t) \tan(\pi\alpha/2)))
When \alpha = 1 the term \tan(\pi \alpha/2) is replaced by -(2/\pi)\log|t|
. There is no explicit solution for the form of p(x) and the library does
not define a corresponding pdf function. For \alpha = 2 the distribution
reduces to a Gaussian distribution with \sigma = \sqrt{2} c and the
skewness parameter has no effect. For \alpha < 1 the tails of the
distribution become extremely wide. The symmetric distribution corresponds
to \beta = 0 .
The algorithm only works for 0 < \alpha \le 2 .
The Levy alpha-stable distributions have the property that if N alpha-stable
variates are drawn from the distribution p(c, \alpha, \beta) then the sum Y = X_1
+ X_2 + \dots + X_N will also be distributed as an alpha-stable variate,
p(N^{1/\alpha} c, \alpha, \beta) .
_images/rand-levyskew.png
The Gamma Distribution[279]¶
double gsl_ran_gamma(const [280]gsl_rng *r, double a, double b)[281]¶
This function returns a random variate from the gamma distribution. The
distribution function is,
p(x) dx = {1 \over \Gamma(a) b^a} x^{a-1} e^{-x/b} dx
for x > 0 .
The gamma distribution with an integer parameter [282]a is known as the
Erlang distribution.
The variates are computed using the Marsaglia-Tsang fast gamma method.
This function for this method was previously called gsl_ran_gamma_mt() and
can still be accessed using this name.
double gsl_ran_gamma_knuth(const [283]gsl_rng *r, double a, double b)[284]¶
This function returns a gamma variate using the algorithms from Knuth (vol
2).
double gsl_ran_gamma_pdf(double x, double a, double b)[285]¶
This function computes the probability density p(x) at [286]x for a gamma
distribution with parameters [287]a and [288]b, using the formula given
above.
_images/rand-gamma.png
double gsl_cdf_gamma_P(double x, double a, double b)[289]¶
double gsl_cdf_gamma_Q(double x, double a, double b)[290]¶
double gsl_cdf_gamma_Pinv(double P, double a, double b)[291]¶
double gsl_cdf_gamma_Qinv(double Q, double a, double b)[292]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the gamma distribution with parameters [293]a and
[294]b.
The Flat (Uniform) Distribution[295]¶
double gsl_ran_flat(const [296]gsl_rng *r, double a, double b)[297]¶
This function returns a random variate from the flat (uniform)
distribution from [298]a to [299]b. The distribution is,
p(x) dx = {1 \over (b-a)} dx
if a \le x < b and 0 otherwise.
double gsl_ran_flat_pdf(double x, double a, double b)[300]¶
This function computes the probability density p(x) at [301]x for a
uniform distribution from [302]a to [303]b, using the formula given above.
_images/rand-flat.png
double gsl_cdf_flat_P(double x, double a, double b)[304]¶
double gsl_cdf_flat_Q(double x, double a, double b)[305]¶
double gsl_cdf_flat_Pinv(double P, double a, double b)[306]¶
double gsl_cdf_flat_Qinv(double Q, double a, double b)[307]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for a uniform distribution from [308]a to [309]b.
The Lognormal Distribution[310]¶
double gsl_ran_lognormal(const [311]gsl_rng *r, double zeta, double sigma)[312]¶
This function returns a random variate from the lognormal distribution.
The distribution function is,
p(x) dx = {1 \over x \sqrt{2 \pi \sigma^2}} \exp(-(\ln(x) - \zeta)^2/2
\sigma^2) dx
for x > 0 .
double gsl_ran_lognormal_pdf(double x, double zeta, double sigma)[313]¶
This function computes the probability density p(x) at [314]x for a
lognormal distribution with parameters [315]zeta and [316]sigma, using the
formula given above.
_images/rand-lognormal.png
double gsl_cdf_lognormal_P(double x, double zeta, double sigma)[317]¶
double gsl_cdf_lognormal_Q(double x, double zeta, double sigma)[318]¶
double gsl_cdf_lognormal_Pinv(double P, double zeta, double sigma)[319]¶
double gsl_cdf_lognormal_Qinv(double Q, double zeta, double sigma)[320]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the lognormal distribution with parameters
[321]zeta and [322]sigma.
The Chi-squared Distribution[323]¶
The chi-squared distribution arises in statistics. If Y_i are n independent
Gaussian random variates with unit variance then the sum-of-squares,
X_i = \sum_i Y_i^2
has a chi-squared distribution with n degrees of freedom.
double gsl_ran_chisq(const [324]gsl_rng *r, double nu)[325]¶
This function returns a random variate from the chi-squared distribution
with [326]nu degrees of freedom. The distribution function is,
p(x) dx = {1 \over 2 \Gamma(\nu/2) } (x/2)^{\nu/2 - 1} \exp(-x/2) dx
for x \ge 0 .
double gsl_ran_chisq_pdf(double x, double nu)[327]¶
This function computes the probability density p(x) at [328]x for a
chi-squared distribution with [329]nu degrees of freedom, using the
formula given above.
_images/rand-chisq.png
double gsl_cdf_chisq_P(double x, double nu)[330]¶
double gsl_cdf_chisq_Q(double x, double nu)[331]¶
double gsl_cdf_chisq_Pinv(double P, double nu)[332]¶
double gsl_cdf_chisq_Qinv(double Q, double nu)[333]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the chi-squared distribution with [334]nu degrees
of freedom.
The F-distribution[335]¶
The F-distribution arises in statistics. If Y_1 and Y_2 are chi-squared deviates
with \nu_1 and \nu_2 degrees of freedom then the ratio,
X = { (Y_1 / \nu_1) \over (Y_2 / \nu_2) }
has an F-distribution F(x;\nu_1,\nu_2) .
double gsl_ran_fdist(const [336]gsl_rng *r, double nu1, double nu2)[337]¶
This function returns a random variate from the F-distribution with
degrees of freedom [338]nu1 and [339]nu2. The distribution function is,
p(x) dx = { \Gamma((\nu_1 + \nu_2)/2) \over \Gamma(\nu_1/2)
\Gamma(\nu_2/2) } \nu_1^{\nu_1/2} \nu_2^{\nu_2/2} x^{\nu_1/2 - 1} (\nu_2 +
\nu_1 x)^{-\nu_1/2 -\nu_2/2}
for x \ge 0 .
double gsl_ran_fdist_pdf(double x, double nu1, double nu2)[340]¶
This function computes the probability density p(x) at [341]x for an
F-distribution with [342]nu1 and [343]nu2 degrees of freedom, using the
formula given above.
_images/rand-fdist.png
double gsl_cdf_fdist_P(double x, double nu1, double nu2)[344]¶
double gsl_cdf_fdist_Q(double x, double nu1, double nu2)[345]¶
double gsl_cdf_fdist_Pinv(double P, double nu1, double nu2)[346]¶
double gsl_cdf_fdist_Qinv(double Q, double nu1, double nu2)[347]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the F-distribution with [348]nu1 and [349]nu2
degrees of freedom.
The t-distribution[350]¶
The t-distribution arises in statistics. If Y_1 has a normal distribution and Y_2
has a chi-squared distribution with \nu degrees of freedom then the ratio,
X = { Y_1 \over \sqrt{Y_2 / \nu} }
has a t-distribution t(x;\nu) with \nu degrees of freedom.
double gsl_ran_tdist(const [351]gsl_rng *r, double nu)[352]¶
This function returns a random variate from the t-distribution. The
distribution function is,
p(x) dx = {\Gamma((\nu + 1)/2) \over \sqrt{\pi \nu} \Gamma(\nu/2)} (1 +
x^2/\nu)^{-(\nu + 1)/2} dx
for -\infty < x < +\infty .
double gsl_ran_tdist_pdf(double x, double nu)[353]¶
This function computes the probability density p(x) at [354]x for a
t-distribution with [355]nu degrees of freedom, using the formula given
above.
_images/rand-tdist.png
double gsl_cdf_tdist_P(double x, double nu)[356]¶
double gsl_cdf_tdist_Q(double x, double nu)[357]¶
double gsl_cdf_tdist_Pinv(double P, double nu)[358]¶
double gsl_cdf_tdist_Qinv(double Q, double nu)[359]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the t-distribution with [360]nu degrees of freedom.
The Beta Distribution[361]¶
double gsl_ran_beta(const [362]gsl_rng *r, double a, double b)[363]¶
This function returns a random variate from the beta distribution. The
distribution function is,
p(x) dx = {\Gamma(a+b) \over \Gamma(a) \Gamma(b)} x^{a-1} (1-x)^{b-1} dx
for 0 \le x \le 1 .
double gsl_ran_beta_pdf(double x, double a, double b)[364]¶
This function computes the probability density p(x) at [365]x for a beta
distribution with parameters [366]a and [367]b, using the formula given
above.
_images/rand-beta.png
double gsl_cdf_beta_P(double x, double a, double b)[368]¶
double gsl_cdf_beta_Q(double x, double a, double b)[369]¶
double gsl_cdf_beta_Pinv(double P, double a, double b)[370]¶
double gsl_cdf_beta_Qinv(double Q, double a, double b)[371]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the beta distribution with parameters [372]a and
[373]b.
The Logistic Distribution[374]¶
double gsl_ran_logistic(const [375]gsl_rng *r, double a)[376]¶
This function returns a random variate from the logistic distribution. The
distribution function is,
p(x) dx = { \exp(-x/a) \over a (1 + \exp(-x/a))^2 } dx
for -\infty < x < +\infty .
double gsl_ran_logistic_pdf(double x, double a)[377]¶
This function computes the probability density p(x) at [378]x for a
logistic distribution with scale parameter [379]a, using the formula given
above.
_images/rand-logistic.png
double gsl_cdf_logistic_P(double x, double a)[380]¶
double gsl_cdf_logistic_Q(double x, double a)[381]¶
double gsl_cdf_logistic_Pinv(double P, double a)[382]¶
double gsl_cdf_logistic_Qinv(double Q, double a)[383]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the logistic distribution with scale parameter
[384]a.
The Pareto Distribution[385]¶
double gsl_ran_pareto(const [386]gsl_rng *r, double a, double b)[387]¶
This function returns a random variate from the Pareto distribution of
order [388]a. The distribution function is,
p(x) dx = (a/b) / (x/b)^{a+1} dx
for x \ge b .
double gsl_ran_pareto_pdf(double x, double a, double b)[389]¶
This function computes the probability density p(x) at [390]x for a Pareto
distribution with exponent [391]a and scale [392]b, using the formula
given above.
_images/rand-pareto.png
double gsl_cdf_pareto_P(double x, double a, double b)[393]¶
double gsl_cdf_pareto_Q(double x, double a, double b)[394]¶
double gsl_cdf_pareto_Pinv(double P, double a, double b)[395]¶
double gsl_cdf_pareto_Qinv(double Q, double a, double b)[396]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the Pareto distribution with exponent [397]a and
scale [398]b.
Spherical Vector Distributions[399]¶
The spherical distributions generate random vectors, located on a spherical
surface. They can be used as random directions, for example in the steps of a
random walk.
void gsl_ran_dir_2d(const [400]gsl_rng *r, double *x, double *y)[401]¶
void gsl_ran_dir_2d_trig_method(const [402]gsl_rng *r, double *x, double
*y)[403]¶
This function returns a random direction vector v = ([404]x, [405]y) in
two dimensions. The vector is normalized such that |v|^2 = x^2 + y^2 = 1 .
The obvious way to do this is to take a uniform random number between 0
and 2\pi and let [406]x and [407]y be the sine and cosine respectively.
Two trig functions would have been expensive in the old days, but with
modern hardware implementations, this is sometimes the fastest way to go.
This is the case for the Pentium (but not the case for the Sun
Sparcstation). One can avoid the trig evaluations by choosing [408]x and
[409]y in the interior of a unit circle (choose them at random from the
interior of the enclosing square, and then reject those that are outside
the unit circle), and then dividing by \sqrt{x^2 + y^2} . A much cleverer
approach, attributed to von Neumann (See Knuth, v2, 3rd ed, p140, exercise
23), requires neither trig nor a square root. In this approach, u and v
are chosen at random from the interior of a unit circle, and then
x=(u^2-v^2)/(u^2+v^2) and y=2uv/(u^2+v^2) .
void gsl_ran_dir_3d(const [410]gsl_rng *r, double *x, double *y, double *z)[411]¶
This function returns a random direction vector v = ([412]x, [413]y,
[414]z) in three dimensions. The vector is normalized such that |v|^2 =
x^2 + y^2 + z^2 = 1 . The method employed is due to Robert E. Knop (CACM
13, 326 (1970)), and explained in Knuth, v2, 3rd ed, p136. It uses the
surprising fact that the distribution projected along any axis is actually
uniform (this is only true for 3 dimensions).
void gsl_ran_dir_nd(const [415]gsl_rng *r, size_t n, double *x)[416]¶
This function returns a random direction vector v = (x_1,x_2,\ldots,x_n)
in [417]n dimensions. The vector is normalized such that |v|^2 = x_1^2 +
x_2^2 + \cdots + x_n^2 = 1 . The method uses the fact that a multivariate
Gaussian distribution is spherically symmetric. Each component is
generated to have a Gaussian distribution, and then the components are
normalized. The method is described by Knuth, v2, 3rd ed, p135-136, and
attributed to G. W. Brown, Modern Mathematics for the Engineer (1956).
The Weibull Distribution[418]¶
double gsl_ran_weibull(const [419]gsl_rng *r, double a, double b)[420]¶
This function returns a random variate from the Weibull distribution. The
distribution function is,
p(x) dx = {b \over a^b} x^{b-1} \exp(-(x/a)^b) dx
for x \ge 0 .
double gsl_ran_weibull_pdf(double x, double a, double b)[421]¶
This function computes the probability density p(x) at [422]x for a
Weibull distribution with scale [423]a and exponent [424]b, using the
formula given above.
_images/rand-weibull.png
double gsl_cdf_weibull_P(double x, double a, double b)[425]¶
double gsl_cdf_weibull_Q(double x, double a, double b)[426]¶
double gsl_cdf_weibull_Pinv(double P, double a, double b)[427]¶
double gsl_cdf_weibull_Qinv(double Q, double a, double b)[428]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the Weibull distribution with scale [429]a and
exponent [430]b.
The Type-1 Gumbel Distribution[431]¶
double gsl_ran_gumbel1(const [432]gsl_rng *r, double a, double b)[433]¶
This function returns a random variate from the Type-1 Gumbel
distribution. The Type-1 Gumbel distribution function is,
p(x) dx = a b \exp(-(b \exp(-ax) + ax)) dx
for -\infty < x < \infty .
double gsl_ran_gumbel1_pdf(double x, double a, double b)[434]¶
This function computes the probability density p(x) at [435]x for a Type-1
Gumbel distribution with parameters [436]a and [437]b, using the formula
given above.
_images/rand-gumbel1.png
double gsl_cdf_gumbel1_P(double x, double a, double b)[438]¶
double gsl_cdf_gumbel1_Q(double x, double a, double b)[439]¶
double gsl_cdf_gumbel1_Pinv(double P, double a, double b)[440]¶
double gsl_cdf_gumbel1_Qinv(double Q, double a, double b)[441]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the Type-1 Gumbel distribution with parameters
[442]a and [443]b.
The Type-2 Gumbel Distribution[444]¶
double gsl_ran_gumbel2(const [445]gsl_rng *r, double a, double b)[446]¶
This function returns a random variate from the Type-2 Gumbel
distribution. The Type-2 Gumbel distribution function is,
p(x) dx = a b x^{-a-1} \exp(-b x^{-a}) dx
for 0 < x < \infty .
double gsl_ran_gumbel2_pdf(double x, double a, double b)[447]¶
This function computes the probability density p(x) at [448]x for a Type-2
Gumbel distribution with parameters [449]a and [450]b, using the formula
given above.
_images/rand-gumbel2.png
double gsl_cdf_gumbel2_P(double x, double a, double b)[451]¶
double gsl_cdf_gumbel2_Q(double x, double a, double b)[452]¶
double gsl_cdf_gumbel2_Pinv(double P, double a, double b)[453]¶
double gsl_cdf_gumbel2_Qinv(double Q, double a, double b)[454]¶
These functions compute the cumulative distribution functions P(x) , Q(x)
and their inverses for the Type-2 Gumbel distribution with parameters
[455]a and [456]b.
The Dirichlet Distribution[457]¶
void gsl_ran_dirichlet(const [458]gsl_rng *r, size_t K, const double alpha[],
double theta[])[459]¶
This function returns an array of [460]K random variates from a Dirichlet
distribution of order [461]K-1. The distribution function is
p(\theta_1,\ldots,\theta_K) \, d\theta_1 \cdots d\theta_K = {1 \over Z}
\prod_{i=1}^{K} \theta_i^{\alpha_i - 1} \; \delta(1 -\sum_{i=1}^K
\theta_i) d\theta_1 \cdots d\theta_K
for \theta_i \ge 0 and \alpha_i > 0 . The delta function ensures that \sum
\theta_i = 1 . The normalization factor Z is
Z = {\prod_{i=1}^K \Gamma(\alpha_i) \over \Gamma( \sum_{i=1}^K \alpha_i)}
The random variates are generated by sampling [462]K values from gamma
distributions with parameters a=\alpha_i$, $b=1 , and renormalizing. See
A.M. Law, W.D. Kelton, Simulation Modeling and Analysis (1991).
double gsl_ran_dirichlet_pdf(size_t K, const double alpha[], const double
theta[])[463]¶
This function computes the probability density p(\theta_1, \ldots ,
\theta_K) at theta[K] for a Dirichlet distribution with parameters
alpha[K], using the formula given above.
double gsl_ran_dirichlet_lnpdf(size_t K, const double alpha[], const double
theta[])[464]¶
This function computes the logarithm of the probability density
p(\theta_1, \ldots , \theta_K) for a Dirichlet distribution with
parameters alpha[K].
General Discrete Distributions[465]¶
Given K discrete events with different probabilities P[k] , produce a random
value k consistent with its probability.
The obvious way to do this is to preprocess the probability list by generating a
cumulative probability array with K + 1 elements:
C[0] & = 0 \\ C[k+1] &= C[k] + P[k]
Note that this construction produces C[K] = 1 . Now choose a uniform deviate u
between 0 and 1, and find the value of k such that C[k] \le u < C[k+1] . Although
this in principle requires of order \log K steps per random number generation,
they are fast steps, and if you use something like \lfloor uK \rfloor as a
starting point, you can often do pretty well.
But faster methods have been devised. Again, the idea is to preprocess the
probability list, and save the result in some form of lookup table; then the
individual calls for a random discrete event can go rapidly. An approach invented
by G. Marsaglia (Generating discrete random variables in a computer, Comm ACM 6,
37-38 (1963)) is very clever, and readers interested in examples of good
algorithm design are directed to this short and well-written paper.
Unfortunately, for large K , Marsaglia's lookup table can be quite large.
A much better approach is due to Alastair J. Walker (An efficient method for
generating discrete random variables with general distributions, ACM Trans on
Mathematical Software 3, 253-256 (1977); see also Knuth, v2, 3rd ed,
p120-121,139). This requires two lookup tables, one floating point and one
integer, but both only of size K . After preprocessing, the random numbers are
generated in O(1) time, even for large K . The preprocessing suggested by Walker
requires O(K^2) effort, but that is not actually necessary, and the
implementation provided here only takes O(K) effort. In general, more
preprocessing leads to faster generation of the individual random numbers, but a
diminishing return is reached pretty early. Knuth points out that the optimal
preprocessing is combinatorially difficult for large K .
This method can be used to speed up some of the discrete random number generators
below, such as the binomial distribution. To use it for something like the
Poisson Distribution, a modification would have to be made, since it only takes a
finite set of K outcomes.
type gsl_ran_discrete_t[466]¶
This structure contains the lookup table for the discrete random number
generator.
[467]gsl_ran_discrete_t *gsl_ran_discrete_preproc(size_t K, const double
*P)[468]¶
This function returns a pointer to a structure that contains the lookup
table for the discrete random number generator. The array [469]P contains
the probabilities of the discrete events; these array elements must all be
positive, but they needn't add up to one (so you can think of them more
generally as "weights")--the preprocessor will normalize appropriately.
This return value is used as an argument for the [470]gsl_ran_discrete()
function below.
size_t gsl_ran_discrete(const [471]gsl_rng *r, const [472]gsl_ran_discrete_t
*g)[473]¶
After the preprocessor, above, has been called, you use this function to
get the discrete random numbers.
double gsl_ran_discrete_pdf(size_t k, const [474]gsl_ran_discrete_t *g)[475]¶
Returns the probability P[k] of observing the variable [476]k. Since P[k]
is not stored as part of the lookup table, it must be recomputed; this
computation takes O(K) , so if K is large and you care about the original
array P[k] used to create the lookup table, then you should just keep this
original array P[k] around.
void gsl_ran_discrete_free([477]gsl_ran_discrete_t *g)[478]¶
De-allocates the lookup table pointed to by [479]g.
The Poisson Distribution[480]¶
unsigned int gsl_ran_poisson(const [481]gsl_rng *r, double mu)[482]¶
This function returns a random integer from the Poisson distribution with
mean [483]mu. The probability distribution for Poisson variates is,
p(k) = {\mu^k \over k!} \exp(-\mu)
for k \ge 0 .
double gsl_ran_poisson_pdf(unsigned int k, double mu)[484]¶
This function computes the probability p(k) of obtaining [485]k from a
Poisson distribution with mean [486]mu, using the formula given above.
_images/rand-poisson.png
double gsl_cdf_poisson_P(unsigned int k, double mu)[487]¶
double gsl_cdf_poisson_Q(unsigned int k, double mu)[488]¶
These functions compute the cumulative distribution functions P(k) , Q(k)
for the Poisson distribution with parameter [489]mu.
The Bernoulli Distribution[490]¶
unsigned int gsl_ran_bernoulli(const [491]gsl_rng *r, double p)[492]¶
This function returns either 0 or 1, the result of a Bernoulli trial with
probability [493]p. The probability distribution for a Bernoulli trial is,
p(0) & = 1 - p \\ p(1) & = p
double gsl_ran_bernoulli_pdf(unsigned int k, double p)[494]¶
This function computes the probability p(k) of obtaining [495]k from a
Bernoulli distribution with probability parameter [496]p, using the
formula given above.
_images/rand-bernoulli.png
The Binomial Distribution[497]¶
unsigned int gsl_ran_binomial(const [498]gsl_rng *r, double p, unsigned int
n)[499]¶
This function returns a random integer from the binomial distribution, the
number of successes in [500]n independent trials with probability [501]p.
The probability distribution for binomial variates is,
p(k) = {n! \over k! (n-k)!} p^k (1-p)^{n-k}
for 0 \le k \le n .
double gsl_ran_binomial_pdf(unsigned int k, double p, unsigned int n)[502]¶
This function computes the probability p(k) of obtaining [503]k from a
binomial distribution with parameters [504]p and [505]n, using the formula
given above.
_images/rand-binomial.png
double gsl_cdf_binomial_P(unsigned int k, double p, unsigned int n)[506]¶
double gsl_cdf_binomial_Q(unsigned int k, double p, unsigned int n)[507]¶
These functions compute the cumulative distribution functions P(k) , Q(k)
for the binomial distribution with parameters [508]p and [509]n.
The Multinomial Distribution[510]¶
void gsl_ran_multinomial(const [511]gsl_rng *r, size_t K, unsigned int N, const
double p[], unsigned int n[])[512]¶
This function computes a random sample [513]n from the multinomial
distribution formed by [514]N trials from an underlying distribution p[K].
The distribution function for [515]n is,
P(n_1, n_2,\cdots, n_K) = {{ N!}\over{n_1 ! n_2 ! \cdots n_K !}} \,
p_1^{n_1} p_2^{n_2} \cdots p_K^{n_K}
where (n_1, n_2, \ldots, n_K) are nonnegative integers with \sum_{k=1}^{K}
n_k = N , and (p_1, p_2, \ldots, p_K) is a probability distribution with
\sum p_i = 1 . If the array p[K] is not normalized then its entries will
be treated as weights and normalized appropriately. The arrays [516]n and
[517]p must both be of length [518]K.
Random variates are generated using the conditional binomial method (see
C.S. Davis, The computer generation of multinomial random variates, Comp.
Stat. Data Anal. 16 (1993) 205-217 for details).
double gsl_ran_multinomial_pdf(size_t K, const double p[], const unsigned int
n[])[519]¶
This function computes the probability P(n_1, n_2, \ldots, n_K) of
sampling n[K] from a multinomial distribution with parameters p[K], using
the formula given above.
double gsl_ran_multinomial_lnpdf(size_t K, const double p[], const unsigned int
n[])[520]¶
This function returns the logarithm of the probability for the multinomial
distribution P(n_1, n_2, \ldots, n_K) with parameters p[K].
The Negative Binomial Distribution[521]¶
unsigned int gsl_ran_negative_binomial(const [522]gsl_rng *r, double p, double
n)[523]¶
This function returns a random integer from the negative binomial
distribution, the number of failures occurring before [524]n successes in
independent trials with probability [525]p of success. The probability
distribution for negative binomial variates is,
p(k) = {\Gamma(n + k) \over \Gamma(k+1) \Gamma(n) } p^n (1-p)^k
Note that n is not required to be an integer.
double gsl_ran_negative_binomial_pdf(unsigned int k, double p, double n)[526]¶
This function computes the probability p(k) of obtaining [527]k from a
negative binomial distribution with parameters [528]p and [529]n, using
the formula given above.
_images/rand-nbinomial.png
double gsl_cdf_negative_binomial_P(unsigned int k, double p, double n)[530]¶
double gsl_cdf_negative_binomial_Q(unsigned int k, double p, double n)[531]¶
These functions compute the cumulative distribution functions P(k) , Q(k)
for the negative binomial distribution with parameters [532]p and [533]n.
The Pascal Distribution[534]¶
unsigned int gsl_ran_pascal(const [535]gsl_rng *r, double p, unsigned int
n)[536]¶
This function returns a random integer from the Pascal distribution. The
Pascal distribution is simply a negative binomial distribution with an
integer value of n .
p(k) = {(n + k - 1)! \over k! (n - 1)! } p^n (1-p)^k
for k \ge 0 .
double gsl_ran_pascal_pdf(unsigned int k, double p, unsigned int n)[537]¶
This function computes the probability p(k) of obtaining [538]k from a
Pascal distribution with parameters [539]p and [540]n, using the formula
given above.
_images/rand-pascal.png
double gsl_cdf_pascal_P(unsigned int k, double p, unsigned int n)[541]¶
double gsl_cdf_pascal_Q(unsigned int k, double p, unsigned int n)[542]¶
These functions compute the cumulative distribution functions P(k) , Q(k)
for the Pascal distribution with parameters [543]p and [544]n.
The Geometric Distribution[545]¶
unsigned int gsl_ran_geometric(const [546]gsl_rng *r, double p)[547]¶
This function returns a random integer from the geometric distribution,
the number of independent trials with probability [548]p until the first
success. The probability distribution for geometric variates is,
p(k) = p (1-p)^{k-1}
for k \ge 1 . Note that the distribution begins with k = 1 with this
definition. There is another convention in which the exponent k - 1 is
replaced by k .
double gsl_ran_geometric_pdf(unsigned int k, double p)[549]¶
This function computes the probability p(k) of obtaining [550]k from a
geometric distribution with probability parameter [551]p, using the
formula given above.
_images/rand-geometric.png
double gsl_cdf_geometric_P(unsigned int k, double p)[552]¶
double gsl_cdf_geometric_Q(unsigned int k, double p)[553]¶
These functions compute the cumulative distribution functions P(k) , Q(k)
for the geometric distribution with parameter [554]p.
The Hypergeometric Distribution[555]¶
unsigned int gsl_ran_hypergeometric(const [556]gsl_rng *r, unsigned int n1,
unsigned int n2, unsigned int t)[557]¶
This function returns a random integer from the hypergeometric
distribution. The probability distribution for hypergeometric random
variates is,
p(k) = C(n_1, k) C(n_2, t - k) / C(n_1 + n_2, t)
where C(a,b) = a!/(b!(a-b)!) and t \leq n_1 + n_2 . The domain of k is
\max(0, t - n_2), \ldots, \min(t, n_1)
If a population contains n_1 elements of "type 1" and n_2 elements of
"type 2" then the hypergeometric distribution gives the probability of
obtaining k elements of "type 1" in t samples from the population without
replacement.
double gsl_ran_hypergeometric_pdf(unsigned int k, unsigned int n1, unsigned int
n2, unsigned int t)[558]¶
This function computes the probability p(k) of obtaining [559]k from a
hypergeometric distribution with parameters [560]n1, [561]n2, [562]t,
using the formula given above.
_images/rand-hypergeometric.png
double gsl_cdf_hypergeometric_P(unsigned int k, unsigned int n1, unsigned int n2,
unsigned int t)[563]¶
double gsl_cdf_hypergeometric_Q(unsigned int k, unsigned int n1, unsigned int n2,
unsigned int t)[564]¶
These functions compute the cumulative distribution functions P(k) , Q(k)
for the hypergeometric distribution with parameters [565]n1, [566]n2 and
[567]t.
The Logarithmic Distribution[568]¶
unsigned int gsl_ran_logarithmic(const [569]gsl_rng *r, double p)[570]¶
This function returns a random integer from the logarithmic distribution.
The probability distribution for logarithmic random variates is,
p(k) = {-1 \over \log(1-p)} {\left( p^k \over k \right)}
for k \ge 1 .
double gsl_ran_logarithmic_pdf(unsigned int k, double p)[571]¶
This function computes the probability p(k) of obtaining [572]k from a
logarithmic distribution with probability parameter [573]p, using the
formula given above.
_images/rand-logarithmic.png
The Wishart Distribution[574]¶
int gsl_ran_wishart(const [575]gsl_rng *r, const double n, const [576]gsl_matrix
*L, [577]gsl_matrix *result, [578]gsl_matrix *work)[579]¶
This function computes a random symmetric p -by- p matrix from the Wishart
distribution. The probability distribution for Wishart random variates is,
p(X) = \frac{|X|^{(n-p-1)/2} e^{-\textrm{tr}\left( V^{-1}
X\right)/2}}{2^{\frac{np}{2}} \left| V \right|^{n/2}
\Gamma_p(\frac{n}{2})}
Here, n > p - 1 is the number of degrees of freedom, V is a symmetric
positive definite p -by- p scale matrix, whose Cholesky factor is
specified by [580]L, and [581]work is p -by- p workspace. The p -by- p
Wishart distributed matrix X is stored in [582]result on output.
int gsl_ran_wishart_pdf(const [583]gsl_matrix *X, const [584]gsl_matrix *L_X,
const double n, const [585]gsl_matrix *L, double *result, [586]gsl_matrix
*work)[587]¶
int gsl_ran_wishart_log_pdf(const [588]gsl_matrix *X, const [589]gsl_matrix *L_X,
const double n, const [590]gsl_matrix *L, double *result, [591]gsl_matrix
*work)[592]¶
These functions compute p(X) or \log{p(X)} for the p -by- p matrix [593]X,
whose Cholesky factor is specified in [594]L_X. The degrees of freedom is
given by [595]n, the Cholesky factor of the scale matrix V is specified in
[596]L, and [597]work is p -by- p workspace. The probably density value is
returned in [598]result.
Shuffling and Sampling[599]¶
The following functions allow the shuffling and sampling of a set of objects. The
algorithms rely on a random number generator as a source of randomness and a poor
quality generator can lead to correlations in the output. In particular it is
important to avoid generators with a short period. For more information see
Knuth, v2, 3rd ed, Section 3.4.2, "Random Sampling and Shuffling".
void gsl_ran_shuffle(const [600]gsl_rng *r, void *base, size_t n, size_t
size)[601]¶
This function randomly shuffles the order of [602]n objects, each of size
[603]size, stored in the array base[0..n-1]. The output of the random
number generator [604]r is used to produce the permutation. The algorithm
generates all possible n! permutations with equal probability, assuming a
perfect source of random numbers.
The following code shows how to shuffle the numbers from 0 to 51:
int a[52];
for (i = 0; i < 52; i++)
{
a[i] = i;
}
gsl_ran_shuffle (r, a, 52, sizeof (int));
int gsl_ran_choose(const [605]gsl_rng *r, void *dest, size_t k, void *src, size_t
n, size_t size)[606]¶
This function fills the array dest[k] with [607]k objects taken randomly
from the [608]n elements of the array src[0..n-1]. The objects are each of
size [609]size. The output of the random number generator [610]r is used
to make the selection. The algorithm ensures all possible samples are
equally likely, assuming a perfect source of randomness.
The objects are sampled without replacement, thus each object can only
appear once in [611]dest. It is required that [612]k be less than or equal
to [613]n. The objects in [614]dest will be in the same relative order as
those in [615]src. You will need to call gsl_ran_shuffle(r, dest, n, size)
if you want to randomize the order.
The following code shows how to select a random sample of three unique
numbers from the set 0 to 99:
double a[3], b[100];
for (i = 0; i < 100; i++)
{
b[i] = (double) i;
}
gsl_ran_choose (r, a, 3, b, 100, sizeof (double));
void gsl_ran_sample(const [616]gsl_rng *r, void *dest, size_t k, void *src,
size_t n, size_t size)[617]¶
This function is like [618]gsl_ran_choose() but samples [619]k items from
the original array of [620]n items [621]src with replacement, so the same
object can appear more than once in the output sequence [622]dest. There
is no requirement that [623]k be less than [624]n in this case.
Examples[625]¶
The following program demonstrates the use of a random number generator to
produce variates from a distribution. It prints 10 samples from the Poisson
distribution with a mean of 3.
#include
#include
#include
int
main (void)
{
const gsl_rng_type * T;
gsl_rng * r;
int i, n = 10;
double mu = 3.0;
/* create a generator chosen by the
environment variable GSL_RNG_TYPE */
gsl_rng_env_setup();
T = gsl_rng_default;
r = gsl_rng_alloc (T);
/* print n random variates chosen from
the poisson distribution with mean
parameter mu */
for (i = 0; i < n; i++)
{
unsigned int k = gsl_ran_poisson (r, mu);
printf (" %u", k);
}
printf ("\n");
gsl_rng_free (r);
return 0;
}
If the library and header files are installed under /usr/local (the default
location) then the program can be compiled with these options:
$ gcc -Wall demo.c -lgsl -lgslcblas -lm
Here is the output of the program,
2 5 5 2 1 0 3 4 1 1
The variates depend on the seed used by the generator. The seed for the default
generator type [626]gsl_rng_default can be changed with the [627]GSL_RNG_SEED
environment variable to produce a different stream of variates:
$ GSL_RNG_SEED=123 ./a.out
giving output
4 5 6 3 3 1 4 2 5 5
The following program generates a random walk in two dimensions.
#include
#include
#include
int
main (void)
{
int i;
double x = 0, y = 0, dx, dy;
const gsl_rng_type * T;
gsl_rng * r;
gsl_rng_env_setup();
T = gsl_rng_default;
r = gsl_rng_alloc (T);
printf ("%g %g\n", x, y);
for (i = 0; i < 10; i++)
{
gsl_ran_dir_2d (r, &dx, &dy);
x += dx; y += dy;
printf ("%g %g\n", x, y);
}
gsl_rng_free (r);
return 0;
}
[628]Fig. 5 shows the output from the program.
[629]_images/random-walk.png
Fig. 5 Four 10-step random walks from the origin.[630]¶
The following program computes the upper and lower cumulative distribution
functions for the standard normal distribution at x = 2 .
#include
#include
int
main (void)
{
double P, Q;
double x = 2.0;
P = gsl_cdf_ugaussian_P (x);
printf ("prob(x < %f) = %f\n", x, P);
Q = gsl_cdf_ugaussian_Q (x);
printf ("prob(x > %f) = %f\n", x, Q);
x = gsl_cdf_ugaussian_Pinv (P);
printf ("Pinv(%f) = %f\n", P, x);
x = gsl_cdf_ugaussian_Qinv (Q);
printf ("Qinv(%f) = %f\n", Q, x);
return 0;
}
Here is the output of the program,
prob(x < 2.000000) = 0.977250
prob(x > 2.000000) = 0.022750
Pinv(0.977250) = 2.000000
Qinv(0.022750) = 2.000000
References and Further Reading[631]¶
For an encyclopaedic coverage of the subject readers are advised to consult the
book "Non-Uniform Random Variate Generation" by Luc Devroye. It covers every
imaginable distribution and provides hundreds of algorithms.
* Luc Devroye, "Non-Uniform Random Variate Generation", Springer-Verlag, ISBN
0-387-96305-7. Available online at
[632]http://cg.scs.carleton.ca/~luc/rnbookindex.html.
The subject of random variate generation is also reviewed by Knuth, who describes
algorithms for all the major distributions.
* Donald E. Knuth, "The Art of Computer Programming: Seminumerical Algorithms"
(Vol 2, 3rd Ed, 1997), Addison-Wesley, ISBN 0201896842.
The Particle Data Group provides a short review of techniques for generating
distributions of random numbers in the "Monte Carlo" section of its Annual Review
of Particle Physics.
* Review of Particle Properties, R.M. Barnett et al., Physical Review D54, 1
(1996) [633]http://pdg.lbl.gov/.
The Review of Particle Physics is available online in postscript and pdf format.
An overview of methods used to compute cumulative distribution functions can be
found in Statistical Computing by W.J. Kennedy and J.E. Gentle. Another general
reference is Elements of Statistical Computing by R.A. Thisted.
* William E. Kennedy and James E. Gentle, Statistical Computing (1980), Marcel
Dekker, ISBN 0-8247-6898-1.
* Ronald A. Thisted, Elements of Statistical Computing (1988), Chapman & Hall,
ISBN 0-412-01371-1.
The cumulative distribution functions for the Gaussian distribution are based on
the following papers,
* Rational Chebyshev Approximations Using Linear Equations, W.J. Cody, W.
Fraser, J.F. Hart. Numerische Mathematik 12, 242-251 (1968).
* Rational Chebyshev Approximations for the Error Function, W.J. Cody.
Mathematics of Computation 23, n107, 631-637 (July 1969).
[634]Next [635]Previous
____________________________________________________________________________
© Copyright 1996-2024 The GSL Team.
Built with [636]Sphinx using a [637]theme provided by [638]Read the Docs.
References
Visible links:
1. http://www.gnu.org/software/gsl/doc/html/genindex.html
2. http://www.gnu.org/software/gsl/doc/html/search.html
3. http://www.gnu.org/software/gsl/doc/html/statistics.html
4. http://www.gnu.org/software/gsl/doc/html/qrng.html
5. http://www.gnu.org/software/gsl/doc/html/index.html
6. http://www.gnu.org/software/gsl/doc/html/intro.html
7. http://www.gnu.org/software/gsl/doc/html/usage.html
8. http://www.gnu.org/software/gsl/doc/html/err.html
9. http://www.gnu.org/software/gsl/doc/html/math.html
10. http://www.gnu.org/software/gsl/doc/html/complex.html
11. http://www.gnu.org/software/gsl/doc/html/poly.html
12. http://www.gnu.org/software/gsl/doc/html/specfunc.html
13. http://www.gnu.org/software/gsl/doc/html/vectors.html
14. http://www.gnu.org/software/gsl/doc/html/permutation.html
15. http://www.gnu.org/software/gsl/doc/html/combination.html
16. http://www.gnu.org/software/gsl/doc/html/multiset.html
17. http://www.gnu.org/software/gsl/doc/html/sort.html
18. http://www.gnu.org/software/gsl/doc/html/blas.html
19. http://www.gnu.org/software/gsl/doc/html/linalg.html
20. http://www.gnu.org/software/gsl/doc/html/eigen.html
21. http://www.gnu.org/software/gsl/doc/html/fft.html
22. http://www.gnu.org/software/gsl/doc/html/integration.html
23. http://www.gnu.org/software/gsl/doc/html/rng.html
24. http://www.gnu.org/software/gsl/doc/html/qrng.html
25. http://www.gnu.org/software/gsl/doc/html/randist.html
26. http://www.gnu.org/software/gsl/doc/html/randist.html#introduction
27. http://www.gnu.org/software/gsl/doc/html/randist.html#the-gaussian-distribution
28. http://www.gnu.org/software/gsl/doc/html/randist.html#the-gaussian-tail-distribution
29. http://www.gnu.org/software/gsl/doc/html/randist.html#the-bivariate-gaussian-distribution
30. http://www.gnu.org/software/gsl/doc/html/randist.html#the-multivariate-gaussian-distribution
31. http://www.gnu.org/software/gsl/doc/html/randist.html#the-exponential-distribution
32. http://www.gnu.org/software/gsl/doc/html/randist.html#the-laplace-distribution
33. http://www.gnu.org/software/gsl/doc/html/randist.html#the-exponential-power-distribution
34. http://www.gnu.org/software/gsl/doc/html/randist.html#the-cauchy-distribution
35. http://www.gnu.org/software/gsl/doc/html/randist.html#the-rayleigh-distribution
36. http://www.gnu.org/software/gsl/doc/html/randist.html#the-rayleigh-tail-distribution
37. http://www.gnu.org/software/gsl/doc/html/randist.html#the-landau-distribution
38. http://www.gnu.org/software/gsl/doc/html/randist.html#the-levy-alpha-stable-distributions
39. http://www.gnu.org/software/gsl/doc/html/randist.html#the-levy-skew-alpha-stable-distribution
40. http://www.gnu.org/software/gsl/doc/html/randist.html#the-gamma-distribution
41. http://www.gnu.org/software/gsl/doc/html/randist.html#the-flat-uniform-distribution
42. http://www.gnu.org/software/gsl/doc/html/randist.html#the-lognormal-distribution
43. http://www.gnu.org/software/gsl/doc/html/randist.html#the-chi-squared-distribution
44. http://www.gnu.org/software/gsl/doc/html/randist.html#the-f-distribution
45. http://www.gnu.org/software/gsl/doc/html/randist.html#the-t-distribution
46. http://www.gnu.org/software/gsl/doc/html/randist.html#the-beta-distribution
47. http://www.gnu.org/software/gsl/doc/html/randist.html#the-logistic-distribution
48. http://www.gnu.org/software/gsl/doc/html/randist.html#the-pareto-distribution
49. http://www.gnu.org/software/gsl/doc/html/randist.html#spherical-vector-distributions
50. http://www.gnu.org/software/gsl/doc/html/randist.html#the-weibull-distribution
51. http://www.gnu.org/software/gsl/doc/html/randist.html#the-type-1-gumbel-distribution
52. http://www.gnu.org/software/gsl/doc/html/randist.html#the-type-2-gumbel-distribution
53. http://www.gnu.org/software/gsl/doc/html/randist.html#the-dirichlet-distribution
54. http://www.gnu.org/software/gsl/doc/html/randist.html#general-discrete-distributions
55. http://www.gnu.org/software/gsl/doc/html/randist.html#the-poisson-distribution
56. http://www.gnu.org/software/gsl/doc/html/randist.html#the-bernoulli-distribution
57. http://www.gnu.org/software/gsl/doc/html/randist.html#the-binomial-distribution
58. http://www.gnu.org/software/gsl/doc/html/randist.html#the-multinomial-distribution
59. http://www.gnu.org/software/gsl/doc/html/randist.html#the-negative-binomial-distribution
60. http://www.gnu.org/software/gsl/doc/html/randist.html#the-pascal-distribution
61. http://www.gnu.org/software/gsl/doc/html/randist.html#the-geometric-distribution
62. http://www.gnu.org/software/gsl/doc/html/randist.html#the-hypergeometric-distribution
63. http://www.gnu.org/software/gsl/doc/html/randist.html#the-logarithmic-distribution
64. http://www.gnu.org/software/gsl/doc/html/randist.html#the-wishart-distribution
65. http://www.gnu.org/software/gsl/doc/html/randist.html#shuffling-and-sampling
66. http://www.gnu.org/software/gsl/doc/html/randist.html#examples
67. http://www.gnu.org/software/gsl/doc/html/randist.html#references-and-further-reading
68. http://www.gnu.org/software/gsl/doc/html/statistics.html
69. http://www.gnu.org/software/gsl/doc/html/rstat.html
70. http://www.gnu.org/software/gsl/doc/html/movstat.html
71. http://www.gnu.org/software/gsl/doc/html/filter.html
72. http://www.gnu.org/software/gsl/doc/html/histogram.html
73. http://www.gnu.org/software/gsl/doc/html/ntuple.html
74. http://www.gnu.org/software/gsl/doc/html/montecarlo.html
75. http://www.gnu.org/software/gsl/doc/html/siman.html
76. http://www.gnu.org/software/gsl/doc/html/ode-initval.html
77. http://www.gnu.org/software/gsl/doc/html/interp.html
78. http://www.gnu.org/software/gsl/doc/html/diff.html
79. http://www.gnu.org/software/gsl/doc/html/cheb.html
80. http://www.gnu.org/software/gsl/doc/html/sum.html
81. http://www.gnu.org/software/gsl/doc/html/dwt.html
82. http://www.gnu.org/software/gsl/doc/html/dht.html
83. http://www.gnu.org/software/gsl/doc/html/roots.html
84. http://www.gnu.org/software/gsl/doc/html/min.html
85. http://www.gnu.org/software/gsl/doc/html/multiroots.html
86. http://www.gnu.org/software/gsl/doc/html/multimin.html
87. http://www.gnu.org/software/gsl/doc/html/lls.html
88. http://www.gnu.org/software/gsl/doc/html/nls.html
89. http://www.gnu.org/software/gsl/doc/html/bspline.html
90. http://www.gnu.org/software/gsl/doc/html/spmatrix.html
91. http://www.gnu.org/software/gsl/doc/html/spblas.html
92. http://www.gnu.org/software/gsl/doc/html/splinalg.html
93. http://www.gnu.org/software/gsl/doc/html/const.html
94. http://www.gnu.org/software/gsl/doc/html/ieee754.html
95. http://www.gnu.org/software/gsl/doc/html/debug.html
96. http://www.gnu.org/software/gsl/doc/html/contrib.html
97. http://www.gnu.org/software/gsl/doc/html/autoconf.html
98. http://www.gnu.org/software/gsl/doc/html/cblas.html
99. http://www.gnu.org/software/gsl/doc/html/gpl.html
100. http://www.gnu.org/software/gsl/doc/html/fdl.html
101. http://www.gnu.org/software/gsl/doc/html/index.html
102. http://www.gnu.org/software/gsl/doc/html/_sources/randist.rst.txt
103. http://www.gnu.org/software/gsl/doc/html/statistics.html
104. http://www.gnu.org/software/gsl/doc/html/qrng.html
105. http://www.gnu.org/software/gsl/doc/html/randist.html#random-number-distributions
106. http://www.gnu.org/software/gsl/doc/html/randist.html#introduction
107. http://www.gnu.org/software/gsl/doc/html/randist.html#the-gaussian-distribution
108. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
109. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian
110. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian
111. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian
112. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian
113. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian_pdf
114. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian_pdf
115. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian_pdf
116. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
117. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian_ziggurat
118. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
119. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian_ratio_method
120. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
121. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_ugaussian
122. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_ugaussian_pdf
123. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
124. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_ugaussian_ratio_method
125. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gaussian_P
126. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gaussian_Q
127. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gaussian_Pinv
128. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gaussian_Qinv
129. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gaussian_Qinv
130. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_ugaussian_P
131. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_ugaussian_Q
132. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_ugaussian_Pinv
133. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_ugaussian_Qinv
134. http://www.gnu.org/software/gsl/doc/html/randist.html#the-gaussian-tail-distribution
135. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
136. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian_tail
137. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian_tail
138. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian_tail
139. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian_tail_pdf
140. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian_tail_pdf
141. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian_tail_pdf
142. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gaussian_tail_pdf
143. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
144. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_ugaussian_tail
145. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_ugaussian_tail_pdf
146. http://www.gnu.org/software/gsl/doc/html/randist.html#the-bivariate-gaussian-distribution
147. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
148. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bivariate_gaussian
149. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bivariate_gaussian
150. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bivariate_gaussian
151. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bivariate_gaussian
152. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bivariate_gaussian
153. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bivariate_gaussian_pdf
154. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bivariate_gaussian_pdf
155. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bivariate_gaussian_pdf
156. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bivariate_gaussian_pdf
157. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bivariate_gaussian_pdf
158. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bivariate_gaussian_pdf
159. http://www.gnu.org/software/gsl/doc/html/randist.html#the-multivariate-gaussian-distribution
160. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
161. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_vector
162. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
163. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_vector
164. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian
165. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian
166. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian
167. http://www.gnu.org/software/gsl/doc/html/linalg.html#c.gsl_linalg_cholesky_decomp
168. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian
169. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_vector
170. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_vector
171. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
172. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_vector
173. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian_pdf
174. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_vector
175. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_vector
176. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
177. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_vector
178. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian_log_pdf
179. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian_log_pdf
180. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian_log_pdf
181. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian_log_pdf
182. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian_log_pdf
183. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
184. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_vector
185. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian_mean
186. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian_mean
187. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian_mean
188. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
189. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
190. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian_vcov
191. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian_vcov
192. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multivariate_gaussian_vcov
193. http://www.gnu.org/software/gsl/doc/html/randist.html#the-exponential-distribution
194. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
195. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_exponential
196. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_exponential
197. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_exponential_pdf
198. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_exponential_pdf
199. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_exponential_pdf
200. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_exponential_P
201. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_exponential_Q
202. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_exponential_Pinv
203. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_exponential_Qinv
204. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_exponential_Qinv
205. http://www.gnu.org/software/gsl/doc/html/randist.html#the-laplace-distribution
206. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
207. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_laplace
208. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_laplace
209. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_laplace_pdf
210. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_laplace_pdf
211. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_laplace_pdf
212. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_laplace_P
213. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_laplace_Q
214. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_laplace_Pinv
215. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_laplace_Qinv
216. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_laplace_Qinv
217. http://www.gnu.org/software/gsl/doc/html/randist.html#the-exponential-power-distribution
218. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
219. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_exppow
220. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_exppow
221. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_exppow
222. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_exppow_pdf
223. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_exppow_pdf
224. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_exppow_pdf
225. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_exppow_pdf
226. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_exppow_P
227. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_exppow_Q
228. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_exppow_Q
229. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_exppow_Q
230. http://www.gnu.org/software/gsl/doc/html/randist.html#the-cauchy-distribution
231. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
232. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_cauchy
233. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_cauchy
234. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_cauchy_pdf
235. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_cauchy_pdf
236. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_cauchy_pdf
237. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_cauchy_P
238. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_cauchy_Q
239. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_cauchy_Pinv
240. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_cauchy_Qinv
241. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_cauchy_Qinv
242. http://www.gnu.org/software/gsl/doc/html/randist.html#the-rayleigh-distribution
243. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
244. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_rayleigh
245. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_rayleigh
246. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_rayleigh_pdf
247. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_rayleigh_pdf
248. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_rayleigh_pdf
249. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_rayleigh_P
250. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_rayleigh_Q
251. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_rayleigh_Pinv
252. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_rayleigh_Qinv
253. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_rayleigh_Qinv
254. http://www.gnu.org/software/gsl/doc/html/randist.html#the-rayleigh-tail-distribution
255. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
256. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_rayleigh_tail
257. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_rayleigh_tail
258. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_rayleigh_tail
259. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_rayleigh_tail_pdf
260. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_rayleigh_tail_pdf
261. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_rayleigh_tail_pdf
262. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_rayleigh_tail_pdf
263. http://www.gnu.org/software/gsl/doc/html/randist.html#the-landau-distribution
264. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
265. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_landau
266. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_landau_pdf
267. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_landau_pdf
268. http://www.gnu.org/software/gsl/doc/html/randist.html#the-levy-alpha-stable-distributions
269. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
270. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_levy
271. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_levy
272. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_levy
273. http://www.gnu.org/software/gsl/doc/html/randist.html#the-levy-skew-alpha-stable-distribution
274. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
275. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_levy_skew
276. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_levy_skew
277. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_levy_skew
278. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_levy_skew
279. http://www.gnu.org/software/gsl/doc/html/randist.html#the-gamma-distribution
280. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
281. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gamma
282. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gamma
283. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
284. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gamma_knuth
285. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gamma_pdf
286. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gamma_pdf
287. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gamma_pdf
288. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gamma_pdf
289. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gamma_P
290. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gamma_Q
291. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gamma_Pinv
292. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gamma_Qinv
293. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gamma_Qinv
294. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gamma_Qinv
295. http://www.gnu.org/software/gsl/doc/html/randist.html#the-flat-uniform-distribution
296. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
297. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_flat
298. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_flat
299. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_flat
300. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_flat_pdf
301. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_flat_pdf
302. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_flat_pdf
303. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_flat_pdf
304. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_flat_P
305. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_flat_Q
306. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_flat_Pinv
307. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_flat_Qinv
308. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_flat_Qinv
309. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_flat_Qinv
310. http://www.gnu.org/software/gsl/doc/html/randist.html#the-lognormal-distribution
311. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
312. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_lognormal
313. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_lognormal_pdf
314. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_lognormal_pdf
315. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_lognormal_pdf
316. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_lognormal_pdf
317. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_lognormal_P
318. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_lognormal_Q
319. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_lognormal_Pinv
320. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_lognormal_Qinv
321. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_lognormal_Qinv
322. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_lognormal_Qinv
323. http://www.gnu.org/software/gsl/doc/html/randist.html#the-chi-squared-distribution
324. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
325. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_chisq
326. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_chisq
327. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_chisq_pdf
328. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_chisq_pdf
329. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_chisq_pdf
330. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_chisq_P
331. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_chisq_Q
332. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_chisq_Pinv
333. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_chisq_Qinv
334. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_chisq_Qinv
335. http://www.gnu.org/software/gsl/doc/html/randist.html#the-f-distribution
336. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
337. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_fdist
338. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_fdist
339. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_fdist
340. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_fdist_pdf
341. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_fdist_pdf
342. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_fdist_pdf
343. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_fdist_pdf
344. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_fdist_P
345. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_fdist_Q
346. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_fdist_Pinv
347. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_fdist_Qinv
348. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_fdist_Qinv
349. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_fdist_Qinv
350. http://www.gnu.org/software/gsl/doc/html/randist.html#the-t-distribution
351. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
352. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_tdist
353. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_tdist_pdf
354. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_tdist_pdf
355. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_tdist_pdf
356. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_tdist_P
357. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_tdist_Q
358. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_tdist_Pinv
359. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_tdist_Qinv
360. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_tdist_Qinv
361. http://www.gnu.org/software/gsl/doc/html/randist.html#the-beta-distribution
362. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
363. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_beta
364. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_beta_pdf
365. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_beta_pdf
366. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_beta_pdf
367. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_beta_pdf
368. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_beta_P
369. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_beta_Q
370. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_beta_Pinv
371. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_beta_Qinv
372. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_beta_Qinv
373. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_beta_Qinv
374. http://www.gnu.org/software/gsl/doc/html/randist.html#the-logistic-distribution
375. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
376. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_logistic
377. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_logistic_pdf
378. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_logistic_pdf
379. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_logistic_pdf
380. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_logistic_P
381. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_logistic_Q
382. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_logistic_Pinv
383. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_logistic_Qinv
384. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_logistic_Qinv
385. http://www.gnu.org/software/gsl/doc/html/randist.html#the-pareto-distribution
386. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
387. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_pareto
388. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_pareto
389. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_pareto_pdf
390. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_pareto_pdf
391. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_pareto_pdf
392. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_pareto_pdf
393. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_pareto_P
394. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_pareto_Q
395. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_pareto_Pinv
396. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_pareto_Qinv
397. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_pareto_Qinv
398. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_pareto_Qinv
399. http://www.gnu.org/software/gsl/doc/html/randist.html#spherical-vector-distributions
400. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
401. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_2d
402. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
403. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_2d_trig_method
404. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_2d_trig_method
405. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_2d_trig_method
406. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_2d_trig_method
407. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_2d_trig_method
408. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_2d_trig_method
409. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_2d_trig_method
410. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
411. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_3d
412. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_3d
413. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_3d
414. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_3d
415. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
416. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_nd
417. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dir_nd
418. http://www.gnu.org/software/gsl/doc/html/randist.html#the-weibull-distribution
419. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
420. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_weibull
421. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_weibull_pdf
422. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_weibull_pdf
423. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_weibull_pdf
424. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_weibull_pdf
425. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_weibull_P
426. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_weibull_Q
427. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_weibull_Pinv
428. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_weibull_Qinv
429. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_weibull_Qinv
430. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_weibull_Qinv
431. http://www.gnu.org/software/gsl/doc/html/randist.html#the-type-1-gumbel-distribution
432. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
433. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gumbel1
434. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gumbel1_pdf
435. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gumbel1_pdf
436. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gumbel1_pdf
437. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gumbel1_pdf
438. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gumbel1_P
439. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gumbel1_Q
440. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gumbel1_Pinv
441. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gumbel1_Qinv
442. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gumbel1_Qinv
443. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gumbel1_Qinv
444. http://www.gnu.org/software/gsl/doc/html/randist.html#the-type-2-gumbel-distribution
445. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
446. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gumbel2
447. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gumbel2_pdf
448. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gumbel2_pdf
449. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gumbel2_pdf
450. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_gumbel2_pdf
451. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gumbel2_P
452. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gumbel2_Q
453. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gumbel2_Pinv
454. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gumbel2_Qinv
455. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gumbel2_Qinv
456. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_gumbel2_Qinv
457. http://www.gnu.org/software/gsl/doc/html/randist.html#the-dirichlet-distribution
458. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
459. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dirichlet
460. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dirichlet
461. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dirichlet
462. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dirichlet
463. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dirichlet_pdf
464. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_dirichlet_lnpdf
465. http://www.gnu.org/software/gsl/doc/html/randist.html#general-discrete-distributions
466. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete_t
467. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete_t
468. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete_preproc
469. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete_preproc
470. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete
471. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
472. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete_t
473. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete
474. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete_t
475. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete_pdf
476. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete_pdf
477. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete_t
478. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete_free
479. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_discrete_free
480. http://www.gnu.org/software/gsl/doc/html/randist.html#the-poisson-distribution
481. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
482. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_poisson
483. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_poisson
484. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_poisson_pdf
485. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_poisson_pdf
486. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_poisson_pdf
487. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_poisson_P
488. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_poisson_Q
489. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_poisson_Q
490. http://www.gnu.org/software/gsl/doc/html/randist.html#the-bernoulli-distribution
491. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
492. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bernoulli
493. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bernoulli
494. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bernoulli_pdf
495. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bernoulli_pdf
496. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_bernoulli_pdf
497. http://www.gnu.org/software/gsl/doc/html/randist.html#the-binomial-distribution
498. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
499. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_binomial
500. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_binomial
501. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_binomial
502. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_binomial_pdf
503. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_binomial_pdf
504. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_binomial_pdf
505. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_binomial_pdf
506. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_binomial_P
507. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_binomial_Q
508. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_binomial_Q
509. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_binomial_Q
510. http://www.gnu.org/software/gsl/doc/html/randist.html#the-multinomial-distribution
511. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
512. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multinomial
513. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multinomial
514. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multinomial
515. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multinomial
516. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multinomial
517. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multinomial
518. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multinomial
519. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multinomial_pdf
520. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_multinomial_lnpdf
521. http://www.gnu.org/software/gsl/doc/html/randist.html#the-negative-binomial-distribution
522. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
523. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_negative_binomial
524. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_negative_binomial
525. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_negative_binomial
526. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_negative_binomial_pdf
527. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_negative_binomial_pdf
528. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_negative_binomial_pdf
529. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_negative_binomial_pdf
530. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_negative_binomial_P
531. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_negative_binomial_Q
532. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_negative_binomial_Q
533. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_negative_binomial_Q
534. http://www.gnu.org/software/gsl/doc/html/randist.html#the-pascal-distribution
535. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
536. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_pascal
537. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_pascal_pdf
538. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_pascal_pdf
539. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_pascal_pdf
540. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_pascal_pdf
541. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_pascal_P
542. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_pascal_Q
543. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_pascal_Q
544. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_pascal_Q
545. http://www.gnu.org/software/gsl/doc/html/randist.html#the-geometric-distribution
546. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
547. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_geometric
548. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_geometric
549. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_geometric_pdf
550. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_geometric_pdf
551. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_geometric_pdf
552. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_geometric_P
553. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_geometric_Q
554. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_geometric_Q
555. http://www.gnu.org/software/gsl/doc/html/randist.html#the-hypergeometric-distribution
556. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
557. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_hypergeometric
558. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_hypergeometric_pdf
559. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_hypergeometric_pdf
560. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_hypergeometric_pdf
561. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_hypergeometric_pdf
562. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_hypergeometric_pdf
563. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_hypergeometric_P
564. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_hypergeometric_Q
565. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_hypergeometric_Q
566. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_hypergeometric_Q
567. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_cdf_hypergeometric_Q
568. http://www.gnu.org/software/gsl/doc/html/randist.html#the-logarithmic-distribution
569. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
570. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_logarithmic
571. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_logarithmic_pdf
572. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_logarithmic_pdf
573. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_logarithmic_pdf
574. http://www.gnu.org/software/gsl/doc/html/randist.html#the-wishart-distribution
575. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
576. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
577. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
578. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
579. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_wishart
580. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_wishart
581. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_wishart
582. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_wishart
583. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
584. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
585. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
586. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
587. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_wishart_pdf
588. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
589. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
590. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
591. http://www.gnu.org/software/gsl/doc/html/vectors.html#c.gsl_matrix
592. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_wishart_log_pdf
593. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_wishart_log_pdf
594. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_wishart_log_pdf
595. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_wishart_log_pdf
596. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_wishart_log_pdf
597. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_wishart_log_pdf
598. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_wishart_log_pdf
599. http://www.gnu.org/software/gsl/doc/html/randist.html#shuffling-and-sampling
600. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
601. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_shuffle
602. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_shuffle
603. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_shuffle
604. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_shuffle
605. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
606. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_choose
607. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_choose
608. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_choose
609. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_choose
610. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_choose
611. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_choose
612. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_choose
613. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_choose
614. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_choose
615. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_choose
616. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng
617. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_sample
618. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_choose
619. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_sample
620. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_sample
621. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_sample
622. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_sample
623. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_sample
624. http://www.gnu.org/software/gsl/doc/html/randist.html#c.gsl_ran_sample
625. http://www.gnu.org/software/gsl/doc/html/randist.html#examples
626. http://www.gnu.org/software/gsl/doc/html/rng.html#c.gsl_rng_default
627. http://www.gnu.org/software/gsl/doc/html/rng.html#c.GSL_RNG_SEED
628. http://www.gnu.org/software/gsl/doc/html/randist.html#fig-rand-walk
629. http://www.gnu.org/software/gsl/doc/html/_images/random-walk.png
630. http://www.gnu.org/software/gsl/doc/html/randist.html#id1
631. http://www.gnu.org/software/gsl/doc/html/randist.html#references-and-further-reading
632. http://cg.scs.carleton.ca/~luc/rnbookindex.html
633. http://pdg.lbl.gov/
634. http://www.gnu.org/software/gsl/doc/html/statistics.html
635. http://www.gnu.org/software/gsl/doc/html/qrng.html
636. https://www.sphinx-doc.org/
637. https://github.com/readthedocs/sphinx_rtd_theme
638. https://readthedocs.org/
Hidden links:
640. http://www.gnu.org/software/gsl/doc/html/index.html
Usage: http://www.kk-software.de/kklynxview/get/URL
e.g. http://www.kk-software.de/kklynxview/get/http://www.kk-software.de
Errormessages are in German, sorry ;-)