Octopus
quickrnd.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2024 M. Marques, A. Castro, A. Rubio, G. Bertsch, X. Andrade, A. Buccheri
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module quickrnd_oct_m
22 use, intrinsic :: iso_fortran_env
23 use debug_oct_m
24 use global_oct_m
27
28 implicit none
29
30 private :: &
31 xorshift64, &
33 public :: &
34 quickrnd, &
35 shiftseed, &
37
38 ! xorshift64 needs a good seed, in particular it cannot be seeded with 0.
39 ! These seeds are generated with the splitmix64 generator and have a
40 ! reasonable amount of entropy.
41 integer(int64), parameter, public :: &
42 splitmix64_123 = int(z'B4DC9BD462DE412B', int64), & ! splitmix64(123)
43 splitmix64_321 = int(z'E95F1C4EFCF85DEE', int64) ! splitmix64(321)
44
45 interface quickrnd
47 end interface quickrnd
48
49contains
50
55 function xorshift64(x) result(res)
56 integer(int64), intent(inout) :: x
57 integer(int64) :: res
58
59 assert(x /= 0_int64)
60
61 x = ieor(x, shiftl(x, 13))
62 x = ieor(x, shiftr(x, 7))
63 x = ieor(x, shiftl(x, 17))
64
65 res = x
66 end function xorshift64
67
76 function to_double(x) result(d)
77 integer(int64), intent(in) :: x
78 real(real64) :: d
79 real(real64), parameter :: two_pow_minus_53 = real(z'3CA0000000000000', real64)
80
81 d = real(shiftr(x, 11), real64) * two_pow_minus_53
82 end function to_double
83
84 ! ---------------------------------------------------------
85
86 subroutine dquickrnd_single(iseed, rnd)
87 integer(int64), intent(inout) :: iseed
88 real(real64), intent(out) :: rnd
89
90 ! no PUSH_SUB, called too often
91
92 rnd = to_double(xorshift64(iseed))
93
94 end subroutine dquickrnd_single
95
96 ! ---------------------------------------------------------
97
98 subroutine dquickrnd_array(iseed, nn, rnd)
99 integer(int64), intent(inout) :: iseed
100 integer, intent(in) :: nn
101 real(real64), intent(inout) :: rnd(:)
102
103 integer :: ii
104
105 push_sub(quickrnd_array)
106
107 do ii = 1, nn
108 rnd(ii) = to_double(xorshift64(iseed))
109 end do
110
111 pop_sub(quickrnd_array)
112
113 end subroutine dquickrnd_array
115 ! ---------------------------------------------------------
116
117 subroutine zquickrnd_array(iseed, nn, rnd)
118 integer(int64), intent(inout) :: iseed
119 integer, intent(in) :: nn
120 complex(real64), intent(inout) :: rnd(:)
121
122 real(real64), parameter :: inv_sqrt2 = 1.0_real64 / sqrt(m_two)
123 real(real64) :: re, im
124 integer :: ii
125
126 push_sub(quickrnd_array)
127
128 do ii = 1, nn
129 re = to_double(xorshift64(iseed))
130 im = to_double(xorshift64(iseed))
131 rnd(ii) = cmplx(re, im, real64) * inv_sqrt2
132 end do
133
134 pop_sub(quickrnd_array)
135
136 end subroutine zquickrnd_array
137
138 ! ---------------------------------------------------------
139
140 subroutine shiftseed(iseed, n)
141 integer(int64), intent(inout) :: iseed
142 integer(int64), intent(in) :: n
143
144 integer(int64) :: discard
145 integer(int64) :: ii
146
147 push_sub(shiftseed)
149 do ii = 1, n
150 discard = xorshift64(iseed)
151 end do
152
153 pop_sub(shiftseed)
154
155 end subroutine shiftseed
156
157
158 ! TODO(Alex) Issue #1074 . When MR !2523 is merged, scrap `xorshift_rnd_32` and `random_real_in_range`
159 ! in favour of functions introduced. Use them to implement `random_number_in_range`, as described in the
160 ! issue.
161 subroutine xorshift_rnd_32_biased(x)
162 integer(int32), intent(inout) :: x
163
164 x = ieor(x, ishft(x, 13_int32)) ! x ^= x << 13
165 x = ieor(x, ishft(x, -17_int32)) ! x ^= x >> 17
166 x = ieor(x, ishft(x, 5_int32)) ! x ^= x << 5
167
168 ! Ensure the highest bit is zero.
169 ! This keeps the result in the range [0, 2^31 - 1] (positive) BUT biases against zero
170 x = iand(x, z'7FFFFFFF')
171
172 end subroutine xorshift_rnd_32_biased
173
174
176 function random_real_in_range(x_min, x_max, seed) result(rand)
177 real(real64), intent(in ) :: x_min, x_max
178 integer(int32), intent(inout) :: seed
179 real(real64) :: rand
180
181 real(real64), parameter :: max_value_xorshift_32 = 2147483647_real64
182 real(real64) :: scale
183
184 call xorshift_rnd_32_biased(seed)
185 scale = real(seed, real64) / max_value_xorshift_32
186 rand = scale * (x_max - x_min) + x_min
187
188 end function random_real_in_range
189
190
209 subroutine reservoir_sampling_aexpj(weight, reservoir, seed_value)
210 real(real64), intent(in ) :: weight(:)
211 integer(int32), intent(out) :: reservoir(:)
212 integer(int32), intent(in ), optional :: seed_value
213
214 integer(int32) :: m
215 integer(int32) :: n
216 real(real64) :: u_i, u_i2
217 real(real64) :: Xw
218 real(real64), allocatable :: key(:)
219
220 integer(int32) :: min_key_index, i, seed_size
221 integer(int32), allocatable :: seed(:)
222 real(real64) :: thres_w
223
224 ! Random number range must be (0, 1), as log(0) would result in a key of -inf, and log(1) = 0
225 ! Set range to [a, b], chosen empirically (may not be optimal)
226 real(real64), parameter :: a = 1.e-10_real64
227 real(real64), parameter :: b = 1._real64 - 1.e-10_real64
228
230
231 m = size(reservoir)
232 n = size(weight)
233 ! Cannot ask for more numbers than are present in the population
234 assert(m <= n)
235
236 if (present(seed_value)) then
237 safe_allocate_source(seed(1:1), seed_value)
238 else
239 ! Fortran intrinsic returns an array of seeds
240 call random_seed(size=seed_size)
241 safe_allocate(seed(1:seed_size))
242 call random_seed(get=seed)
243 end if
244
245 ! Initialise the reservoir with the first m items of the population
246 min_key_index = 1
247 safe_allocate(key(1:m))
248
249 do i = 1, m
250 reservoir(i) = i
251 u_i = random_real_in_range(a, b, seed(1))
252 key(i) = u_i ** (1._real64 / weight(i))
253 if (key(i) < key(min_key_index)) min_key_index = i
254 enddo
255
256 u_i = random_real_in_range(a, b, seed(1))
257 xw = log(u_i) / log(key(min_key_index))
258
259 ! Perform swaps at jump-intervals, until the population is exhausted
260 i = m + 1
261 do while(i <= n)
262 xw = xw - weight(i)
263 if (xw <= 0._real64) then
264 reservoir(min_key_index) = i
265 thres_w = key(min_key_index)**weight(i)
266 ! U__{i2} \in (t_w, 1)
267 u_i2 = random_real_in_range(thres_w + a, b, seed(1))
268 key(min_key_index) = u_i2 ** (1._real64 / weight(i))
269 min_key_index = minloc(key, dim=1)
270 u_i = random_real_in_range(a, b, seed(1))
271 xw = log(u_i) / log(key(min_key_index))
272 endif
273 i = i + 1
274 enddo
275
276 safe_deallocate_a(key)
277 safe_deallocate_a(seed)
278
280
281 end subroutine reservoir_sampling_aexpj
282
283end module quickrnd_oct_m
284
285!! Local Variables:
286!! mode: f90
287!! coding: utf-8
288!! End:
double log(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_two
Definition: global.F90:189
subroutine dquickrnd_array(iseed, nn, rnd)
Definition: quickrnd.F90:192
integer(int64), parameter, public splitmix64_321
Definition: quickrnd.F90:134
subroutine, public reservoir_sampling_aexpj(weight, reservoir, seed_value)
Algorithm "A-ExpJ" weighted reservoir sampling without replacement.
Definition: quickrnd.F90:303
subroutine, public shiftseed(iseed, n)
Definition: quickrnd.F90:234
real(real64) function random_real_in_range(x_min, x_max, seed)
Generate a random real64 in the range .
Definition: quickrnd.F90:270
real(real64) function, private to_double(x)
Generating uniform doubles in the unit interval.
Definition: quickrnd.F90:170
subroutine xorshift_rnd_32_biased(x)
Definition: quickrnd.F90:255
subroutine dquickrnd_single(iseed, rnd)
Definition: quickrnd.F90:180
subroutine zquickrnd_array(iseed, nn, rnd)
Definition: quickrnd.F90:211
integer(int64) function, private xorshift64(x)
xorshift64 pseudorandom number generator
Definition: quickrnd.F90:149