This is the derived type for declaring and generating objects of type splitmix64_type containing a unique instance of an splitmix64 random number generator (RNG).
See also the documentation of splitmix64_typer for information on the constructor of this type.
Splitmix64 is a pseudo-random number generator algorithm that originated from the Java programming language and is used in many other programming languages.
It is a fairly simple algorithm that is fast but unsuitable for cryptographic purposes and comparable tasks.
The Splitmix64 algorithm is frequently used to calculate random initial states (seed) of other more complex pseudo-random number generators.
The conventional splitmix64 holds one 64bit state (seed) variable and returns 64bits of random binary data upon subsequent calls.
Splitmix64 is comparatively fast; It requires only 9 64-bit arithmetic/logical operations per 64 bits of random binary stream generation.
A conventional linear RNG provides a generate method that returns one pseudorandom value and updates the state of the RNG.
But the splitable RNG also has a second operation, split, that replaces the original RNG with two (seemingly) independent RNG.
This is done by creating and returning a new such RNG and updating the state of the original object.
Applications
Splitable RNG make it easy to organize the use of pseudorandom numbers in multithreaded programs structured using forkjoin parallelism.
No locking or synchronization is required (other than the usual memory fence immediately after RNG creation).
Because the generate method has no loops or conditionals, it is also suitable for SIMD or GPU implementation.
Usage instructions This derived type contains only the most recently updated state and random bit stream of the RNG.
To generate random values of arbitrary intrinsic kinds (character
, integer
, logical
, complex
, real
) the user must,
-
declare an object of type splitmix64_type and initialize the object via the type constructor splitmix64_typer (see below for the possible calling interfaces),
-
pass the generated RNG instance to the desired random number generating routines,
-
getUnifRand,
-
setUnifRand,
to generate the desired scalar or sequence of random values.
- Parameters
-
[in] | seed | : The input scalar of type integer of kind IK64, containing an integer that serves as the starting point to generate the full RNG seed.
Specify this input argument if you wish to make random simulations reproducible and deterministic, even between multiple independent runs of the program compiled by the same compiler.
(optional. If missing, it is set to a processor-dependent value based on the current date and time.) |
[in] | imageID | : The input positive scalar integer of default kind IK containing the ID of the current image/thread/process.
This can be,
-
The Coarray image ID as returned by Fortran intrinsic
this_image() within a global team of Coarray images.
-
The MPI rank of the processor (plus one) as returned by the MPI library intrinsic
mpi_comm_rank() .
-
The OpenMP thread number as returned by the OpenMP library intrinsic
omp_get_thread_num() .
-
Any (positive) integer that uniquely identifies the current processor from other processes.
The image/process/thread ID can be readily obtained by calling getImageID.
This number will be used to set the RNG seed uniquely on each processor.
(optional. If missing, the RNG seed will be set as if it is called in a serial application (or called on the first image).) |
- Returns
rng
: The output scalar object (or array of objects, of the same rank and shape as the input array-like arguments) of type splitmix64_type containing an instance of a splitmix64 random number generator.
Possible calling interfaces ⛓
type(splitmix64_type) :: rng
print *, rng%state
print *, rng%stream
!
This module contains classes and procedures for computing various statistical quantities related to t...
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
This is the derived type for declaring and generating objects of type splitmix64_type containing a un...
- Warning
- The condition
0 < imageID
must hold for the corresponding input arguments.
This condition is verified only if the library is built with the preprocessor macro CHECK_ENABLED=1
.
-
Although the components of this derived type are
public
, they are theoretically protected
.
The end user must not change the values of the components of an object of type splitmix64_type at anytime during the random value generation.
-
The Splitmix64 algorithm is not necessarily the best option for generating random values, particularly, for parallel applications.
Use xoshiro256** algorithm instead.
- See also
- rngf
isHead
getUnifCDF
getUnifRand
setUnifRand
getUnifRandState
setUnifRandState
rngu_type
rngf_type
splitmix64_type
xoshiro256ssw_type
getUnifRandStateSize
Example usage ⛓
14 integer :: itry, ntry
= 5
15 type(splitmix64_type) :: rng
16 type(display_type) :: disp
21 call disp%show(
"rng = splitmix64_type()")
26 integer :: rand, lb, ub
28 call disp%show(
"lb = -3_IKD; ub = 5_IKD")
29 lb
= -3_IKD; ub
= 5_IKD
30 call disp%show(
"rand = getUnifRand(rng, lb, ub)")
34 call disp%show(
.and."call setAsserted(lb <= rand rand <= ub)")
42 integer(IKG) :: rand, lb, ub
46 call disp%show(
"lb = -3_IKG; ub = 5_IKG")
47 lb
= -3_IKG; ub
= 5_IKG
48 call disp%show(
"rand = getUnifRand(rng, lb, ub)")
52 call disp%show(
.and."call setAsserted(lb <= rand rand <= ub)")
60 integer(IKG) :: rand, lb, ub
64 call disp%show(
"lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG")
65 lb
= -huge(
0_IKG); ub
= huge(
0_IKG)
/ 2_IKG
66 call disp%show(
"rand = getUnifRand(rng, lb, ub)")
70 call disp%show(
.and."call setAsserted(lb <= rand rand <= ub)")
77 character(
2) :: rand, lb, ub
79 call disp%show(
"lb = 'ai'; ub = 'by'")
81 call disp%show(
"rand = getUnifRand(rng, lb, ub)")
85 call disp%show(
.and."call setAsserted(lb <= rand rand <= ub)")
93 logical :: rand, lb, ub
95 call disp%show(
"lb = .false.; ub = .true.")
96 lb
= .false.; ub
= .true.
97 call disp%show(
"rand = getUnifRand(rng, lb, ub)")
101 call disp%show(
.and."call setAsserted(lb <= rand rand <= ub)")
109 complex :: rand, lb, ub
111 call disp%show(
"lb = (-1., +1.); ub = (1., +2.)")
112 lb
= (
-1.,
+1.); ub
= (
1.,
+2.)
113 call disp%show(
"rand = getUnifRand(rng, lb, ub)")
117 call disp%show(
.and."call setAsserted(lb <= rand rand < ub)")
125 real(RKG) :: rand, lb, ub
126 call disp%show(
"lb = 2._RKG; ub = lb + spacing(lb)")
127 lb
= 2._RKG; ub
= lb
+ spacing(lb)
129 call disp%show(
"call setUnifRand(rng, rand, lb, ub)")
131 call disp%show(
"[lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK)")
132 call disp%show( [lb, rand, ub],
format = getFormat(width
= 42_IK, ndigit
= 35_IK) )
133 call disp%show(
.and."call setAsserted(lb <= rand rand < ub)")
144 call disp%show(
"rand = getUnifRand(rng, lb, ub)")
148 call disp%show(
.and."call setAsserted(lb <= rand rand < ub)")
159 integer :: rand(
5000)
161 if (
0 /= getErrTableWrite(SK_
"splitmix64_type.IK.txt", rand))
error stop "Table writing failed."
165 complex :: rand(
5000)
167 if (
0 /= getErrTableWrite(SK_
"splitmix64_type.CK.txt", rand))
error stop "Table writing failed."
173 if (
0 /= getErrTableWrite(SK_
"splitmix64_type.RK.txt", rand))
error stop "Table writing failed."
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
Return a uniform random scalar or contiguous array of arbitrary rank of randomly uniformly distribute...
Verify the input assertion holds and if it does not, print the (optional) input message on stdout and...
Generate and return an object of type stop_type with the user-specified input attributes.
Generate and return the iostat code resulting from writing the input table of rank 1 or 2 to the spec...
This is a generic method of the derived type display_type with pass attribute.
This module contains procedures and generic interfaces for checking if both of the corresponding real...
This module contains classes and procedures for reporting and handling errors.
This module contains classes and procedures for input/output (IO) or generic display operations on st...
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
integer, parameter IKS
The single-precision integer kind in Fortran mode. On most platforms, this is a 32-bit integer kind.
integer, parameter IKL
The scalar integer constant of intrinsic default kind, representing the lowest range integer kind typ...
integer, parameter IKD
The double precision integer kind in Fortran mode. On most platforms, this is a 64-bit integer kind.
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
integer, parameter RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
This module contains procedures and generic interfaces for performing a variety of logical comparison...
Generate and return an object of type display_type.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example output ⛓
13lb
= -3_IKD; ub
= 5_IKD
18lb
= -3_IKD; ub
= 5_IKD
23lb
= -3_IKD; ub
= 5_IKD
31lb
= -3_IKG; ub
= 5_IKG
38lb
= -3_IKG; ub
= 5_IKG
45lb
= -3_IKG; ub
= 5_IKG
52lb
= -3_IKG; ub
= 5_IKG
59lb
= -3_IKG; ub
= 5_IKG
67lb
= -huge(
0_IKG); ub
= huge(
0_IKG)
/ 2_IKG
74lb
= -huge(
0_IKG); ub
= huge(
0_IKG)
/ 2_IKG
81lb
= -huge(
0_IKG); ub
= huge(
0_IKG)
/ 2_IKG
88lb
= -huge(
0_IKG); ub
= huge(
0_IKG)
/ 2_IKG
95lb
= -huge(
0_IKG); ub
= huge(
0_IKG)
/ 2_IKG
127lb
= .false.; ub
= .true.
132lb
= .false.; ub
= .true.
137lb
= .false.; ub
= .true.
142lb
= .false.; ub
= .true.
147lb
= .false.; ub
= .true.
152lb
= .false.; ub
= .true.
157lb
= .false.; ub
= .true.
162lb
= .false.; ub
= .true.
167lb
= .false.; ub
= .true.
172lb
= .false.; ub
= .true.
178lb
= (
-1.,
+1.); ub
= (
1.,
+2.)
181(
-1.00000000,
+1.00000000), (
+0.391060114E-1,
+1.99536943), (
+1.00000000,
+2.00000000)
183lb
= (
-1.,
+1.); ub
= (
1.,
+2.)
186(
-1.00000000,
+1.00000000), (
-0.726765394E-1,
+1.40922379), (
+1.00000000,
+2.00000000)
188lb
= (
-1.,
+1.); ub
= (
1.,
+2.)
191(
-1.00000000,
+1.00000000), (
-0.234985471,
+1.64044619), (
+1.00000000,
+2.00000000)
193lb
= (
-1.,
+1.); ub
= (
1.,
+2.)
196(
-1.00000000,
+1.00000000), (
-0.747754216,
+1.98392427), (
+1.00000000,
+2.00000000)
198lb
= (
-1.,
+1.); ub
= (
1.,
+2.)
201(
-1.00000000,
+1.00000000), (
-0.245621204E-1,
+1.83962965), (
+1.00000000,
+2.00000000)
204lb
= 2._RKG; ub
= lb
+ spacing(lb)
206[lb, rand, ub],
format = getFormat(width
= 42_IK, ndigit
= 35_IK)
207 2.0000000000000000000000000000000000 ,
2.0000000000000000000000000000000000 ,
2.0000000000000000000000000000000004
210[lb, rand, ub],
format = getFormat(width
= 42_IK, ndigit
= 35_IK)
211 2.0000000000000000000000000000000000 ,
2.0000000000000000000000000000000000 ,
2.0000000000000000000000000000000004
214[lb, rand, ub],
format = getFormat(width
= 42_IK, ndigit
= 35_IK)
215 2.0000000000000000000000000000000000 ,
2.0000000000000000000000000000000000 ,
2.0000000000000000000000000000000004
218[lb, rand, ub],
format = getFormat(width
= 42_IK, ndigit
= 35_IK)
219 2.0000000000000000000000000000000000 ,
2.0000000000000000000000000000000000 ,
2.0000000000000000000000000000000004
222[lb, rand, ub],
format = getFormat(width
= 42_IK, ndigit
= 35_IK)
223 2.0000000000000000000000000000000000 ,
2.0000000000000000000000000000000000 ,
2.0000000000000000000000000000000004
229-3.00000000,
+2.26409197,
+5.00000000
234-3.00000000,
+0.109087944,
+5.00000000
239-3.00000000,
+4.85142183,
+5.00000000
244-3.00000000,
+4.53214359,
+5.00000000
249-3.00000000,
+3.14565611,
+5.00000000
Postprocessing of the example output ⛓
3import matplotlib.pyplot
as plt
16xlab = {
"CK" :
"Uniform Random Value ( real/imaginary components )"
17 ,
"IK" :
"Uniform Random Value ( integer-valued )"
18 ,
"RK" :
"Uniform Random Value ( real-valued )"
25for kind
in [
"IK",
"CK",
"RK"]:
27 pattern =
"*." + kind +
".txt"
28 fileList = glob.glob(pattern)
29 if len(fileList) == 1:
31 df = pd.read_csv(fileList[0], delimiter =
",", header =
None)
33 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
36 for j
in range(len(df.values[0,:])):
38 plt.hist( df.values[:,j]
39 , histtype =
"stepfilled"
44 plt.hist( df.values[:,j]
45 , histtype =
"stepfilled"
52 plt.xticks(fontsize = fontsize - 2)
53 plt.yticks(fontsize = fontsize - 2)
54 ax.set_xlabel(xlab[kind], fontsize = 17)
55 ax.set_ylabel(
"Count", fontsize = 17)
56 ax.set_title(
"Histograms of {} Uniform random values".format(len(df.values[:, 0])), fontsize = 17)
58 plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
59 ax.tick_params(axis =
"y", which =
"minor")
60 ax.tick_params(axis =
"x", which =
"minor")
62 plt.savefig(fileList[0].replace(
".txt",
".png"))
64 elif len(fileList) > 1:
66 sys.exit(
"Ambiguous file list exists.")
Visualization of the example output ⛓
- Test:
- test_pm_distUnif
- Todo:
- High Priority: An illustration of the distribution of the probability of individual bits being
0
or 1
in the mantissa of real
-valued random numbers and integer
random numbers must be added to the example.
Final Remarks ⛓
If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.
-
If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
-
If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.
This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.
- Copyright
- Computational Data Science Lab
- Author:
- Amir Shahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas Austin
Definition at line 3593 of file pm_distUnif.F90.