ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_distUnif::splitmix64_type Type Reference

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). More...

Inheritance diagram for pm_distUnif::splitmix64_type:
Collaboration diagram for pm_distUnif::splitmix64_type:

Public Attributes

integer(IK64) stream
 The scalar of type integer of kind IK64, containing the most recently generated random 64-bit stream. More...
 
integer(IK64) state = 324108011427370141_IK64
 The scalar of type integer of kind IK64, containing the most recent RNG state.
More...
 

Detailed Description

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,

  1. declare an object of type splitmix64_type and initialize the object via the type constructor splitmix64_typer (see below for the possible calling interfaces),
  2. pass the generated RNG instance to the desired random number generating routines,
    1. getUnifRand,
    2. 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,
  1. The Coarray image ID as returned by Fortran intrinsic this_image() within a global team of Coarray images.
  2. The MPI rank of the processor (plus one) as returned by the MPI library intrinsic mpi_comm_rank().
  3. The OpenMP thread number as returned by the OpenMP library intrinsic omp_get_thread_num().
  4. 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

use pm_kind, only: IK
type(splitmix64_type) :: rng
rng = splitmix64_type(seed = seed, imageID = imageID)
print *, rng%state ! The current RNG state.
print *, rng%stream ! The current 64bit integer random number.
!
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...
Definition: pm_kind.F90:268
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
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

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_err, only: setAsserted
5 use pm_io, only: display_type
6 use pm_io, only: getFormat
7 use pm_distUnif, only: getUnifRand
8 use pm_distUnif, only: setUnifRand
10 use pm_io, only: getErrTableWrite
11
12 implicit none
13
14 integer :: itry, ntry = 5
15 type(splitmix64_type) :: rng
16 type(display_type) :: disp
17
18 disp = display_type(file = "main.out.F90")
19
20 call disp%skip()
21 call disp%show("rng = splitmix64_type()")
22 rng = splitmix64_type()
23
24 block
25 use pm_kind, only: IKD
26 integer :: rand, lb, ub
27 do itry = 1, ntry
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)")
31 rand = getUnifRand(rng, lb, ub)
32 call disp%show("[lb, rand, ub]")
33 call disp%show( [lb, rand, ub] )
34 call disp%show(.and."call setAsserted(lb <= rand rand <= ub)")
35 call setAsserted(lb <= rand .and. rand <= ub)
36 end do
37 call disp%skip()
38 end block
39
40 block
41 use pm_kind, only: IKG => IKS
42 integer(IKG) :: rand, lb, ub
43 do itry = 1, ntry
44 call disp%show("digits(0_IKG)")
45 call disp%show( digits(0_IKG) )
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)")
49 rand = getUnifRand(rng, lb, ub)
50 call disp%show("[lb, rand, ub]")
51 call disp%show( [lb, rand, ub] )
52 call disp%show(.and."call setAsserted(lb <= rand rand <= ub)")
53 call setAsserted(lb <= rand .and. rand <= ub)
54 end do
55 call disp%skip()
56 end block
57
58 block
59 use pm_kind, only: IKG => IKL
60 integer(IKG) :: rand, lb, ub
61 do itry = 1, ntry
62 call disp%show("digits(0_IKG)")
63 call disp%show( digits(0_IKG) )
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)")
67 rand = getUnifRand(rng, lb, ub)
68 call disp%show("[lb, rand, ub]")
69 call disp%show( [lb, rand, ub] )
70 call disp%show(.and."call setAsserted(lb <= rand rand <= ub)")
71 call setAsserted(lb <= rand .and. rand <= ub)
72 end do
73 call disp%skip()
74 end block
75
76 block
77 character(2) :: rand, lb, ub
78 do itry = 1, ntry
79 call disp%show("lb = 'ai'; ub = 'by'")
80 lb = 'ai'; ub = 'by'
81 call disp%show("rand = getUnifRand(rng, lb, ub)")
82 rand = getUnifRand(rng, lb, ub)
83 call disp%show("[lb, rand, ub]")
84 call disp%show( [lb, rand, ub] )
85 call disp%show(.and."call setAsserted(lb <= rand rand <= ub)")
86 call setAsserted(lb <= rand .and. rand <= ub)
87 end do
88 call disp%skip()
89 end block
90
91 block
92 use pm_logicalCompare, only: operator(<=)
93 logical :: rand, lb, ub
94 do itry = 1, 10
95 call disp%show("lb = .false.; ub = .true.")
96 lb = .false.; ub = .true.
97 call disp%show("rand = getUnifRand(rng, lb, ub)")
98 rand = getUnifRand(rng, lb, ub)
99 call disp%show("[lb, rand, ub]")
100 call disp%show( [lb, rand, ub] )
101 call disp%show(.and."call setAsserted(lb <= rand rand <= ub)")
102 call setAsserted(lb <= rand .and. rand <= ub)
103 end do
104 call disp%skip()
105 end block
106
107 block
108 use pm_complexCompareAll, only: operator(<=), operator(<)
109 complex :: rand, lb, ub
110 do itry = 1, ntry
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)")
114 rand = getUnifRand(rng, lb, ub)
115 call disp%show("[lb, rand, ub]")
116 call disp%show( [lb, rand, ub] )
117 call disp%show(.and."call setAsserted(lb <= rand rand < ub)")
118 call setAsserted(lb <= rand .and. rand < ub)
119 end do
120 call disp%skip()
121 end block
122
123 block
124 use pm_kind, only: RKG => RKH
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)
128 do itry = 1, ntry
129 call disp%show("call setUnifRand(rng, rand, lb, ub)")
130 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)")
134 call setAsserted(lb <= rand .and. rand < ub)
135 end do
136 call disp%skip()
137 end block
138
139 block
140 real :: rand, lb, ub
141 do itry = 1, ntry
142 call disp%show("lb = -3.; ub = 5.")
143 lb = -3.; ub = 5.
144 call disp%show("rand = getUnifRand(rng, lb, ub)")
145 rand = getUnifRand(rng, lb, ub)
146 call disp%show("[lb, rand, ub]")
147 call disp%show( [lb, rand, ub] )
148 call disp%show(.and."call setAsserted(lb <= rand rand < ub)")
149 call setAsserted(lb <= rand .and. rand < ub)
150 end do
151 call disp%skip()
152 end block
153
154 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
155 ! Output an example rand array for visualization.
156 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
157
158 block
159 integer :: rand(5000)
160 call setUnifRand(rng, rand, -2, 3)
161 if (0 /= getErrTableWrite(SK_"splitmix64_type.IK.txt", rand)) error stop "Table writing failed."
162 end block
163
164 block
165 complex :: rand(5000)
166 call setUnifRand(rng, rand, (-2., +2.), (3., 5.))
167 if (0 /= getErrTableWrite(SK_"splitmix64_type.CK.txt", rand)) error stop "Table writing failed."
168 end block
169
170 block
171 real :: rand(5000)
172 call setUnifRand(rng, rand)
173 if (0 /= getErrTableWrite(SK_"splitmix64_type.RK.txt", rand)) error stop "Table writing failed."
174 end block
175
176end program example
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...
Definition: pm_err.F90:735
Generate and return an object of type stop_type with the user-specified input attributes.
Definition: pm_err.F90:1618
Generate and return the iostat code resulting from writing the input table of rank 1 or 2 to the spec...
Definition: pm_io.F90:5940
Generate and return a generic or type/kind-specific IO format with the requested specifications that ...
Definition: pm_io.F90:18485
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
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.
Definition: pm_err.F90:52
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter IKS
The single-precision integer kind in Fortran mode. On most platforms, this is a 32-bit integer kind.
Definition: pm_kind.F90:563
integer, parameter IKL
The scalar integer constant of intrinsic default kind, representing the lowest range integer kind typ...
Definition: pm_kind.F90:749
integer, parameter IKD
The double precision integer kind in Fortran mode. On most platforms, this is a 64-bit integer kind.
Definition: pm_kind.F90:564
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
integer, parameter RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
Definition: pm_kind.F90:858
This module contains procedures and generic interfaces for performing a variety of logical comparison...
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2rng = splitmix64_type()
3lb = -3_IKD; ub = 5_IKD
4rand = getUnifRand(rng, lb, ub)
5[lb, rand, ub]
6-3, -3, +5
7call setAsserted(lb <= rand .and. rand <= ub)
8lb = -3_IKD; ub = 5_IKD
9rand = getUnifRand(rng, lb, ub)
10[lb, rand, ub]
11-3, -1, +5
12call setAsserted(lb <= rand .and. rand <= ub)
13lb = -3_IKD; ub = 5_IKD
14rand = getUnifRand(rng, lb, ub)
15[lb, rand, ub]
16-3, +5, +5
17call setAsserted(lb <= rand .and. rand <= ub)
18lb = -3_IKD; ub = 5_IKD
19rand = getUnifRand(rng, lb, ub)
20[lb, rand, ub]
21-3, +4, +5
22call setAsserted(lb <= rand .and. rand <= ub)
23lb = -3_IKD; ub = 5_IKD
24rand = getUnifRand(rng, lb, ub)
25[lb, rand, ub]
26-3, +3, +5
27call setAsserted(lb <= rand .and. rand <= ub)
28
29digits(0_IKG)
30+31
31lb = -3_IKG; ub = 5_IKG
32rand = getUnifRand(rng, lb, ub)
33[lb, rand, ub]
34-3, +1, +5
35call setAsserted(lb <= rand .and. rand <= ub)
36digits(0_IKG)
37+31
38lb = -3_IKG; ub = 5_IKG
39rand = getUnifRand(rng, lb, ub)
40[lb, rand, ub]
41-3, +5, +5
42call setAsserted(lb <= rand .and. rand <= ub)
43digits(0_IKG)
44+31
45lb = -3_IKG; ub = 5_IKG
46rand = getUnifRand(rng, lb, ub)
47[lb, rand, ub]
48-3, +2, +5
49call setAsserted(lb <= rand .and. rand <= ub)
50digits(0_IKG)
51+31
52lb = -3_IKG; ub = 5_IKG
53rand = getUnifRand(rng, lb, ub)
54[lb, rand, ub]
55-3, +4, +5
56call setAsserted(lb <= rand .and. rand <= ub)
57digits(0_IKG)
58+31
59lb = -3_IKG; ub = 5_IKG
60rand = getUnifRand(rng, lb, ub)
61[lb, rand, ub]
62-3, -1, +5
63call setAsserted(lb <= rand .and. rand <= ub)
64
65digits(0_IKG)
66+7
67lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG
68rand = getUnifRand(rng, lb, ub)
69[lb, rand, ub]
70-127, +1, +63
71call setAsserted(lb <= rand .and. rand <= ub)
72digits(0_IKG)
73+7
74lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG
75rand = getUnifRand(rng, lb, ub)
76[lb, rand, ub]
77-127, +50, +63
78call setAsserted(lb <= rand .and. rand <= ub)
79digits(0_IKG)
80+7
81lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG
82rand = getUnifRand(rng, lb, ub)
83[lb, rand, ub]
84-127, +60, +63
85call setAsserted(lb <= rand .and. rand <= ub)
86digits(0_IKG)
87+7
88lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG
89rand = getUnifRand(rng, lb, ub)
90[lb, rand, ub]
91-127, +63, +63
92call setAsserted(lb <= rand .and. rand <= ub)
93digits(0_IKG)
94+7
95lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG
96rand = getUnifRand(rng, lb, ub)
97[lb, rand, ub]
98-127, +45, +63
99call setAsserted(lb <= rand .and. rand <= ub)
100
101lb = 'ai'; ub = 'by'
102rand = getUnifRand(rng, lb, ub)
103[lb, rand, ub]
104ai, bn, by
105call setAsserted(lb <= rand .and. rand <= ub)
106lb = 'ai'; ub = 'by'
107rand = getUnifRand(rng, lb, ub)
108[lb, rand, ub]
109ai, ai, by
110call setAsserted(lb <= rand .and. rand <= ub)
111lb = 'ai'; ub = 'by'
112rand = getUnifRand(rng, lb, ub)
113[lb, rand, ub]
114ai, bx, by
115call setAsserted(lb <= rand .and. rand <= ub)
116lb = 'ai'; ub = 'by'
117rand = getUnifRand(rng, lb, ub)
118[lb, rand, ub]
119ai, am, by
120call setAsserted(lb <= rand .and. rand <= ub)
121lb = 'ai'; ub = 'by'
122rand = getUnifRand(rng, lb, ub)
123[lb, rand, ub]
124ai, bq, by
125call setAsserted(lb <= rand .and. rand <= ub)
126
127lb = .false.; ub = .true.
128rand = getUnifRand(rng, lb, ub)
129[lb, rand, ub]
130F, F, T
131call setAsserted(lb <= rand .and. rand <= ub)
132lb = .false.; ub = .true.
133rand = getUnifRand(rng, lb, ub)
134[lb, rand, ub]
135F, T, T
136call setAsserted(lb <= rand .and. rand <= ub)
137lb = .false.; ub = .true.
138rand = getUnifRand(rng, lb, ub)
139[lb, rand, ub]
140F, F, T
141call setAsserted(lb <= rand .and. rand <= ub)
142lb = .false.; ub = .true.
143rand = getUnifRand(rng, lb, ub)
144[lb, rand, ub]
145F, T, T
146call setAsserted(lb <= rand .and. rand <= ub)
147lb = .false.; ub = .true.
148rand = getUnifRand(rng, lb, ub)
149[lb, rand, ub]
150F, T, T
151call setAsserted(lb <= rand .and. rand <= ub)
152lb = .false.; ub = .true.
153rand = getUnifRand(rng, lb, ub)
154[lb, rand, ub]
155F, T, T
156call setAsserted(lb <= rand .and. rand <= ub)
157lb = .false.; ub = .true.
158rand = getUnifRand(rng, lb, ub)
159[lb, rand, ub]
160F, T, T
161call setAsserted(lb <= rand .and. rand <= ub)
162lb = .false.; ub = .true.
163rand = getUnifRand(rng, lb, ub)
164[lb, rand, ub]
165F, F, T
166call setAsserted(lb <= rand .and. rand <= ub)
167lb = .false.; ub = .true.
168rand = getUnifRand(rng, lb, ub)
169[lb, rand, ub]
170F, T, T
171call setAsserted(lb <= rand .and. rand <= ub)
172lb = .false.; ub = .true.
173rand = getUnifRand(rng, lb, ub)
174[lb, rand, ub]
175F, F, T
176call setAsserted(lb <= rand .and. rand <= ub)
177
178lb = (-1., +1.); ub = (1., +2.)
179rand = getUnifRand(rng, lb, ub)
180[lb, rand, ub]
181(-1.00000000, +1.00000000), (-0.626914859, +1.19532084), (+1.00000000, +2.00000000)
182call setAsserted(lb <= rand .and. rand < ub)
183lb = (-1., +1.); ub = (1., +2.)
184rand = getUnifRand(rng, lb, ub)
185[lb, rand, ub]
186(-1.00000000, +1.00000000), (+0.613393784E-1, +1.88830113), (+1.00000000, +2.00000000)
187call setAsserted(lb <= rand .and. rand < ub)
188lb = (-1., +1.); ub = (1., +2.)
189rand = getUnifRand(rng, lb, ub)
190[lb, rand, ub]
191(-1.00000000, +1.00000000), (+0.405627251, +1.67114520), (+1.00000000, +2.00000000)
192call setAsserted(lb <= rand .and. rand < ub)
193lb = (-1., +1.); ub = (1., +2.)
194rand = getUnifRand(rng, lb, ub)
195[lb, rand, ub]
196(-1.00000000, +1.00000000), (-0.466741443, +1.79003263), (+1.00000000, +2.00000000)
197call setAsserted(lb <= rand .and. rand < ub)
198lb = (-1., +1.); ub = (1., +2.)
199rand = getUnifRand(rng, lb, ub)
200[lb, rand, ub]
201(-1.00000000, +1.00000000), (-0.753254175, +1.98474681), (+1.00000000, +2.00000000)
202call setAsserted(lb <= rand .and. rand < ub)
203
204lb = 2._RKG; ub = lb + spacing(lb)
205call setUnifRand(rng, rand, lb, ub)
206[lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK)
207 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000004
208call setAsserted(lb <= rand .and. rand < ub)
209call setUnifRand(rng, rand, lb, ub)
210[lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK)
211 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000004
212call setAsserted(lb <= rand .and. rand < ub)
213call setUnifRand(rng, rand, lb, ub)
214[lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK)
215 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000004
216call setAsserted(lb <= rand .and. rand < ub)
217call setUnifRand(rng, rand, lb, ub)
218[lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK)
219 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000004
220call setAsserted(lb <= rand .and. rand < ub)
221call setUnifRand(rng, rand, lb, ub)
222[lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK)
223 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000004
224call setAsserted(lb <= rand .and. rand < ub)
225
226lb = -3.; ub = 5.
227rand = getUnifRand(rng, lb, ub)
228[lb, rand, ub]
229-3.00000000, +2.10108042, +5.00000000
230call setAsserted(lb <= rand .and. rand < ub)
231lb = -3.; ub = 5.
232rand = getUnifRand(rng, lb, ub)
233[lb, rand, ub]
234-3.00000000, +1.36400688, +5.00000000
235call setAsserted(lb <= rand .and. rand < ub)
236lb = -3.; ub = 5.
237rand = getUnifRand(rng, lb, ub)
238[lb, rand, ub]
239-3.00000000, +1.66196489, +5.00000000
240call setAsserted(lb <= rand .and. rand < ub)
241lb = -3.; ub = 5.
242rand = getUnifRand(rng, lb, ub)
243[lb, rand, ub]
244-3.00000000, +2.49997282, +5.00000000
245call setAsserted(lb <= rand .and. rand < ub)
246lb = -3.; ub = 5.
247rand = getUnifRand(rng, lb, ub)
248[lb, rand, ub]
249-3.00000000, -0.377151489, +5.00000000
250call setAsserted(lb <= rand .and. rand < ub)
251
252

