PartMC  2.1.5
rand.F90
Go to the documentation of this file.
00001 ! Copyright (C) 2007-2011 Matthew West
00002 ! Licensed under the GNU General Public License version 2 or (at your
00003 ! option) any later version. See the file COPYING for details.
00004 
00005 !> \file
00006 !> The pmc_rand module.
00007 
00008 !> Random number generators.
00009 module pmc_rand
00010   
00011   use pmc_util
00012   use pmc_constants
00013   use pmc_mpi
00014 #ifdef PMC_USE_GSL
00015   use iso_c_binding
00016 #endif
00017 
00018   !> Length of a UUID string.
00019   integer, parameter :: PMC_UUID_LEN = 36
00020 
00021   !> Result code indicating successful completion.
00022   integer, parameter :: PMC_RAND_GSL_SUCCESS      = 0
00023   !> Result code indicating initialization failure.
00024   integer, parameter :: PMC_RAND_GSL_INIT_FAIL    = 1
00025   !> Result code indicating the generator was not initialized when it
00026   !> should have been.
00027   integer, parameter :: PMC_RAND_GSL_NOT_INIT     = 2
00028   !> Result code indicating the generator was already initialized when
00029   !> an initialization was attempted.
00030   integer, parameter :: PMC_RAND_GSL_ALREADY_INIT = 3
00031   
00032 contains
00033 
00034 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00035 
00036 #ifdef PMC_USE_GSL
00037   !> Check the return value of a call to one of the GSL RNG functions.
00038   subroutine rand_check_gsl(code, value)
00039 
00040     !> Error code.
00041     integer :: code
00042     !> Return value.
00043     integer(kind=c_int) :: value
00044 
00045     if (value == PMC_RAND_GSL_SUCCESS) then
00046        return
00047     elseif (value == PMC_RAND_GSL_INIT_FAIL) then
00048        call die_msg(code, "GSL RNG initialization failed")
00049     elseif (value == PMC_RAND_GSL_NOT_INIT) then
00050        call die_msg(code, "GSL RNG has not been successfully initialized")
00051     elseif (value == PMC_RAND_GSL_ALREADY_INIT) then
00052        call die_msg(code, "GSL RNG was already initialized")
00053     else
00054        call die_msg(code, "Unknown GSL RNG interface return value: " &
00055             // trim(integer_to_string(value)))
00056     end if
00057 
00058   end subroutine rand_check_gsl
00059 #endif
00060 
00061 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00062 
00063   !> Initializes the random number generator to the state defined by
00064   !> the given seed plus offset. If the seed is 0 then a seed is
00065   !> auto-generated from the current time plus offset.
00066   subroutine pmc_srand(seed, offset)
00067 
00068     !> Random number generator seed.
00069     integer, intent(in) :: seed
00070     !> Random number generator offset.
00071     integer, intent(in) :: offset
00072 
00073     integer :: i, n, clock
00074     integer, allocatable :: seed_vec(:)
00075 #ifdef PMC_USE_GSL
00076     integer(kind=c_int) :: c_clock
00077 #endif
00078 
00079 #ifdef PMC_USE_GSL
00080 #ifndef DOXYGEN_SKIP_DOC
00081     interface
00082        integer(kind=c_int) function pmc_srand_gsl(seed) bind(c)
00083          use iso_c_binding
00084          integer(kind=c_int), value :: seed
00085        end function pmc_srand_gsl
00086     end interface
00087 #endif
00088 #endif
00089 
00090     if (seed == 0) then
00091        if (pmc_mpi_rank() == 0) then
00092           call system_clock(count = clock)
00093        end if
00094        ! ensure all nodes use exactly the same seed base, to avoid
00095        ! accidental correlations
00096        call pmc_mpi_bcast_integer(clock)
00097     else
00098        clock = seed
00099     end if
00100     clock = clock + 67 * offset
00101 #ifdef PMC_USE_GSL
00102     c_clock = int(clock, kind=c_int)
00103     call rand_check_gsl(100489590, pmc_srand_gsl(c_clock))
00104 #else
00105     call random_seed(size = n)
00106     allocate(seed_vec(n))
00107     i = 0 ! HACK to shut up gfortran warning
00108     seed_vec = clock + 37 * (/ (i - 1, i = 1, n) /)
00109     call random_seed(put = seed_vec)
00110     deallocate(seed_vec)
00111 #endif
00112 
00113   end subroutine pmc_srand
00114 
00115 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00116 
00117   !> Cleanup the random number generator.
00118   subroutine pmc_rand_finalize()
00119 
00120 #ifdef PMC_USE_GSL
00121 
00122 #ifndef DOXYGEN_SKIP_DOC
00123     interface
00124        integer(kind=c_int) function pmc_rand_finalize_gsl() bind(c)
00125          use iso_c_binding
00126        end function pmc_rand_finalize_gsl
00127     end interface
00128 #endif
00129 
00130     call rand_check_gsl(489538382, pmc_rand_finalize_gsl())
00131 #endif
00132 
00133   end subroutine pmc_rand_finalize
00134 
00135 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00136 
00137   !> Returns a random number between 0 and 1.
00138   real(kind=dp) function pmc_random()
00139 
00140 #ifdef PMC_USE_GSL
00141     real(kind=c_double), target :: rnd
00142     type(c_ptr) :: rnd_ptr
00143 #else
00144     real(kind=dp) :: rnd
00145 #endif
00146 
00147 #ifdef PMC_USE_GSL
00148 #ifndef DOXYGEN_SKIP_DOC
00149     interface
00150        integer(kind=c_int) function pmc_rand_gsl(harvest) bind(c)
00151          use iso_c_binding
00152          type(c_ptr), value :: harvest
00153        end function pmc_rand_gsl
00154     end interface
00155 #endif
00156 #endif
00157 
00158 #ifdef PMC_USE_GSL
00159     rnd_ptr = c_loc(rnd)
00160     call rand_check_gsl(843777138, pmc_rand_gsl(rnd_ptr))
00161     pmc_random = real(rnd, kind=dp)
00162 #else
00163     call random_number(rnd)
00164     pmc_random = rnd
00165 #endif
00166 
00167   end function pmc_random
00168   
00169 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00170 
00171   !> Returns a random integer between 1 and n.
00172   integer function pmc_rand_int(n)
00173 
00174     !> Maximum random number to generate.
00175     integer, intent(in) :: n
00176 
00177 #ifdef PMC_USE_GSL
00178     integer(kind=c_int) :: n_c
00179     integer(kind=c_int), target :: harvest
00180     type(c_ptr) :: harvest_ptr
00181 #endif
00182 
00183 #ifdef PMC_USE_GSL
00184 #ifndef DOXYGEN_SKIP_DOC
00185     interface
00186        integer(kind=c_int) function pmc_rand_int_gsl(n, harvest) bind(c)
00187          use iso_c_binding
00188          integer(kind=c_int), value :: n
00189          type(c_ptr), value :: harvest
00190        end function pmc_rand_int_gsl
00191     end interface
00192 #endif
00193 #endif
00194 
00195     call assert(669532625, n >= 1)
00196 #ifdef PMC_USE_GSL
00197     n_c = int(n, kind=c_int)
00198     harvest_ptr = c_loc(harvest)
00199     call rand_check_gsl(388234845, pmc_rand_int_gsl(n_c, harvest_ptr))
00200     pmc_rand_int = int(harvest)
00201 #else
00202     pmc_rand_int = mod(int(pmc_random() * real(n, kind=dp)), n) + 1
00203 #endif
00204     call assert(515838689, pmc_rand_int >= 1)
00205     call assert(802560153, pmc_rand_int <= n)
00206 
00207   end function pmc_rand_int
00208 
00209 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00210 
00211   !> Round val to \c floor(val) or \c ceiling(val) with probability
00212   !> proportional to the relative distance from \c val. That is,
00213   !> Prob(prob_round(val) == floor(val)) = ceil(val) - val.
00214   integer function prob_round(val)
00215 
00216     !> Value to round.
00217     real(kind=dp), intent(in) :: val
00218 
00219     ! FIXME: can replace this with:
00220     ! prob_round = floor(val + pmc_random())
00221     if (pmc_random() < real(ceiling(val), kind=dp) - val) then
00222        prob_round = floor(val)
00223     else
00224        prob_round = ceiling(val)
00225     endif
00226 
00227   end function prob_round
00228 
00229 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00230 
00231   !> Generate a Poisson-distributed random number with the given
00232   !> mean.
00233   !!
00234   !! See http://en.wikipedia.org/wiki/Poisson_distribution for
00235   !! information on the method. The method used at present is rather
00236   !! inefficient and inaccurate (brute force for mean below 10 and
00237   !! normal approximation above that point).
00238   !!
00239   !! The best known method appears to be due to Ahrens and Dieter (ACM
00240   !! Trans. Math. Software, 8(2), 163-179, 1982) and is available (in
00241   !! various forms) from:
00242   !!     - http://www.netlib.org/toms/599
00243   !!     - http://www.netlib.org/random/ranlib.f.tar.gz
00244   !!     - http://users.bigpond.net.au/amiller/random/random.f90
00245   !!     - http://www.netlib.org/random/random.f90
00246   !!
00247   !! Unfortunately the above code is under the non-free license:
00248   !!     - http://www.acm.org/pubs/copyright_policy/softwareCRnotice.html
00249   !!
00250   !! For other reasonable methods see L. Devroye, "Non-Uniform Random
00251   !! Variate Generation", Springer-Verlag, 1986.
00252   integer function rand_poisson(mean)
00253 
00254     !> Mean of the distribution.
00255     real(kind=dp), intent(in) :: mean
00256 
00257 #ifdef PMC_USE_GSL
00258     real(kind=c_double) :: mean_c
00259     integer(kind=c_int), target :: harvest
00260     type(c_ptr) :: harvest_ptr
00261 #else
00262     real(kind=dp) :: L, p
00263     integer :: k
00264 #endif
00265 
00266 #ifdef PMC_USE_GSL
00267 #ifndef DOXYGEN_SKIP_DOC
00268     interface
00269        integer(kind=c_int) function pmc_rand_poisson_gsl(mean, harvest) &
00270             bind(c)
00271          use iso_c_binding
00272          real(kind=c_double), value :: mean
00273          type(c_ptr), value :: harvest
00274        end function pmc_rand_poisson_gsl
00275     end interface
00276 #endif
00277 #endif
00278 
00279     call assert(368397056, mean >= 0d0)
00280 #ifdef PMC_USE_GSL
00281     mean_c = real(mean, kind=c_double)
00282     harvest_ptr = c_loc(harvest)
00283     call rand_check_gsl(208869397, &
00284          pmc_rand_poisson_gsl(mean_c, harvest_ptr))
00285     rand_poisson = int(harvest)
00286 #else
00287     if (mean <= 10d0) then
00288        ! exact method due to Knuth
00289        L = exp(-mean)
00290        k = 0
00291        p = 1d0
00292        do
00293           k = k + 1
00294           p = p * pmc_random()
00295           if (p < L) exit
00296        end do
00297        rand_poisson = k - 1
00298     else
00299        ! normal approximation with a continuity correction
00300        k = nint(rand_normal(mean, sqrt(mean)))
00301        rand_poisson = max(k, 0)
00302     end if
00303 #endif
00304 
00305   end function rand_poisson
00306 
00307 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00308 
00309   !> Generate a Binomial-distributed random number with the given
00310   !> parameters.
00311   !!
00312   !! See http://en.wikipedia.org/wiki/Binomial_distribution for
00313   !! information on the method. The method used at present is rather
00314   !! inefficient and inaccurate (brute force for \f$np \le 10\f$ and
00315   !! \f$n(1-p) \le 10\f$, otherwise normal approximation).
00316   !!
00317   !! For better methods, see L. Devroye, "Non-Uniform Random Variate
00318   !! Generation", Springer-Verlag, 1986.
00319   integer function rand_binomial(n, p)
00320 
00321     !> Sample size.
00322     integer, intent(in) :: n
00323     !> Sample probability.
00324     real(kind=dp), intent(in) :: p
00325 
00326 #ifdef PMC_USE_GSL
00327     integer(kind=c_int) :: n_c
00328     real(kind=c_double) :: p_c
00329     integer(kind=c_int), target :: harvest
00330     type(c_ptr) :: harvest_ptr
00331 #else
00332     real(kind=dp) :: np, nomp, q, G_real
00333     integer :: k, G, sum
00334 #endif
00335 
00336 #ifdef PMC_USE_GSL
00337 #ifndef DOXYGEN_SKIP_DOC
00338     interface
00339        integer(kind=c_int) function pmc_rand_binomial_gsl(n, p, harvest) &
00340             bind(c)
00341          use iso_c_binding
00342          integer(kind=c_int), value :: n
00343          real(kind=c_double), value :: p
00344          type(c_ptr), value :: harvest
00345        end function pmc_rand_binomial_gsl
00346     end interface
00347 #endif
00348 #endif
00349 
00350     call assert(130699849, n >= 0)
00351     call assert(754379796, p >= 0d0)
00352     call assert(678506029, p <= 1d0)
00353 #ifdef PMC_USE_GSL
00354     n_c = int(n, kind=c_int)
00355     p_c = real(p, kind=c_double)
00356     harvest_ptr = c_loc(harvest)
00357     call rand_check_gsl(208869397, &
00358          pmc_rand_binomial_gsl(n_c, p_c, harvest_ptr))
00359     rand_binomial = int(harvest)
00360 #else
00361     np = real(n, kind=dp) * p
00362     nomp = real(n, kind=dp) * (1d0 - p)
00363     if ((np >= 10d0) .and. (nomp >= 10d0)) then
00364        ! normal approximation with continuity correction
00365        k = nint(rand_normal(np, sqrt(np * (1d0 - p))))
00366        rand_binomial = min(max(k, 0), n)
00367     elseif (np < 1d-200) then
00368        rand_binomial = 0
00369     elseif (nomp < 1d-200) then
00370        rand_binomial = n
00371     else
00372        ! First waiting time method of Devroye (p. 525).
00373        ! We want p small, so if p > 1/2 then we compute with 1 - p and
00374        ! take n - k at the end.
00375        if (p <= 0.5d0) then
00376           q = p
00377        else
00378           q = 1d0 - p
00379        end if
00380        k = 0
00381        sum = 0
00382        do
00383           ! G is geometric(q)
00384           G_real = log(pmc_random()) / log(1d0 - q)
00385           ! early bailout for cases to avoid integer overflow
00386           if (G_real > real(n - sum, kind=dp)) exit
00387           G = ceiling(G_real)
00388           sum = sum + G
00389           if (sum > n) exit
00390           k = k + 1
00391        end do
00392        if (p <= 0.5d0) then
00393           rand_binomial = k
00394        else
00395           rand_binomial = n - k
00396        end if
00397        call assert(359087410, rand_binomial <= n)
00398     end if
00399 #endif
00400 
00401   end function rand_binomial
00402 
00403 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00404 
00405   !> Generates a normally distributed random number with the given
00406   !> mean and standard deviation.
00407   real(kind=dp) function rand_normal(mean, stddev)
00408 
00409     !> Mean of distribution.
00410     real(kind=dp), intent(in) :: mean
00411     !> Standard deviation of distribution.
00412     real(kind=dp), intent(in) :: stddev
00413 
00414 #ifdef PMC_USE_GSL
00415     real(kind=c_double) :: mean_c, stddev_c
00416     real(kind=c_double), target :: harvest
00417     type(c_ptr) :: harvest_ptr
00418 #else
00419     real(kind=dp) :: u1, u2, r, theta, z0, z1
00420 #endif
00421 
00422 #ifdef PMC_USE_GSL
00423 #ifndef DOXYGEN_SKIP_DOC
00424     interface
00425        integer(kind=c_int) function pmc_rand_normal_gsl(mean, stddev, &
00426             harvest) bind(c)
00427          use iso_c_binding
00428          real(kind=c_double), value :: mean
00429          real(kind=c_double), value :: stddev
00430          type(c_ptr), value :: harvest
00431        end function pmc_rand_normal_gsl
00432     end interface
00433 #endif
00434 #endif
00435 
00436     call assert(898978929, stddev >= 0d0)
00437 #ifdef PMC_USE_GSL
00438     mean_c = real(mean, kind=c_double)
00439     stddev_c = real(stddev, kind=c_double)
00440     harvest_ptr = c_loc(harvest)
00441     call rand_check_gsl(102078576, &
00442          pmc_rand_normal_gsl(mean_c, stddev_c, harvest_ptr))
00443     rand_normal = real(harvest, kind=dp)
00444 #else
00445     ! Uses the Box-Muller transform
00446     ! http://en.wikipedia.org/wiki/Box-Muller_transform
00447     u1 = pmc_random()
00448     u2 = pmc_random()
00449     r = sqrt(-2d0 * log(u1))
00450     theta = 2d0 * const%pi * u2
00451     z0 = r * cos(theta)
00452     z1 = r * sin(theta)
00453     ! z0 and z1 are now independent N(0,1) random variables
00454     ! We throw away z1, but we could use a SAVE variable to only do
00455     ! the computation on every second call of this function.
00456     rand_normal = stddev * z0 + mean
00457 #endif
00458 
00459   end function rand_normal
00460 
00461 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00462 
00463   !> Sample the given continuous probability density function.
00464   !!
00465   !! That is, return a number k = 1,...,n such that prob(k) = pdf(k) /
00466   !! sum(pdf). Uses accept-reject.
00467   integer function sample_cts_pdf(pdf)
00468 
00469     !> Probability density function (not normalized).
00470     real(kind=dp), intent(in) :: pdf(:)
00471 
00472     real(kind=dp) :: pdf_max
00473     integer :: k
00474     logical :: found
00475 
00476     ! use accept-reject
00477     pdf_max = maxval(pdf)
00478     if (minval(pdf) < 0d0) then
00479        call die_msg(121838078, 'pdf contains negative values')
00480     end if
00481     if (pdf_max <= 0d0) then
00482        call die_msg(119208863, 'pdf is not positive')
00483     end if
00484     found = .false.
00485     do while (.not. found)
00486        k = pmc_rand_int(size(pdf))
00487        if (pmc_random() < pdf(k) / pdf_max) then
00488           found = .true.
00489        end if
00490     end do
00491     sample_cts_pdf = k
00492 
00493   end function sample_cts_pdf
00494 
00495 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00496 
00497   !> Sample the given discrete probability density function.
00498   !!
00499   !! That is, return a number k = 1,...,n such that prob(k) = pdf(k) /
00500   !! sum(pdf). Uses accept-reject.
00501   integer function sample_disc_pdf(pdf)
00502 
00503     !> Probability density function.
00504     integer, intent(in) :: pdf(:)
00505 
00506     integer :: pdf_max, k
00507     logical :: found
00508 
00509     ! use accept-reject
00510     pdf_max = maxval(pdf)
00511     if (minval(pdf) < 0) then
00512        call die_msg(598024763, 'pdf contains negative values')
00513     end if
00514     if (pdf_max <= 0) then
00515        call die_msg(109961454, 'pdf is not positive')
00516     end if
00517     found = .false.
00518     do while (.not. found)
00519        k = pmc_rand_int(size(pdf))
00520        if (pmc_random() < real(pdf(k), kind=dp) / real(pdf_max, kind=dp)) then
00521           found = .true.
00522        end if
00523     end do
00524     sample_disc_pdf = k
00525 
00526   end function sample_disc_pdf
00527   
00528 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00529 
00530   !> Convert a real-valued vector into an integer-valued vector by
00531   !> sampling.
00532   !!
00533   !! Use n_samp samples. Each discrete entry is sampled with a PDF
00534   !! given by vec_cts. This is very slow for large n_samp or large n.
00535   subroutine sample_vec_cts_to_disc(vec_cts, n_samp, vec_disc)
00536     
00537     !> Continuous vector.
00538     real(kind=dp), intent(in) :: vec_cts(:)
00539     !> Number of discrete samples to use.
00540     integer, intent(in) :: n_samp
00541     !> Discretized vector.
00542     integer, intent(out) :: vec_disc(size(vec_cts))
00543 
00544     integer :: i_samp, k
00545 
00546     call assert(617770167, n_samp >= 0)
00547     vec_disc = 0
00548     do i_samp = 1,n_samp
00549        k = sample_cts_pdf(vec_cts)
00550        vec_disc(k) = vec_disc(k) + 1
00551     end do
00552     
00553   end subroutine sample_vec_cts_to_disc
00554   
00555 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00556 
00557   !> Generate a random hexadecimal character.
00558   character function rand_hex_char()
00559 
00560     integer :: i
00561 
00562     i = pmc_rand_int(16)
00563     if (i <= 10) then
00564        rand_hex_char = achar(iachar('0') + i - 1)
00565     else
00566        rand_hex_char = achar(iachar('A') + i - 11)
00567     end if
00568 
00569   end function rand_hex_char
00570 
00571 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00572 
00573   !> Generate a version 4 UUID as a string.
00574   !!
00575   !! See http://en.wikipedia.org/wiki/Universally_Unique_Identifier
00576   !! for format details.
00577   subroutine uuid4_str(uuid)
00578 
00579     !> The newly generated UUID string.
00580     character(len=PMC_UUID_LEN), intent(out) :: uuid
00581 
00582     integer :: i
00583 
00584     do i = 1,8
00585        uuid(i:i) = rand_hex_char()
00586     end do
00587     uuid(9:9) = '-'
00588     do i = 1,4
00589        uuid((i + 9):(i + 9)) = rand_hex_char()
00590     end do
00591     uuid(14:14) = '-'
00592     do i = 1,4
00593        uuid((i + 14):(i + 14)) = rand_hex_char()
00594     end do
00595     uuid(19:19) = '-'
00596     do i = 1,4
00597        uuid((i + 19):(i + 19)) = rand_hex_char()
00598     end do
00599     uuid(24:24) = '-'
00600     do i = 1,12
00601        uuid((i + 24):(i + 24)) = rand_hex_char()
00602     end do
00603 
00604     uuid(15:15) = '4'
00605 
00606     i = pmc_rand_int(4)
00607     if (i <= 2) then
00608        uuid(20:20) = achar(iachar('8') + i - 1)
00609     else
00610        uuid(20:20) = achar(iachar('A') + i - 3)
00611     end if
00612 
00613   end subroutine uuid4_str
00614 
00615 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00616 
00617 end module pmc_rand