| 1 | ! radiation_random_numbers.F90 - Generate random numbers for McICA solver |
|---|
| 2 | ! |
|---|
| 3 | ! (C) Copyright 2020- ECMWF. |
|---|
| 4 | ! |
|---|
| 5 | ! This software is licensed under the terms of the Apache Licence Version 2.0 |
|---|
| 6 | ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0. |
|---|
| 7 | ! |
|---|
| 8 | ! In applying this licence, ECMWF does not waive the privileges and immunities |
|---|
| 9 | ! granted to it by virtue of its status as an intergovernmental organisation |
|---|
| 10 | ! nor does it submit to any jurisdiction. |
|---|
| 11 | ! |
|---|
| 12 | ! Author: Robin Hogan |
|---|
| 13 | ! Email: r.j.hogan@ecmwf.int |
|---|
| 14 | ! |
|---|
| 15 | ! The derived type "rng_type" is a random number generator that uses |
|---|
| 16 | ! either (1) Fortran's built-in random_number function, or (2) a |
|---|
| 17 | ! vectorized version of the MINSTD linear congruential generator. In |
|---|
| 18 | ! the case of (2), an rng_type object is initialized with a seed that |
|---|
| 19 | ! is used to fill up a state of "nmaxstreams" elements using the C++ |
|---|
| 20 | ! minstd_rand0 version of the MINSTD linear congruential generator |
|---|
| 21 | ! (LNG), which has the form istate[i+1] = mod(istate[i]*A0, M) from |
|---|
| 22 | ! i=1 to i=nmaxstreams. Subsequent requests for blocks of nmaxstreams |
|---|
| 23 | ! of random numbers use the C++ minstd_ran algorithm in a vectorizable |
|---|
| 24 | ! form, which modifies the state elements via istate[i] <- |
|---|
| 25 | ! mod(istate[i]*A, M). Uniform deviates are returned that normalize |
|---|
| 26 | ! the state elements to the range 0-1. |
|---|
| 27 | ! |
|---|
| 28 | ! The MINSTD generator was coded because the random_numbers_mix |
|---|
| 29 | ! generator in the IFS was found not to vectorize well on some |
|---|
| 30 | ! hardware. I am no expert on random number generators, so my |
|---|
| 31 | ! implementation should really be looked at and improved by someone |
|---|
| 32 | ! who knows what they are doing. |
|---|
| 33 | ! |
|---|
| 34 | ! Reference for MINSTD: Park, Stephen K.; Miller, Keith |
|---|
| 35 | ! W. (1988). "Random Number Generators: Good Ones Are Hard To Find" |
|---|
| 36 | ! (PDF). Communications of the ACM. 31 (10): |
|---|
| 37 | ! 1192-1201. doi:10.1145/63039.63042 |
|---|
| 38 | ! |
|---|
| 39 | ! Modifications |
|---|
| 40 | ! 2022-12-01 R. Hogan Fixed zeroed state in single precision |
|---|
| 41 | |
|---|
| 42 | module radiation_random_numbers |
|---|
| 43 | |
|---|
| 44 | use parkind1, only : jprb, jprd, jpim, jpib |
|---|
| 45 | |
|---|
| 46 | implicit none |
|---|
| 47 | |
|---|
| 48 | public :: rng_type, IRngMinstdVector, IRngNative |
|---|
| 49 | |
|---|
| 50 | enum, bind(c) |
|---|
| 51 | enumerator IRngNative, & ! Built-in Fortran-90 RNG |
|---|
| 52 | & IRngMinstdVector ! Vector MINSTD algorithm |
|---|
| 53 | end enum |
|---|
| 54 | |
|---|
| 55 | ! Maximum number of random numbers that can be computed in one call |
|---|
| 56 | ! - this can be increased |
|---|
| 57 | integer(kind=jpim), parameter :: NMaxStreams = 512 |
|---|
| 58 | |
|---|
| 59 | ! A requirement of the generator is that the operation mod(A*X,M) is |
|---|
| 60 | ! performed with no loss of precision, so type used for A and X must |
|---|
| 61 | ! be able to hold the largest possible value of A*X without |
|---|
| 62 | ! overflowing, going negative or losing precision. The largest |
|---|
| 63 | ! possible value is 48271*2147483647 = 103661183124337. This number |
|---|
| 64 | ! can be held in either a double-precision real number, or an 8-byte |
|---|
| 65 | ! integer. Either may be used, but on some hardwares it has been |
|---|
| 66 | ! found that operations on double-precision reals are faster. Select |
|---|
| 67 | ! which you prefer by defining USE_REAL_RNG_STATE for double |
|---|
| 68 | ! precision, or undefining it for an 8-byte integer. |
|---|
| 69 | #define USE_REAL_RNG_STATE 1 |
|---|
| 70 | |
|---|
| 71 | ! Define RNG_STATE_TYPE based on USE_REAL_RNG_STATE, where jprd |
|---|
| 72 | ! refers to a double-precision number regardless of the working |
|---|
| 73 | ! precision described by jprb, while jpib describes an 8-byte |
|---|
| 74 | ! integer |
|---|
| 75 | #ifdef USE_REAL_RNG_STATE |
|---|
| 76 | #define RNG_STATE_TYPE real(kind=jprd) |
|---|
| 77 | #else |
|---|
| 78 | #define RNG_STATE_TYPE integer(kind=jpib) |
|---|
| 79 | #endif |
|---|
| 80 | |
|---|
| 81 | ! The constants used in the main random number generator |
|---|
| 82 | RNG_STATE_TYPE , parameter :: IMinstdA = 48271 |
|---|
| 83 | RNG_STATE_TYPE , parameter :: IMinstdM = 2147483647 |
|---|
| 84 | |
|---|
| 85 | ! An alternative value of A that can be used to initialize the |
|---|
| 86 | ! members of the state from a single seed |
|---|
| 87 | RNG_STATE_TYPE , parameter :: IMinstdA0 = 16807 |
|---|
| 88 | |
|---|
| 89 | ! Scaling to convert the state to a uniform deviate in the range 0 |
|---|
| 90 | ! to 1 in working precision |
|---|
| 91 | real(kind=jprb), parameter :: IMinstdScale = 1.0_jprb / real(IMinstdM,jprb) |
|---|
| 92 | |
|---|
| 93 | !--------------------------------------------------------------------- |
|---|
| 94 | ! A random number generator type: after being initialized with a |
|---|
| 95 | ! seed, type and optionally a number of vector streams, subsequent |
|---|
| 96 | ! calls to "uniform_distribution" are used to fill 1D or 2D arrays |
|---|
| 97 | ! with random numbers in a way that ought to be fast. |
|---|
| 98 | type rng_type |
|---|
| 99 | |
|---|
| 100 | integer(kind=jpim) :: itype = IRngNative |
|---|
| 101 | RNG_STATE_TYPE :: istate(NMaxStreams) |
|---|
| 102 | integer(kind=jpim) :: nmaxstreams = NMaxStreams |
|---|
| 103 | integer(kind=jpim) :: iseed |
|---|
| 104 | |
|---|
| 105 | contains |
|---|
| 106 | procedure :: initialize |
|---|
| 107 | procedure :: uniform_distribution_1d, & |
|---|
| 108 | & uniform_distribution_2d, & |
|---|
| 109 | & uniform_distribution_2d_masked |
|---|
| 110 | generic :: uniform_distribution => uniform_distribution_1d, & |
|---|
| 111 | & uniform_distribution_2d, & |
|---|
| 112 | & uniform_distribution_2d_masked |
|---|
| 113 | |
|---|
| 114 | end type rng_type |
|---|
| 115 | |
|---|
| 116 | contains |
|---|
| 117 | |
|---|
| 118 | !--------------------------------------------------------------------- |
|---|
| 119 | ! Initialize a random number generator, where "itype" may be either |
|---|
| 120 | ! IRngNative, indicating to use Fortran's native random_number |
|---|
| 121 | ! subroutine, or IRngMinstdVector, indicating to use the MINSTD |
|---|
| 122 | ! linear congruential generator (LCG). In the latter case |
|---|
| 123 | ! "nmaxstreams" should be provided indicating that random numbers |
|---|
| 124 | ! will be requested in blocks of this length. The generator is |
|---|
| 125 | ! seeded with "iseed". |
|---|
| 126 | subroutine initialize(this, itype, iseed, nmaxstreams) |
|---|
| 127 | |
|---|
| 128 | class(rng_type), intent(inout) :: this |
|---|
| 129 | integer(kind=jpim), intent(in), optional :: itype |
|---|
| 130 | integer(kind=jpim), intent(in), optional :: iseed |
|---|
| 131 | integer(kind=jpim), intent(in), optional :: nmaxstreams |
|---|
| 132 | |
|---|
| 133 | integer, allocatable :: iseednative(:) |
|---|
| 134 | integer :: nseed, jseed, jstr |
|---|
| 135 | real(jprd) :: rseed ! Note this must be in double precision |
|---|
| 136 | |
|---|
| 137 | if (present(itype)) then |
|---|
| 138 | this%itype = itype |
|---|
| 139 | else |
|---|
| 140 | this%itype = IRngNative |
|---|
| 141 | end if |
|---|
| 142 | |
|---|
| 143 | if (present(iseed)) then |
|---|
| 144 | this%iseed = iseed |
|---|
| 145 | else |
|---|
| 146 | this%iseed = 1 |
|---|
| 147 | end if |
|---|
| 148 | |
|---|
| 149 | if (present(nmaxstreams)) then |
|---|
| 150 | this%nmaxstreams = nmaxstreams |
|---|
| 151 | else |
|---|
| 152 | this%nmaxstreams = NMaxStreams |
|---|
| 153 | end if |
|---|
| 154 | |
|---|
| 155 | if (this%itype == IRngMinstdVector) then |
|---|
| 156 | ! ! OPTION 1: Use the C++ minstd_rand0 algorithm to populate the |
|---|
| 157 | ! ! state: this loop is not vectorizable because the state in |
|---|
| 158 | ! ! one stream depends on the one in the previous stream. |
|---|
| 159 | ! this%istate(1) = this%iseed |
|---|
| 160 | ! do jseed = 2,this%nmaxstreams |
|---|
| 161 | ! this%istate(jseed) = mod(IMinstdA0 * this%istate(jseed-1), IMinstdM) |
|---|
| 162 | ! end do |
|---|
| 163 | |
|---|
| 164 | ! OPTION 2: Use a modified (and vectorized) C++ minstd_rand0 algorithm to |
|---|
| 165 | ! populate the state |
|---|
| 166 | rseed = real(abs(this%iseed),jprd) |
|---|
| 167 | do jstr = 1,this%nmaxstreams |
|---|
| 168 | ! Note that nint returns an integer of type jpib (8-byte) |
|---|
| 169 | ! which may be converted to double if that is the type of |
|---|
| 170 | ! istate |
|---|
| 171 | this%istate(jstr) = nint(mod(rseed*jstr*(1.0_jprd-0.05_jprd*jstr & |
|---|
| 172 | & +0.005_jprd*jstr**2)*IMinstdA0, real(IMinstdM,jprd)),kind=jpib) |
|---|
| 173 | end do |
|---|
| 174 | |
|---|
| 175 | ! One warmup of the C++ minstd_rand algorithm |
|---|
| 176 | do jstr = 1,this%nmaxstreams |
|---|
| 177 | this%istate(jstr) = mod(IMinstdA * this%istate(jstr), IMinstdM) |
|---|
| 178 | end do |
|---|
| 179 | |
|---|
| 180 | else |
|---|
| 181 | ! Native generator by default |
|---|
| 182 | call random_seed(size=nseed) |
|---|
| 183 | allocate(iseednative(nseed)) |
|---|
| 184 | do jseed = 1,nseed |
|---|
| 185 | iseednative(jseed) = this%iseed + jseed - 1 |
|---|
| 186 | end do |
|---|
| 187 | call random_seed(put=iseednative) |
|---|
| 188 | deallocate(iseednative) |
|---|
| 189 | end if |
|---|
| 190 | |
|---|
| 191 | end subroutine initialize |
|---|
| 192 | |
|---|
| 193 | !--------------------------------------------------------------------- |
|---|
| 194 | ! Populate vector "randnum" with pseudo-random numbers; if rannum is |
|---|
| 195 | ! of length greater than nmaxstreams (specified when the generator |
|---|
| 196 | ! was initialized) then only the first nmaxstreams elements will be |
|---|
| 197 | ! assigned. |
|---|
| 198 | subroutine uniform_distribution_1d(this, randnum) |
|---|
| 199 | |
|---|
| 200 | class(rng_type), intent(inout) :: this |
|---|
| 201 | real(kind=jprb), intent(out) :: randnum(:) |
|---|
| 202 | |
|---|
| 203 | integer :: imax, i |
|---|
| 204 | |
|---|
| 205 | if (this%itype == IRngMinstdVector) then |
|---|
| 206 | |
|---|
| 207 | imax = min(this%nmaxstreams, size(randnum)) |
|---|
| 208 | |
|---|
| 209 | ! C++ minstd_rand algorithm |
|---|
| 210 | do i = 1, imax |
|---|
| 211 | ! The following calculation is computed entirely with 8-byte |
|---|
| 212 | ! numbers (whether real or integer) |
|---|
| 213 | this%istate(i) = mod(IMinstdA * this%istate(i), IMinstdM) |
|---|
| 214 | ! Scale the current state to a number in working precision |
|---|
| 215 | ! (jprb) between 0 and 1 |
|---|
| 216 | randnum(i) = IMinstdScale * this%istate(i) |
|---|
| 217 | end do |
|---|
| 218 | |
|---|
| 219 | else |
|---|
| 220 | |
|---|
| 221 | call random_number(randnum) |
|---|
| 222 | |
|---|
| 223 | end if |
|---|
| 224 | |
|---|
| 225 | end subroutine uniform_distribution_1d |
|---|
| 226 | |
|---|
| 227 | |
|---|
| 228 | !--------------------------------------------------------------------- |
|---|
| 229 | ! Populate matrix "randnum" with pseudo-random numbers; if the inner |
|---|
| 230 | ! dimension of rannum is of length greater than nmaxstreams |
|---|
| 231 | ! (specified when the generator was initialized) then only the first |
|---|
| 232 | ! nmaxstreams elements along this dimension will be assigned. |
|---|
| 233 | subroutine uniform_distribution_2d(this, randnum) |
|---|
| 234 | |
|---|
| 235 | class(rng_type), intent(inout) :: this |
|---|
| 236 | real(kind=jprb), intent(out) :: randnum(:,:) |
|---|
| 237 | |
|---|
| 238 | integer :: imax, jblock, i |
|---|
| 239 | |
|---|
| 240 | if (this%itype == IRngMinstdVector) then |
|---|
| 241 | |
|---|
| 242 | imax = min(this%nmaxstreams, size(randnum,1)) |
|---|
| 243 | |
|---|
| 244 | ! C++ minstd_ran algorithm |
|---|
| 245 | do jblock = 1,size(randnum,2) |
|---|
| 246 | ! These lines should be vectorizable |
|---|
| 247 | do i = 1, imax |
|---|
| 248 | this%istate(i) = mod(IMinstdA * this%istate(i), IMinstdM) |
|---|
| 249 | randnum(i,jblock) = IMinstdScale * this%istate(i) |
|---|
| 250 | end do |
|---|
| 251 | end do |
|---|
| 252 | |
|---|
| 253 | else |
|---|
| 254 | |
|---|
| 255 | call random_number(randnum) |
|---|
| 256 | |
|---|
| 257 | end if |
|---|
| 258 | |
|---|
| 259 | end subroutine uniform_distribution_2d |
|---|
| 260 | |
|---|
| 261 | !--------------------------------------------------------------------- |
|---|
| 262 | ! Populate matrix "randnum" with pseudo-random numbers; if the inner |
|---|
| 263 | ! dimension of rannum is of length greater than nmaxstreams |
|---|
| 264 | ! (specified when the generator was initialized) then only the first |
|---|
| 265 | ! nmaxstreams elements along this dimension will be assigned. This |
|---|
| 266 | ! version only operates on outer dimensions for which "mask" is true. |
|---|
| 267 | subroutine uniform_distribution_2d_masked(this, randnum, mask) |
|---|
| 268 | |
|---|
| 269 | class(rng_type), intent(inout) :: this |
|---|
| 270 | real(kind=jprb), intent(inout) :: randnum(:,:) |
|---|
| 271 | logical, intent(in) :: mask(:) |
|---|
| 272 | |
|---|
| 273 | integer :: imax, jblock, i |
|---|
| 274 | |
|---|
| 275 | if (this%itype == IRngMinstdVector) then |
|---|
| 276 | |
|---|
| 277 | imax = min(this%nmaxstreams, size(randnum,1)) |
|---|
| 278 | |
|---|
| 279 | ! C++ minstd_ran algorithm |
|---|
| 280 | do jblock = 1,size(randnum,2) |
|---|
| 281 | if (mask(jblock)) then |
|---|
| 282 | ! These lines should be vectorizable |
|---|
| 283 | do i = 1, imax |
|---|
| 284 | this%istate(i) = mod(IMinstdA * this%istate(i), IMinstdM) |
|---|
| 285 | randnum(i,jblock) = IMinstdScale * this%istate(i) |
|---|
| 286 | end do |
|---|
| 287 | end if |
|---|
| 288 | end do |
|---|
| 289 | |
|---|
| 290 | else |
|---|
| 291 | |
|---|
| 292 | do jblock = 1,size(randnum,2) |
|---|
| 293 | if (mask(jblock)) then |
|---|
| 294 | call random_number(randnum(:,jblock)) |
|---|
| 295 | end if |
|---|
| 296 | end do |
|---|
| 297 | |
|---|
| 298 | end if |
|---|
| 299 | |
|---|
| 300 | end subroutine uniform_distribution_2d_masked |
|---|
| 301 | |
|---|
| 302 | |
|---|
| 303 | end module radiation_random_numbers |
|---|
| 304 | |
|---|