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

This is the derived type for declaring and generating objects of type xoshiro256ssg_type containing a unique instance of a greedy Xoshiro256** random number generator (RNG). More...

Inheritance diagram for pm_distUnif::xoshiro256ssg_type:
Collaboration diagram for pm_distUnif::xoshiro256ssg_type:

Public Attributes

integer(IK) pos = 0_IK
 The scalar of type integer of default kind IK, containing the position of the first unused bit of the stream component of the RNG object. More...
 
- Public Attributes inherited from pm_distUnif::xoshiro256ss_type
integer(IK64), dimension(xoshiro256ssStateSizestate = [ -5952639272145821898_IK64 , -2790978430781836137_IK64 , -4872796757339724681_IK64 , -6638731986642513151_IK64 ]
 The vector of size xoshiro256ssStateSize of type integer of kind IK64, containing the most recent RNG state.
More...
 
integer(IK64) stream
 The scalar of type integer of kind IK64, containing the most recently generated random 64-bit stream. More...
 

Detailed Description

This is the derived type for declaring and generating objects of type xoshiro256ssg_type containing a unique instance of a greedy Xoshiro256** random number generator (RNG).

Unlike the Xoshiro256** algorithm as implemented by the derived type xoshiro256ssw_type, the greedy version of the algorithm here does not waste any of the randomly generated 64 bits in each update the RNG state.
See also the documentation of xoshiro256ssg_typer for information on the constructor of this type.

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 deterministic 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).)
[in]jump: The input vector of size xoshiro256ssStateSize of type integer of kind IK64, whose value sets the jump size of the random number generator.
It can be,
  1. xoshiro256ssJump128, corresponding to a jump size of imageID * 2**128.
    This jump can be used to generate up to 2**128 unique RNG sequences in parallel, each with length 2**128.
  2. xoshiro256ssJump192, corresponding to a jump size of imageID * 2**192.
    This jump can be used to generate up to 2**64 unique RNG sequences in parallel, each with length 2**192.
(optional. default = xoshiro256ssJump128)
Returns
rng : The output scalar object (or array of objects, of the same rank and shape as the input array-like arguments) of type xoshiro256ssg_type containing an instance of a splitmix64 random number generator.


Possible calling interfaces