Postprocessing of the example output
1#!/usr/bin/env python
2
3import matplotlib.pyplot as plt
4import pandas as pd
5import numpy as np
6import glob
7import sys
8
9linewidth = 2
10fontsize = 17
11
12marker ={ "CK" : "-"
13 , "IK" : "."
14 , "RK" : "-"
15 }
16xlab = { "CK" : "Uniform Random Value ( real/imaginary components )"
17 , "IK" : "Uniform Random Value ( integer-valued )"
18 , "RK" : "Uniform Random Value ( real-valued )"
19 }
20#legends = [ r"$lb = 0., ub = 1.$"
21# , r"$lb = 0., ub = 1.$"
22# , r"$lb = 0., ub = 1.$"
23# ]
24
25for kind in ["IK", "CK", "RK"]:
26
27 pattern = "*." + kind + ".txt"
28 fileList = glob.glob(pattern)
29 if len(fileList) == 1:
30
31 df = pd.read_csv(fileList[0], delimiter = ",", header = None)
32
33 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
34 ax = plt.subplot()
35
36 for j in range(len(df.values[0,:])):
37 if kind == "CK":
38 plt.hist( df.values[:,j]
39 , histtype = "stepfilled"
40 , alpha = 0.5
41 , bins = 75
42 )
43 else:
44 plt.hist( df.values[:,j]
45 , histtype = "stepfilled"
46 , alpha = 0.5
47 , bins = 75
48 )
49 #ax.legend ( legends
50 # , fontsize = fontsize
51 # )
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)
57
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")
61
62 plt.savefig(fileList[0].replace(".txt",".png"))
63
64 elif len(fileList) > 1:
65
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.

  1. 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.
  2. 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.

Author:
Amir Shahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin

Definition at line 3593 of file pm_distUnif.F90.

Member Data Documentation

◆ state

integer(IK64) pm_distUnif::splitmix64_type::state = 324108011427370141_IK64

The scalar of type integer of kind IK64, containing the most recent RNG state.

Definition at line 3599 of file pm_distUnif.F90.

◆ stream

integer(IK64) pm_distUnif::splitmix64_type::stream

The scalar of type integer of kind IK64, containing the most recently generated random 64-bit stream.

Definition at line 3596 of file pm_distUnif.F90.


The documentation for this type was generated from the following file: