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

Generate and return an allocatable array of rank 1 containing the state vector of the Fortran default random number generator (RNG) or, optionally set the RNG state based on a reference input scalar seed, optionally distinctly on each processor. More...

Detailed Description

Generate and return an allocatable array of rank 1 containing the state vector of the Fortran default random number generator (RNG) or, optionally set the RNG state based on a reference input scalar seed, optionally distinctly on each processor.

The procedures of this generic interface are merely wrappers around the procedures of the generic interface setUnifRandState that additionally return the current RNG state.
If both seed and imageID input arguments are missing, the procedures simply return the current RNG state.
Otherwise, the behavior is identical to that of setUnifRandState.

Parameters
[in]seed: See the documentation of the corresponding input argument to setUnifRandState.
[in]imageID: See the documentation of the corresponding input argument to setUnifRandState.
Returns
unifRandState : The output allocatable vector of rank 1 whose length equals the size of the random state vector of the default Fortran RNG as returned by the Fortran intrinsic random_seed(size = size).


Possible calling interfaces

use pm_lkind, only: IK
integer(IK), allocatable :: unifRandState(:)
unifRandState = getUnifRandState(seed = seed, imageID = imageID)
Generate and return an allocatable array of rank 1 containing the state vector of the Fortran default...
This module contains classes and procedures for computing various statistical quantities related to t...
Remarks
The procedures under discussion are impure.
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_io, only: display_type
7
8 implicit none
9
10 integer(IK), allocatable :: unifRandState(:)
11
12 type(display_type) :: disp
13
14 disp = display_type(file = "main.out.F90")
15
16 call disp%skip()
17 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%")
18 call disp%show("! Get the current RNG seed.")
19 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%")
20 call disp%skip()
21
22 call disp%skip()
23 call disp%show("unifRandState = getUnifRandState()")
24 unifRandState = getUnifRandState()
25 call disp%show("unifRandState")
26 call disp%show( unifRandState )
27 call disp%skip()
28
29 call disp%skip()
30 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
31 call disp%show("! Set the current RNG seed based on the input scalar `seed` and return a copy of seed vector.")
32 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
33 call disp%skip()
34
35 call disp%skip()
36 call disp%show("unifRandState = getUnifRandState(seed = 12345_IK)")
37 unifRandState = getUnifRandState(seed = 12345_IK)
38 call disp%show("unifRandState")
39 call disp%show( unifRandState )
40 call disp%skip()
41
42 call disp%skip()
43 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
44 call disp%show("! Set the current RNG seed on each image distinctly based on the current date and time and and the image ID and return a copy of seed vector.")
45 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
46 call disp%skip()
47
48 call disp%skip()
49 call disp%show("unifRandState = getUnifRandState(imageID = getImageID())")
50 unifRandState = getUnifRandState(imageID = getImageID())
51 call disp%show("unifRandState")
52 call disp%show( unifRandState )
53 call disp%skip()
54
55 call disp%skip()
56 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
57 call disp%show("! Set the current RNG seed on each image distinctly based on the input scalar `seed` and and the image ID and return a copy of seed vector.")
58 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
59 call disp%skip()
60
61 call disp%skip()
62 call disp%show("unifRandState = getUnifRandState(seed = 12345_IK, imageID = getImageID())")
63 unifRandState = getUnifRandState(seed = 12345_IK, imageID = getImageID())
64 call disp%show("unifRandState")
65 call disp%show( unifRandState )
66 call disp%skip()
67
68end program example
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
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
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
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
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
This module contains procedures and generic interfaces for facilitating parallel computations or comp...
integer(IK) function getImageID()
Generate and return the ID of the current Coarray image / MPI process / OpenMP thread,...
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
2!%%%%%%%%%%%%%%%%%%%%%%%%%%
3! Get the current RNG seed.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7unifRandState = getUnifRandState()
8unifRandState
9-750981303, +2070083405, +1733008626, -1262172988, -902083335, -1985651255, -1658240627, +1738057454
10
11
12!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
13! Set the current RNG seed based on the input scalar `seed` and return a copy of seed vector.
14!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
15
16
17unifRandState = getUnifRandState(seed = 12345_IK)
18unifRandState
19-930397057, -1271367333, -839688912, +1099379960, +310102611, +1239277010, +1840798532, +493908163
20
21
22!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
23! Set the current RNG seed on each image distinctly based on the current date and time and and the image ID and return a copy of seed vector.
24!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
25
26
27unifRandState = getUnifRandState(imageID = getImageID())
28unifRandState
29-1093144448, -1671167646, +1116058510, +206243038, -1983460214, -474484444, +649781684, -1219745346
30
31
32!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
33! Set the current RNG seed on each image distinctly based on the input scalar `seed` and and the image ID and return a copy of seed vector.
34!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
35
36
37unifRandState = getUnifRandState(seed = 12345_IK, imageID = getImageID())
38unifRandState
39-930397057, -1271367333, -839688912, +1099379960, +310102611, +1239277010, +1840798532, +493908163
40
41
Test:
test_pm_distUnif


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 3902 of file pm_distUnif.F90.


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