use pm_kind, only: IK
type(xoshiro256ssg_type) :: rng
rng = xoshiro256ssg_type(seed = seed, imageID = imageID, jump = jump)
rand = getUnifRand(rng, rand, lb, ub)
call setUnifRand(rng, rand, lb, ub)
call setUnifRand(rng, rand)
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...
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 xoshiro256ssg_type containing a...
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 users must not manipulate the component values at any stages of the random number generation.
Note
Without initializing objects of this derived type, the generated RNGs will always be deterministic, always yielding identical sequences.
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(xoshiro256ssg_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 = xoshiro256ssg_type()")
22 rng = xoshiro256ssg_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("call setUnifRand(rng, rand)")
47 call setUnifRand(rng, rand)
48 call disp%show("[-huge(rand) - 1_IKG, rand, huge(rand)]")
49 call disp%show( [-huge(rand) - 1_IKG, rand, huge(rand)] )
50 end do
51 call disp%skip()
52 end block
53
54 block
55 use pm_kind, only: IKG => IKL
56 integer(IKG) :: rand, lb, ub
57 do itry = 1, ntry
58 call disp%show("digits(0_IKG)")
59 call disp%show( digits(0_IKG) )
60 call disp%show("call setUnifRand(rng, rand)")
61 call setUnifRand(rng, rand)
62 call disp%show("[-huge(rand) - 1_IKG, rand, huge(rand)]")
63 call disp%show( [-huge(rand) - 1_IKG, rand, huge(rand)] )
64 end do
65 call disp%skip()
66 end block
67
68 block
69 use pm_kind, only: IKG => IKS
70 integer(IKG) :: rand, lb, ub
71 do itry = 1, ntry
72 call disp%show("digits(0_IKG)")
73 call disp%show( digits(0_IKG) )
74 call disp%show("lb = -3_IKG; ub = 5_IKG")
75 lb = -3_IKG; ub = 5_IKG
76 call disp%show("rand = getUnifRand(rng, lb, ub)")
77 rand = getUnifRand(rng, lb, ub)
78 call disp%show("[lb, rand, ub]")
79 call disp%show( [lb, rand, ub] )
80 call disp%show(.and."call setAsserted(lb <= rand rand <= ub)")
81 call setAsserted(lb <= rand .and. rand <= ub)
82 end do
83 call disp%skip()
84 end block
85
86 block
87 use pm_kind, only: IKG => IKL
88 integer(IKG) :: rand, lb, ub
89 do itry = 1, ntry
90 call disp%show("digits(0_IKG)")
91 call disp%show( digits(0_IKG) )
92 call disp%show("lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG")
93 lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG
94 call disp%show("rand = getUnifRand(rng, lb, ub)")
95 rand = getUnifRand(rng, lb, ub)
96 call disp%show("[lb, rand, ub]")
97 call disp%show( [lb, rand, ub] )
98 call disp%show(.and."call setAsserted(lb <= rand rand <= ub)")
99 call setAsserted(lb <= rand .and. rand <= ub)
100 end do
101 call disp%skip()
102 end block
103
104 block
105 character(2) :: rand, lb, ub
106 do itry = 1, ntry
107 call disp%show("lb = 'ai'; ub = 'by'")
108 lb = 'ai'; ub = 'by'
109 call disp%show("rand = getUnifRand(rng, lb, ub)")
110 rand = getUnifRand(rng, lb, ub)
111 call disp%show("[lb, rand, ub]")
112 call disp%show( [lb, rand, ub] )
113 call disp%show(.and."call setAsserted(lb <= rand rand <= ub)")
114 call setAsserted(lb <= rand .and. rand <= ub)
115 end do
116 call disp%skip()
117 end block
118
119 block
120 use pm_logicalCompare, only: operator(<=)
121 logical :: rand, lb, ub
122 do itry = 1, 10
123 call disp%show("lb = .false.; ub = .true.")
124 lb = .false.; ub = .true.
125 call disp%show("rand = getUnifRand(rng, lb, ub)")
126 rand = getUnifRand(rng, lb, ub)
127 call disp%show("[lb, rand, ub]")
128 call disp%show( [lb, rand, ub] )
129 call disp%show(.and."call setAsserted(lb <= rand rand <= ub)")
130 call setAsserted(lb <= rand .and. rand <= ub)
131 end do
132 call disp%skip()
133 end block
134
135 block
136 use pm_complexCompareAll, only: operator(<=), operator(<)
137 complex :: rand, lb, ub
138 do itry = 1, ntry
139 call disp%show("lb = (-1., +1.); ub = (1., +2.)")
140 lb = (-1., +1.); ub = (1., +2.)
141 call disp%show("rand = getUnifRand(rng, lb, ub)")
142 rand = getUnifRand(rng, lb, ub)
143 call disp%show("[lb, rand, ub]")
144 call disp%show( [lb, rand, ub] )
145 call disp%show(.and."call setAsserted(lb <= rand rand < ub)")
146 call setAsserted(lb <= rand .and. rand < ub)
147 end do
148 call disp%skip()
149 end block
150
151 block
152 use pm_kind, only: RKG => RKH
153 real(RKG) :: rand, lb, ub
154 call disp%show("lb = 2._RKG; ub = lb + spacing(lb)")
155 lb = 2._RKG; ub = lb + spacing(lb)
156 do itry = 1, ntry !* 100000
157 call disp%show("call setUnifRand(rng, rand, lb, ub)")
158 call setUnifRand(rng, rand, lb, ub)
159 call disp%show("[lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK)")
160 call disp%show( [lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK) )
161 call disp%show(.and."call setAsserted(lb <= rand rand < ub)")
162 call setAsserted(lb <= rand .and. rand < ub)
163 end do
164 call disp%skip()
165 end block
166
167 block
168 real :: rand, lb, ub
169 do itry = 1, ntry !* 100000
170 call disp%show("lb = -3.; ub = 5.")
171 lb = -3.; ub = 5.
172 call disp%show("rand = getUnifRand(rng, lb, ub)")
173 rand = getUnifRand(rng, lb, ub)
174 call disp%show("[lb, rand, ub]")
175 call disp%show( [lb, rand, ub] )
176 call disp%show(.and."call setAsserted(lb <= rand rand < ub)")
177 call setAsserted(lb <= rand .and. rand < ub)
178 end do
179 call disp%skip()
180 end block
181
182 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
183 ! Output an example rand array for visualization.
184 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
185
186 block
187 integer :: rand(5000)
188 call setUnifRand(rng, rand, -2, 3)
189 if (0 /= getErrTableWrite(SK_"xoshiro256ssg_type.IK.txt", rand)) error stop "Table writing failed."
190 end block
191
192 block
193 complex :: rand(5000)
194 call setUnifRand(rng, rand, (-2., +2.), (3., 5.))
195 if (0 /= getErrTableWrite(SK_"xoshiro256ssg_type.CK.txt", rand)) error stop "Table writing failed."
196 end block
197
198 block
199 real :: rand(5000)
200 call setUnifRand(rng, rand)
201 if (0 /= getErrTableWrite(SK_"xoshiro256ssg_type.RK.txt", rand)) error stop "Table writing failed."
202 end block
203
204end program example
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
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, +2, +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, -1, +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, +5, +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, -2, +5
27call setAsserted(lb <= rand .and. rand <= ub)
28
29digits(0_IKG)
30+31
31call setUnifRand(rng, rand)
32[-huge(rand) - 1_IKG, rand, huge(rand)]
33-2147483648, +1251149422, +2147483647
34digits(0_IKG)
35+31
36call setUnifRand(rng, rand)
37[-huge(rand) - 1_IKG, rand, huge(rand)]
38-2147483648, -656359233, +2147483647
39digits(0_IKG)
40+31
41call setUnifRand(rng, rand)
42[-huge(rand) - 1_IKG, rand, huge(rand)]
43-2147483648, +1642178385, +2147483647
44digits(0_IKG)
45+31
46call setUnifRand(rng, rand)
47[-huge(rand) - 1_IKG, rand, huge(rand)]
48-2147483648, +1134633695, +2147483647
49digits(0_IKG)
50+31
51call setUnifRand(rng, rand)
52[-huge(rand) - 1_IKG, rand, huge(rand)]
53-2147483648, -59977553, +2147483647
54
55digits(0_IKG)
56+7
57call setUnifRand(rng, rand)
58[-huge(rand) - 1_IKG, rand, huge(rand)]
59-128, +7, +127
60digits(0_IKG)
61+7
62call setUnifRand(rng, rand)
63[-huge(rand) - 1_IKG, rand, huge(rand)]
64-128, +116, +127
65digits(0_IKG)
66+7
67call setUnifRand(rng, rand)
68[-huge(rand) - 1_IKG, rand, huge(rand)]
69-128, -29, +127
70digits(0_IKG)
71+7
72call setUnifRand(rng, rand)
73[-huge(rand) - 1_IKG, rand, huge(rand)]
74-128, -11, +127
75digits(0_IKG)
76+7
77call setUnifRand(rng, rand)
78[-huge(rand) - 1_IKG, rand, huge(rand)]
79-128, -83, +127
80
81digits(0_IKG)
82+31
83lb = -3_IKG; ub = 5_IKG
84rand = getUnifRand(rng, lb, ub)
85[lb, rand, ub]
86-3, +2, +5
87call setAsserted(lb <= rand .and. rand <= ub)
88digits(0_IKG)
89+31
90lb = -3_IKG; ub = 5_IKG
91rand = getUnifRand(rng, lb, ub)
92[lb, rand, ub]
93-3, +4, +5
94call setAsserted(lb <= rand .and. rand <= ub)
95digits(0_IKG)
96+31
97lb = -3_IKG; ub = 5_IKG
98rand = getUnifRand(rng, lb, ub)
99[lb, rand, ub]
100-3, -2, +5
101call setAsserted(lb <= rand .and. rand <= ub)
102digits(0_IKG)
103+31
104lb = -3_IKG; ub = 5_IKG
105rand = getUnifRand(rng, lb, ub)
106[lb, rand, ub]
107-3, +4, +5
108call setAsserted(lb <= rand .and. rand <= ub)
109digits(0_IKG)
110+31
111lb = -3_IKG; ub = 5_IKG
112rand = getUnifRand(rng, lb, ub)
113[lb, rand, ub]
114-3, +0, +5
115call setAsserted(lb <= rand .and. rand <= ub)
116
117digits(0_IKG)
118+7
119lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG
120rand = getUnifRand(rng, lb, ub)
121[lb, rand, ub]
122-127, +57, +63
123call setAsserted(lb <= rand .and. rand <= ub)
124digits(0_IKG)
125+7
126lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG
127rand = getUnifRand(rng, lb, ub)
128[lb, rand, ub]
129-127, +52, +63
130call setAsserted(lb <= rand .and. rand <= ub)
131digits(0_IKG)
132+7
133lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG
134rand = getUnifRand(rng, lb, ub)
135[lb, rand, ub]
136-127, +40, +63
137call setAsserted(lb <= rand .and. rand <= ub)
138digits(0_IKG)
139+7
140lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG
141rand = getUnifRand(rng, lb, ub)
142[lb, rand, ub]
143-127, +38, +63
144call setAsserted(lb <= rand .and. rand <= ub)
145digits(0_IKG)
146+7
147lb = -huge(0_IKG); ub = huge(0_IKG) / 2_IKG
148rand = getUnifRand(rng, lb, ub)
149[lb, rand, ub]
150-127, +59, +63
151call setAsserted(lb <= rand .and. rand <= ub)
152
153lb = 'ai'; ub = 'by'
154rand = getUnifRand(rng, lb, ub)
155[lb, rand, ub]
156ai, bt, by
157call setAsserted(lb <= rand .and. rand <= ub)
158lb = 'ai'; ub = 'by'
159rand = getUnifRand(rng, lb, ub)
160[lb, rand, ub]
161ai, bj, by
162call setAsserted(lb <= rand .and. rand <= ub)
163lb = 'ai'; ub = 'by'
164rand = getUnifRand(rng, lb, ub)
165[lb, rand, ub]
166ai, aq, by
167call setAsserted(lb <= rand .and. rand <= ub)
168lb = 'ai'; ub = 'by'
169rand = getUnifRand(rng, lb, ub)
170[lb, rand, ub]
171ai, bm, by
172call setAsserted(lb <= rand .and. rand <= ub)
173lb = 'ai'; ub = 'by'
174rand = getUnifRand(rng, lb, ub)
175[lb, rand, ub]
176ai, bj, by
177call setAsserted(lb <= rand .and. rand <= ub)
178
179lb = .false.; ub = .true.
180rand = getUnifRand(rng, lb, ub)
181[lb, rand, ub]
182F, T, T
183call setAsserted(lb <= rand .and. rand <= ub)
184lb = .false.; ub = .true.
185rand = getUnifRand(rng, lb, ub)
186[lb, rand, ub]
187F, T, T
188call setAsserted(lb <= rand .and. rand <= ub)
189lb = .false.; ub = .true.
190rand = getUnifRand(rng, lb, ub)
191[lb, rand, ub]
192F, T, T
193call setAsserted(lb <= rand .and. rand <= ub)
194lb = .false.; ub = .true.
195rand = getUnifRand(rng, lb, ub)
196[lb, rand, ub]
197F, F, T
198call setAsserted(lb <= rand .and. rand <= ub)
199lb = .false.; ub = .true.
200rand = getUnifRand(rng, lb, ub)
201[lb, rand, ub]
202F, F, T
203call setAsserted(lb <= rand .and. rand <= ub)
204lb = .false.; ub = .true.
205rand = getUnifRand(rng, lb, ub)
206[lb, rand, ub]
207F, F, T
208call setAsserted(lb <= rand .and. rand <= ub)
209lb = .false.; ub = .true.
210rand = getUnifRand(rng, lb, ub)
211[lb, rand, ub]
212F, T, T
213call setAsserted(lb <= rand .and. rand <= ub)
214lb = .false.; ub = .true.
215rand = getUnifRand(rng, lb, ub)
216[lb, rand, ub]
217F, T, T
218call setAsserted(lb <= rand .and. rand <= ub)
219lb = .false.; ub = .true.
220rand = getUnifRand(rng, lb, ub)
221[lb, rand, ub]
222F, F, T
223call setAsserted(lb <= rand .and. rand <= ub)
224lb = .false.; ub = .true.
225rand = getUnifRand(rng, lb, ub)
226[lb, rand, ub]
227F, T, T
228call setAsserted(lb <= rand .and. rand <= ub)
229
230lb = (-1., +1.); ub = (1., +2.)
231rand = getUnifRand(rng, lb, ub)
232[lb, rand, ub]
233(-1.00000000, +1.00000000), (+0.366286516, +1.05448020), (+1.00000000, +2.00000000)
234call setAsserted(lb <= rand .and. rand < ub)
235lb = (-1., +1.); ub = (1., +2.)
236rand = getUnifRand(rng, lb, ub)
237[lb, rand, ub]
238(-1.00000000, +1.00000000), (-0.341614723, +1.39216471), (+1.00000000, +2.00000000)
239call setAsserted(lb <= rand .and. rand < ub)
240lb = (-1., +1.); ub = (1., +2.)
241rand = getUnifRand(rng, lb, ub)
242[lb, rand, ub]
243(-1.00000000, +1.00000000), (+0.316033840, +1.70363772), (+1.00000000, +2.00000000)
244call setAsserted(lb <= rand .and. rand < ub)
245lb = (-1., +1.); ub = (1., +2.)
246rand = getUnifRand(rng, lb, ub)
247[lb, rand, ub]
248(-1.00000000, +1.00000000), (+0.414571047, +1.95823407), (+1.00000000, +2.00000000)
249call setAsserted(lb <= rand .and. rand < ub)
250lb = (-1., +1.); ub = (1., +2.)
251rand = getUnifRand(rng, lb, ub)
252[lb, rand, ub]
253(-1.00000000, +1.00000000), (+0.594375849, +1.44850373), (+1.00000000, +2.00000000)
254call setAsserted(lb <= rand .and. rand < ub)
255
256lb = 2._RKG; ub = lb + spacing(lb)
257call setUnifRand(rng, rand, lb, ub)
258[lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK)
259 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000004
260call setAsserted(lb <= rand .and. rand < ub)
261call setUnifRand(rng, rand, lb, ub)
262[lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK)
263 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000004
264call setAsserted(lb <= rand .and. rand < ub)
265call setUnifRand(rng, rand, lb, ub)
266[lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK)
267 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000004
268call setAsserted(lb <= rand .and. rand < ub)
269call setUnifRand(rng, rand, lb, ub)
270[lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK)
271 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000004
272call setAsserted(lb <= rand .and. rand < ub)
273call setUnifRand(rng, rand, lb, ub)
274[lb, rand, ub], format = getFormat(width = 42_IK, ndigit = 35_IK)
275 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000000 , 2.0000000000000000000000000000000004
276call setAsserted(lb <= rand .and. rand < ub)
277
278lb = -3.; ub = 5.
279rand = getUnifRand(rng, lb, ub)
280[lb, rand, ub]
281-3.00000000, -2.09693289, +5.00000000
282call setAsserted(lb <= rand .and. rand < ub)
283lb = -3.; ub = 5.
284rand = getUnifRand(rng, lb, ub)
285[lb, rand, ub]
286-3.00000000, +2.65852451, +5.00000000
287call setAsserted(lb <= rand .and. rand < ub)
288lb = -3.; ub = 5.
289rand = getUnifRand(rng, lb, ub)
290[lb, rand, ub]
291-3.00000000, -1.81541717, +5.00000000
292call setAsserted(lb <= rand .and. rand < ub)
293lb = -3.; ub = 5.
294rand = getUnifRand(rng, lb, ub)
295[lb, rand, ub]
296-3.00000000, +3.02706671, +5.00000000
297call setAsserted(lb <= rand .and. rand < ub)
298lb = -3.; ub = 5.
299rand = getUnifRand(rng, lb, ub)
300[lb, rand, ub]
301-3.00000000, +4.98884773, +5.00000000
302call setAsserted(lb <= rand .and. rand < ub)
303
304

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:
Fatemeh Bagheri, Wednesday 12:20 AM, October 13, 2021, Dallas, TX

Definition at line 3362 of file pm_distUnif.F90.

Member Data Documentation

◆ pos

integer(IK) pm_distUnif::xoshiro256ssg_type::pos = 0_IK

The scalar of type integer of default kind IK, containing the position of the first unused bit of the stream component of the RNG object.

Definition at line 3366 of file pm_distUnif.F90.


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