PartMC  2.6.1
rand.F90
Go to the documentation of this file.
1 ! Copyright (C) 2007-2012 Matthew West
2 ! Licensed under the GNU General Public License version 2 or (at your
3 ! option) any later version. See the file COPYING for details.
4 
5 !> \file
6 !> The pmc_rand module.
7 
8 !> Random number generators.
9 module pmc_rand
10 
11  use pmc_util
12  use pmc_constants
13  use pmc_mpi
14 #ifdef PMC_USE_GSL
15  use iso_c_binding
16 #endif
17 
18  !> Length of a UUID string.
19  integer, parameter :: pmc_uuid_len = 36
20 
21  !> Result code indicating successful completion.
22  integer, parameter :: pmc_rand_gsl_success = 0
23  !> Result code indicating initialization failure.
24  integer, parameter :: pmc_rand_gsl_init_fail = 1
25  !> Result code indicating the generator was not initialized when it
26  !> should have been.
27  integer, parameter :: pmc_rand_gsl_not_init = 2
28  !> Result code indicating the generator was already initialized when
29  !> an initialization was attempted.
30  integer, parameter :: pmc_rand_gsl_already_init = 3
31 
32 contains
33 
34 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
35 
36 #ifdef PMC_USE_GSL
37  !> Check the return value of a call to one of the GSL RNG functions.
38  subroutine rand_check_gsl(code, value)
39 
40  !> Error code.
41  integer :: code
42  !> Return value.
43  integer(kind=c_int) :: value
44 
45  if (value == pmc_rand_gsl_success) then
46  return
47  elseif (value == pmc_rand_gsl_init_fail) then
48  call die_msg(code, "GSL RNG initialization failed")
49  elseif (value == pmc_rand_gsl_not_init) then
50  call die_msg(code, "GSL RNG has not been successfully initialized")
51  elseif (value == pmc_rand_gsl_already_init) then
52  call die_msg(code, "GSL RNG was already initialized")
53  else
54  call die_msg(code, "Unknown GSL RNG interface return value: " &
55  // trim(integer_to_string(value)))
56  end if
57 
58  end subroutine rand_check_gsl
59 #endif
60 
61 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
62 
63  !> Initializes the random number generator to the state defined by
64  !> the given seed plus offset. If the seed is 0 then a seed is
65  !> auto-generated from the current time plus offset.
66  subroutine pmc_srand(seed, offset)
67 
68  !> Random number generator seed.
69  integer, intent(in) :: seed
70  !> Random number generator offset.
71  integer, intent(in) :: offset
72 
73  integer :: i, n, clock
74  integer, allocatable :: seed_vec(:)
75 #ifdef PMC_USE_GSL
76  integer(kind=c_int) :: c_clock
77 #endif
78 
79 #ifdef PMC_USE_GSL
80 #ifndef DOXYGEN_SKIP_DOC
81  interface
82  integer(kind=c_int) function pmc_srand_gsl(seed) bind(c)
83  use iso_c_binding
84  integer(kind=c_int), value :: seed
85  end function pmc_srand_gsl
86  end interface
87 #endif
88 #endif
89 
90  if (seed == 0) then
91  if (pmc_mpi_rank() == 0) then
92  call system_clock(count = clock)
93  end if
94  ! ensure all nodes use exactly the same seed base, to avoid
95  ! accidental correlations
96  call pmc_mpi_bcast_integer(clock)
97  else
98  clock = seed
99  end if
100  clock = clock + 67 * offset
101 #ifdef PMC_USE_GSL
102  c_clock = int(clock, kind=c_int)
103  call rand_check_gsl(100489590, pmc_srand_gsl(c_clock))
104 #else
105  call random_seed(size = n)
106  allocate(seed_vec(n))
107  i = 0 ! HACK to shut up gfortran warning
108  seed_vec = clock + 37 * (/ (i - 1, i = 1, n) /)
109  call random_seed(put = seed_vec)
110  deallocate(seed_vec)
111 #endif
112 
113  end subroutine pmc_srand
114 
115 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
116 
117  !> Cleanup the random number generator.
118  subroutine pmc_rand_finalize()
119 
120 #ifdef PMC_USE_GSL
121 
122 #ifndef DOXYGEN_SKIP_DOC
123  interface
124  integer(kind=c_int) function pmc_rand_finalize_gsl() bind(c)
125  use iso_c_binding
126  end function pmc_rand_finalize_gsl
127  end interface
128 #endif
129 
130  call rand_check_gsl(489538382, pmc_rand_finalize_gsl())
131 #endif
132 
133  end subroutine pmc_rand_finalize
134 
135 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
136 
137  !> Returns a random number between 0 and 1.
138  real(kind=dp) function pmc_random()
139 
140 #ifdef PMC_USE_GSL
141  real(kind=c_double), target :: rnd
142  type(c_ptr) :: rnd_ptr
143 #else
144  real(kind=dp) :: rnd
145 #endif
146 
147 #ifdef PMC_USE_GSL
148 #ifndef DOXYGEN_SKIP_DOC
149  interface
150  integer(kind=c_int) function pmc_rand_gsl(harvest) bind(c)
151  use iso_c_binding
152  type(c_ptr), value :: harvest
153  end function pmc_rand_gsl
154  end interface
155 #endif
156 #endif
157 
158 #ifdef PMC_USE_GSL
159  rnd_ptr = c_loc(rnd)
160  call rand_check_gsl(843777138, pmc_rand_gsl(rnd_ptr))
161  pmc_random = real(rnd, kind=dp)
162 #else
163  call random_number(rnd)
164  pmc_random = rnd
165 #endif
166 
167  end function pmc_random
168 
169 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
170 
171  !> Returns a random integer between 1 and n.
172  integer function pmc_rand_int(n)
173 
174  !> Maximum random number to generate.
175  integer, intent(in) :: n
176 
177 #ifdef PMC_USE_GSL
178  integer(kind=c_int) :: n_c
179  integer(kind=c_int), target :: harvest
180  type(c_ptr) :: harvest_ptr
181 #endif
182 
183 #ifdef PMC_USE_GSL
184 #ifndef DOXYGEN_SKIP_DOC
185  interface
186  integer(kind=c_int) function pmc_rand_int_gsl(n, harvest) bind(c)
187  use iso_c_binding
188  integer(kind=c_int), value :: n
189  type(c_ptr), value :: harvest
190  end function pmc_rand_int_gsl
191  end interface
192 #endif
193 #endif
194 
195  call assert(669532625, n >= 1)
196 #ifdef PMC_USE_GSL
197  n_c = int(n, kind=c_int)
198  harvest_ptr = c_loc(harvest)
199  call rand_check_gsl(388234845, pmc_rand_int_gsl(n_c, harvest_ptr))
200  pmc_rand_int = int(harvest)
201 #else
202  pmc_rand_int = mod(int(pmc_random() * real(n, kind=dp)), n) + 1
203 #endif
204  call assert(515838689, pmc_rand_int >= 1)
205  call assert(802560153, pmc_rand_int <= n)
206 
207  end function pmc_rand_int
208 
209 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
210 
211  !> Round val to \c floor(val) or \c ceiling(val) with probability
212  !> proportional to the relative distance from \c val. That is,
213  !> Prob(prob_round(val) == floor(val)) = ceil(val) - val.
214  integer function prob_round(val)
215 
216  !> Value to round.
217  real(kind=dp), intent(in) :: val
218 
219  ! FIXME: can replace this with:
220  ! prob_round = floor(val + pmc_random())
221  if (pmc_random() < real(ceiling(val), kind=dp) - val) then
222  prob_round = floor(val)
223  else
224  prob_round = ceiling(val)
225  endif
226 
227  end function prob_round
228 
229 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
230 
231  !> Generate a Poisson-distributed random number with the given
232  !> mean.
233  !!
234  !! See http://en.wikipedia.org/wiki/Poisson_distribution for
235  !! information on the method. The method used at present is rather
236  !! inefficient and inaccurate (brute force for mean below 10 and
237  !! normal approximation above that point).
238  !!
239  !! The best known method appears to be due to Ahrens and Dieter (ACM
240  !! Trans. Math. Software, 8(2), 163-179, 1982) and is available (in
241  !! various forms) from:
242  !! - http://www.netlib.org/toms/599
243  !! - http://www.netlib.org/random/ranlib.f.tar.gz
244  !! - http://users.bigpond.net.au/amiller/random/random.f90
245  !! - http://www.netlib.org/random/random.f90
246  !!
247  !! Unfortunately the above code is under the non-free license:
248  !! - http://www.acm.org/pubs/copyright_policy/softwareCRnotice.html
249  !!
250  !! For other reasonable methods see L. Devroye, "Non-Uniform Random
251  !! Variate Generation", Springer-Verlag, 1986.
252  integer function rand_poisson(mean)
253 
254  !> Mean of the distribution.
255  real(kind=dp), intent(in) :: mean
256 
257 #ifdef PMC_USE_GSL
258  real(kind=c_double) :: mean_c
259  integer(kind=c_int), target :: harvest
260  type(c_ptr) :: harvest_ptr
261 #else
262  real(kind=dp) :: l, p
263  integer :: k
264 #endif
265 
266 #ifdef PMC_USE_GSL
267 #ifndef DOXYGEN_SKIP_DOC
268  interface
269  integer(kind=c_int) function pmc_rand_poisson_gsl(mean, harvest) &
270  bind(c)
271  use iso_c_binding
272  real(kind=c_double), value :: mean
273  type(c_ptr), value :: harvest
274  end function pmc_rand_poisson_gsl
275  end interface
276 #endif
277 #endif
278 
279  call assert(368397056, mean >= 0d0)
280 #ifdef PMC_USE_GSL
281  mean_c = real(mean, kind=c_double)
282  harvest_ptr = c_loc(harvest)
283  call rand_check_gsl(353483140, &
284  pmc_rand_poisson_gsl(mean_c, harvest_ptr))
285  rand_poisson = int(harvest)
286 #else
287  if (mean <= 10d0) then
288  ! exact method due to Knuth
289  l = exp(-mean)
290  k = 0
291  p = 1d0
292  do
293  k = k + 1
294  p = p * pmc_random()
295  if (p < l) exit
296  end do
297  rand_poisson = k - 1
298  else
299  ! normal approximation with a continuity correction
300  k = nint(rand_normal(mean, sqrt(mean)))
301  rand_poisson = max(k, 0)
302  end if
303 #endif
304 
305  end function rand_poisson
306 
307 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
308 
309  !> Generate a Binomial-distributed random number with the given
310  !> parameters.
311  !!
312  !! See http://en.wikipedia.org/wiki/Binomial_distribution for
313  !! information on the method. The method used at present is rather
314  !! inefficient and inaccurate (brute force for \f$np \le 10\f$ and
315  !! \f$n(1-p) \le 10\f$, otherwise normal approximation).
316  !!
317  !! For better methods, see L. Devroye, "Non-Uniform Random Variate
318  !! Generation", Springer-Verlag, 1986.
319  integer function rand_binomial(n, p)
320 
321  !> Sample size.
322  integer, intent(in) :: n
323  !> Sample probability.
324  real(kind=dp), intent(in) :: p
325 
326 #ifdef PMC_USE_GSL
327  integer(kind=c_int) :: n_c
328  real(kind=c_double) :: p_c
329  integer(kind=c_int), target :: harvest
330  type(c_ptr) :: harvest_ptr
331 #else
332  real(kind=dp) :: np, nomp, q, g_real
333  integer :: k, g, sum
334 #endif
335 
336 #ifdef PMC_USE_GSL
337 #ifndef DOXYGEN_SKIP_DOC
338  interface
339  integer(kind=c_int) function pmc_rand_binomial_gsl(n, p, harvest) &
340  bind(c)
341  use iso_c_binding
342  integer(kind=c_int), value :: n
343  real(kind=c_double), value :: p
344  type(c_ptr), value :: harvest
345  end function pmc_rand_binomial_gsl
346  end interface
347 #endif
348 #endif
349 
350  call assert(130699849, n >= 0)
351  call assert(754379796, p >= 0d0)
352  call assert(678506029, p <= 1d0)
353 #ifdef PMC_USE_GSL
354  n_c = int(n, kind=c_int)
355  p_c = real(p, kind=c_double)
356  harvest_ptr = c_loc(harvest)
357  call rand_check_gsl(208869397, &
358  pmc_rand_binomial_gsl(n_c, p_c, harvest_ptr))
359  rand_binomial = int(harvest)
360 #else
361  np = real(n, kind=dp) * p
362  nomp = real(n, kind=dp) * (1d0 - p)
363  if ((np >= 10d0) .and. (nomp >= 10d0)) then
364  ! normal approximation with continuity correction
365  k = nint(rand_normal(np, sqrt(np * (1d0 - p))))
366  rand_binomial = min(max(k, 0), n)
367  elseif (np < 1d-200) then
368  rand_binomial = 0
369  elseif (nomp < 1d-200) then
370  rand_binomial = n
371  else
372  ! First waiting time method of Devroye (p. 525).
373  ! We want p small, so if p > 1/2 then we compute with 1 - p and
374  ! take n - k at the end.
375  if (p <= 0.5d0) then
376  q = p
377  else
378  q = 1d0 - p
379  end if
380  k = 0
381  sum = 0
382  do
383  ! G is geometric(q)
384  g_real = log(pmc_random()) / log(1d0 - q)
385  ! early bailout for cases to avoid integer overflow
386  if (g_real > real(n - sum, kind=dp)) exit
387  g = ceiling(g_real)
388  sum = sum + g
389  if (sum > n) exit
390  k = k + 1
391  end do
392  if (p <= 0.5d0) then
393  rand_binomial = k
394  else
395  rand_binomial = n - k
396  end if
397  call assert(359087410, rand_binomial <= n)
398  end if
399 #endif
400 
401  end function rand_binomial
402 
403 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
404 
405  !> Generates a normally distributed random number with the given
406  !> mean and standard deviation.
407  real(kind=dp) function rand_normal(mean, stddev)
408 
409  !> Mean of distribution.
410  real(kind=dp), intent(in) :: mean
411  !> Standard deviation of distribution.
412  real(kind=dp), intent(in) :: stddev
413 
414 #ifdef PMC_USE_GSL
415  real(kind=c_double) :: mean_c, stddev_c
416  real(kind=c_double), target :: harvest
417  type(c_ptr) :: harvest_ptr
418 #else
419  real(kind=dp) :: u1, u2, r, theta, z0, z1
420 #endif
421 
422 #ifdef PMC_USE_GSL
423 #ifndef DOXYGEN_SKIP_DOC
424  interface
425  integer(kind=c_int) function pmc_rand_normal_gsl(mean, stddev, &
426  harvest) bind(c)
427  use iso_c_binding
428  real(kind=c_double), value :: mean
429  real(kind=c_double), value :: stddev
430  type(c_ptr), value :: harvest
431  end function pmc_rand_normal_gsl
432  end interface
433 #endif
434 #endif
435 
436  call assert(898978929, stddev >= 0d0)
437 #ifdef PMC_USE_GSL
438  mean_c = real(mean, kind=c_double)
439  stddev_c = real(stddev, kind=c_double)
440  harvest_ptr = c_loc(harvest)
441  call rand_check_gsl(102078576, &
442  pmc_rand_normal_gsl(mean_c, stddev_c, harvest_ptr))
443  rand_normal = real(harvest, kind=dp)
444 #else
445  ! Uses the Box-Muller transform
446  ! http://en.wikipedia.org/wiki/Box-Muller_transform
447  u1 = pmc_random()
448  u2 = pmc_random()
449  r = sqrt(-2d0 * log(u1))
450  theta = 2d0 * const%pi * u2
451  z0 = r * cos(theta)
452  z1 = r * sin(theta)
453  ! z0 and z1 are now independent N(0,1) random variables
454  ! We throw away z1, but we could use a SAVE variable to only do
455  ! the computation on every second call of this function.
456  rand_normal = stddev * z0 + mean
457 #endif
458 
459  end function rand_normal
460 
461 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
462 
463  !> Generates a vector of normally distributed random numbers with
464  !> the given means and standard deviations. This is a set of
465  !> normally distributed scalars, not a normally distributed vector.
466  subroutine rand_normal_array_1d(mean, stddev, val)
467 
468  !> Array of means.
469  real(kind=dp), intent(in) :: mean(:)
470  !> Array of standard deviations.
471  real(kind=dp), intent(in) :: stddev(size(mean))
472  !> Array of sampled normal random variables.
473  real(kind=dp), intent(out) :: val(size(mean))
474 
475  integer :: i
476 
477  do i = 1,size(mean)
478  val(i) = rand_normal(mean(i), stddev(i))
479  end do
480 
481  end subroutine rand_normal_array_1d
482 
483 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
484 
485  !> Sample the given continuous probability density function.
486  !!
487  !! That is, return a number k = 1,...,n such that prob(k) = pdf(k) /
488  !! sum(pdf). Uses accept-reject.
489  integer function sample_cts_pdf(pdf)
490 
491  !> Probability density function (not normalized).
492  real(kind=dp), intent(in) :: pdf(:)
493 
494  real(kind=dp) :: pdf_max
495  integer :: k
496  logical :: found
497 
498  ! use accept-reject
499  pdf_max = maxval(pdf)
500  if (minval(pdf) < 0d0) then
501  call die_msg(121838078, 'pdf contains negative values')
502  end if
503  if (pdf_max <= 0d0) then
504  call die_msg(119208863, 'pdf is not positive')
505  end if
506  found = .false.
507  do while (.not. found)
508  k = pmc_rand_int(size(pdf))
509  if (pmc_random() < pdf(k) / pdf_max) then
510  found = .true.
511  end if
512  end do
513  sample_cts_pdf = k
514 
515  end function sample_cts_pdf
516 
517 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
518 
519  !> Sample the given discrete probability density function.
520  !!
521  !! That is, return a number k = 1,...,n such that prob(k) = pdf(k) /
522  !! sum(pdf). Uses accept-reject.
523  integer function sample_disc_pdf(pdf)
524 
525  !> Probability density function.
526  integer, intent(in) :: pdf(:)
527 
528  integer :: pdf_max, k
529  logical :: found
530 
531  ! use accept-reject
532  pdf_max = maxval(pdf)
533  if (minval(pdf) < 0) then
534  call die_msg(598024763, 'pdf contains negative values')
535  end if
536  if (pdf_max <= 0) then
537  call die_msg(109961454, 'pdf is not positive')
538  end if
539  found = .false.
540  do while (.not. found)
541  k = pmc_rand_int(size(pdf))
542  if (pmc_random() < real(pdf(k), kind=dp) / real(pdf_max, kind=dp)) then
543  found = .true.
544  end if
545  end do
546  sample_disc_pdf = k
547 
548  end function sample_disc_pdf
549 
550 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
551 
552  !> Convert a real-valued vector into an integer-valued vector by
553  !> sampling.
554  !!
555  !! Use n_samp samples. Each discrete entry is sampled with a PDF
556  !! given by vec_cts. This is very slow for large n_samp or large n.
557  subroutine sample_vec_cts_to_disc(vec_cts, n_samp, vec_disc)
558 
559  !> Continuous vector.
560  real(kind=dp), intent(in) :: vec_cts(:)
561  !> Number of discrete samples to use.
562  integer, intent(in) :: n_samp
563  !> Discretized vector.
564  integer, intent(out) :: vec_disc(size(vec_cts))
565 
566  integer :: i_samp, k
567 
568  call assert(617770167, n_samp >= 0)
569  vec_disc = 0
570  do i_samp = 1,n_samp
571  k = sample_cts_pdf(vec_cts)
572  vec_disc(k) = vec_disc(k) + 1
573  end do
574 
575  end subroutine sample_vec_cts_to_disc
576 
577 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
578 
579  !> Generate a random hexadecimal character.
580  character function rand_hex_char()
581 
582  integer :: i
583 
584  i = pmc_rand_int(16)
585  if (i <= 10) then
586  rand_hex_char = achar(iachar('0') + i - 1)
587  else
588  rand_hex_char = achar(iachar('A') + i - 11)
589  end if
590 
591  end function rand_hex_char
592 
593 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
594 
595  !> Generate a version 4 UUID as a string.
596  !!
597  !! See http://en.wikipedia.org/wiki/Universally_Unique_Identifier
598  !! for format details.
599  subroutine uuid4_str(uuid)
600 
601  !> The newly generated UUID string.
602  character(len=PMC_UUID_LEN), intent(out) :: uuid
603 
604  integer :: i
605 
606  do i = 1,8
607  uuid(i:i) = rand_hex_char()
608  end do
609  uuid(9:9) = '-'
610  do i = 1,4
611  uuid((i + 9):(i + 9)) = rand_hex_char()
612  end do
613  uuid(14:14) = '-'
614  do i = 1,4
615  uuid((i + 14):(i + 14)) = rand_hex_char()
616  end do
617  uuid(19:19) = '-'
618  do i = 1,4
619  uuid((i + 19):(i + 19)) = rand_hex_char()
620  end do
621  uuid(24:24) = '-'
622  do i = 1,12
623  uuid((i + 24):(i + 24)) = rand_hex_char()
624  end do
625 
626  uuid(15:15) = '4'
627 
628  i = pmc_rand_int(4)
629  if (i <= 2) then
630  uuid(20:20) = achar(iachar('8') + i - 1)
631  else
632  uuid(20:20) = achar(iachar('A') + i - 3)
633  end if
634 
635  end subroutine uuid4_str
636 
637 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
638 
639 end module pmc_rand
pmc_rand::prob_round
integer function prob_round(val)
Round val to floor(val) or ceiling(val) with probability proportional to the relative distance from v...
Definition: rand.F90:215
pmc_rand::pmc_rand_gsl_init_fail
integer, parameter pmc_rand_gsl_init_fail
Result code indicating initialization failure.
Definition: rand.F90:24
pmc_mpi
Wrapper functions for MPI.
Definition: mpi.F90:13
pmc_srand_gsl
int pmc_srand_gsl(int seed)
Initialize the random number generator with the given seed.
Definition: rand_gsl.c:45
pmc_rand::sample_vec_cts_to_disc
subroutine sample_vec_cts_to_disc(vec_cts, n_samp, vec_disc)
Convert a real-valued vector into an integer-valued vector by sampling.
Definition: rand.F90:558
pmc_rand::rand_poisson
integer function rand_poisson(mean)
Generate a Poisson-distributed random number with the given mean.
Definition: rand.F90:253
pmc_rand::rand_normal
real(kind=dp) function rand_normal(mean, stddev)
Generates a normally distributed random number with the given mean and standard deviation.
Definition: rand.F90:408
pmc_util::die_msg
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:134
pmc_mpi::pmc_mpi_rank
integer function pmc_mpi_rank()
Returns the rank of the current process.
Definition: mpi.F90:117
pmc_constants::dp
integer, parameter dp
Kind of a double precision real number.
Definition: constants.F90:12
pmc_rand::pmc_rand_gsl_not_init
integer, parameter pmc_rand_gsl_not_init
Result code indicating the generator was not initialized when it should have been.
Definition: rand.F90:27
pmc_util::assert
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:103
pmc_rand_poisson_gsl
int pmc_rand_poisson_gsl(double mean, int *harvest)
Generate a Poisson-distributed random integer.
Definition: rand_gsl.c:126
pmc_rand::pmc_rand_gsl_already_init
integer, parameter pmc_rand_gsl_already_init
Result code indicating the generator was already initialized when an initialization was attempted.
Definition: rand.F90:30
pmc_rand::pmc_rand_int
integer function pmc_rand_int(n)
Returns a random integer between 1 and n.
Definition: rand.F90:173
pmc_rand::pmc_srand
subroutine pmc_srand(seed, offset)
Initializes the random number generator to the state defined by the given seed plus offset....
Definition: rand.F90:67
pmc_mpi::pmc_mpi_bcast_integer
subroutine pmc_mpi_bcast_integer(val)
Broadcast the given value from process 0 to all other processes.
Definition: mpi.F90:288
pmc_util::integer_to_string
character(len=pmc_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition: util.F90:766
pmc_constants::const
type(const_t), save const
Fixed variable for accessing the constant's values.
Definition: constants.F90:73
pmc_rand::pmc_rand_gsl_success
integer, parameter pmc_rand_gsl_success
Result code indicating successful completion.
Definition: rand.F90:22
pmc_rand::sample_cts_pdf
integer function sample_cts_pdf(pdf)
Sample the given continuous probability density function.
Definition: rand.F90:490
pmc_rand_finalize_gsl
int pmc_rand_finalize_gsl()
Cleanup and deallocate the random number generator.
Definition: rand_gsl.c:65
pmc_rand_int_gsl
int pmc_rand_int_gsl(int n, int *harvest)
Generate a uniform random integer in .
Definition: rand_gsl.c:95
pmc_rand
Random number generators.
Definition: rand.F90:9
pmc_rand::rand_normal_array_1d
subroutine rand_normal_array_1d(mean, stddev, val)
Generates a vector of normally distributed random numbers with the given means and standard deviation...
Definition: rand.F90:467
pmc_rand::uuid4_str
subroutine uuid4_str(uuid)
Generate a version 4 UUID as a string.
Definition: rand.F90:600
pmc_constants
Physical constants.
Definition: constants.F90:9
pmc_rand::pmc_uuid_len
integer, parameter pmc_uuid_len
Length of a UUID string.
Definition: rand.F90:19
pmc_util
Common utility subroutines.
Definition: util.F90:9
pmc_rand_gsl
int pmc_rand_gsl(double *harvest)
Generate a uniform random number in .
Definition: rand_gsl.c:80
pmc_rand::pmc_rand_finalize
subroutine pmc_rand_finalize()
Cleanup the random number generator.
Definition: rand.F90:119
pmc_rand::rand_binomial
integer function rand_binomial(n, p)
Generate a Binomial-distributed random number with the given parameters.
Definition: rand.F90:320
pmc_rand::rand_hex_char
character function rand_hex_char()
Generate a random hexadecimal character.
Definition: rand.F90:581
pmc_rand::sample_disc_pdf
integer function sample_disc_pdf(pdf)
Sample the given discrete probability density function.
Definition: rand.F90:524
pmc_rand::pmc_random
real(kind=dp) function pmc_random()
Returns a random number between 0 and 1.
Definition: rand.F90:139
pmc_rand_binomial_gsl
int pmc_rand_binomial_gsl(int n, double p, int *harvest)
Generate a Binomial-distributed random integer.
Definition: rand_gsl.c:142
pmc_rand_normal_gsl
int pmc_rand_normal_gsl(double mean, double stddev, double *harvest)
Generate a normally-distributed random number.
Definition: rand_gsl.c:111