|
subroutine | pmc_srand (seed, offset) |
| Initializes the random number generator to the state defined by the given seed plus offset. If the seed is 0 then a seed is auto-generated from the current time plus offset. More...
|
|
subroutine | pmc_rand_finalize () |
| Cleanup the random number generator. More...
|
|
real(kind=dp) function | pmc_random () |
| Returns a random number between 0 and 1. More...
|
|
integer function | pmc_rand_int (n) |
| Returns a random integer between 1 and n. More...
|
|
integer function | prob_round (val) |
| Round val to floor(val) or ceiling(val) with probability proportional to the relative distance from val . That is, Prob(prob_round(val) == floor(val)) = ceil(val) - val. More...
|
|
integer function | rand_poisson (mean) |
| Generate a Poisson-distributed random number with the given mean. More...
|
|
integer function | rand_binomial (n, p) |
| Generate a Binomial-distributed random number with the given parameters. More...
|
|
real(kind=dp) function | rand_normal (mean, stddev) |
| Generates a normally distributed random number with the given mean and standard deviation. More...
|
|
subroutine | rand_normal_array_1d (mean, stddev, val) |
| Generates a vector of normally distributed random numbers with the given means and standard deviations. This is a set of normally distributed scalars, not a normally distributed vector. More...
|
|
integer function | sample_cts_pdf (pdf) |
| Sample the given continuous probability density function. More...
|
|
integer function | sample_disc_pdf (pdf) |
| Sample the given discrete probability density function. More...
|
|
subroutine | sample_vec_cts_to_disc (vec_cts, n_samp, vec_disc) |
| Convert a real-valued vector into an integer-valued vector by sampling. More...
|
|
character function | rand_hex_char () |
| Generate a random hexadecimal character. More...
|
|
subroutine | uuid4_str (uuid) |
| Generate a version 4 UUID as a string. More...
|
|
Random number generators.
integer function pmc_rand::rand_binomial |
( |
integer, intent(in) |
n, |
|
|
real(kind=dp), intent(in) |
p |
|
) |
| |
Generate a Binomial-distributed random number with the given parameters.
See http://en.wikipedia.org/wiki/Binomial_distribution for information on the method. The method used at present is rather inefficient and inaccurate (brute force for and , otherwise normal approximation).
For better methods, see L. Devroye, "Non-Uniform Random Variate
Generation", Springer-Verlag, 1986.
- Parameters
-
[in] | n | Sample size. |
[in] | p | Sample probability. |
Definition at line 319 of file rand.F90.
integer function pmc_rand::rand_poisson |
( |
real(kind=dp), intent(in) |
mean | ) |
|
Generate a Poisson-distributed random number with the given mean.
See http://en.wikipedia.org/wiki/Poisson_distribution for information on the method. The method used at present is rather inefficient and inaccurate (brute force for mean below 10 and normal approximation above that point).
The best known method appears to be due to Ahrens and Dieter (ACM Trans. Math. Software, 8(2), 163-179, 1982) and is available (in various forms) from:
Unfortunately the above code is under the non-free license:
For other reasonable methods see L. Devroye, "Non-Uniform Random
Variate Generation", Springer-Verlag, 1986.
- Parameters
-
[in] | mean | Mean of the distribution. |
Definition at line 252 of file rand.F90.
integer function pmc_rand::sample_cts_pdf |
( |
real(kind=dp), dimension(:), intent(in) |
pdf | ) |
|
Sample the given continuous probability density function.
That is, return a number k = 1,...,n such that prob(k) = pdf(k) / sum(pdf). Uses accept-reject.
- Parameters
-
[in] | pdf | Probability density function (not normalized). |
Definition at line 489 of file rand.F90.
integer function pmc_rand::sample_disc_pdf |
( |
integer, dimension(:), intent(in) |
pdf | ) |
|
Sample the given discrete probability density function.
That is, return a number k = 1,...,n such that prob(k) = pdf(k) / sum(pdf). Uses accept-reject.
- Parameters
-
[in] | pdf | Probability density function. |
Definition at line 523 of file rand.F90.