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

This is a concrete derived type whose instances can be used to define/request the default uniform random number generator (RNG) of the Fortran standard.
More...

Inheritance diagram for pm_distUnif::rngf_type:
Collaboration diagram for pm_distUnif::rngf_type:

Detailed Description

This is a concrete derived type whose instances can be used to define/request the default uniform random number generator (RNG) of the Fortran standard.

This type does not currently hold any components.
It is merely used to signal/request the use of the intrinsic Fortran RNG.
The Fortran programming language offers a single mechanism through the intrinsic procedure random_seed(size = size, put = put, get = get) to set the seed of the intrinsic uniform random number generator of Fortran random_number(harvest).
However, the internal implementation of the random number generator can vary from compiler to another.
Furthermore, the size and requirements for the initial state (seed) can be also different across compilers.
As such, Fortran programming language also offers a new intrinsic procedure random_init(repeatable, image_distinct).

  1. If repeatable = .true., the seed is set to the same processor-dependent value each time random_init(repeatable = .true., ...) is called from the same image.
    The program is, therefore, guaranteed to generate the same random number sequence after each call to random_init().
    The term same image means a single instance of program execution.
    The sequence of random numbers is different for repeated execution of the program.
    If it is .false., the seed is set to a processor-dependent value.
  2. If image_distinct = .true., the seed is set to a processor-dependent value that is distinct from th seed set by a call to random_init() in another (Coarray) image.
    If it is image_distinct = .false., the seed is set to a value that does not depend which image called random_init().

The Fortran intrinsic random_init() apparently does guarantee the generation of the same random sequence in multiple independent executions of the program.
This requires setting the RNG seed explicitly at the beginning of each program execution.

The rngf_type aims to facilitate the above goal by offering a single interface that combines random_init() and random_seed() intrinsic routines.

Parameters
[in]seed: The input scalar of type integer of default kind IK, containing a positive 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, even between multiple independent runs of the program compiled by the same compiler.
(optional. If missing, it is set to a value determined by 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 identically on all images.)
Returns
rng : The output scalar object of type rngf_type.


Possible calling interfaces

use pm_kind, only: IK
type(rngf_type) :: rng
rng = rngf_type(seed = seed, imageID = imageID)
print *, getUnifRandState()
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...
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 a concrete derived type whose instances can be used to define/request the default uniform ran...
Note
This generic interface does not aim to replicate the behavior of the Fortran intrinsic random_init() but rather to extend it.
If you wish to reset the RNG seed to a reproducible seed without specifying a scalar reference value for it via this generic interface, simply call random_init(repeatable = .true., ...).
However, unlike random_init(), this generic interface can guarantee reproducibility of RNG sequences within and between multiple program runs in serial or parallel Coarray/MPI/OpenMP, all compiled by the same compiler.
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
6 use pm_distUnif, only: getUnifRand
8 use pm_distUnif, only: rngf_type
9
10 implicit none
11
12 type(rngf_type) :: rng
13
14 type(display_type) :: disp
15
16 disp = display_type(file = "main.out.F90")
17
18 call disp%skip()
19 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%")
20 call disp%show("! Set a random RNG seed.")
21 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%")
22 call disp%skip()
23
24 call disp%skip()
25 call disp%show("rng = rngf_type()")
26 rng = rngf_type()
27 call disp%show("getUnifRandState()")
28 call disp%show( getUnifRandState() )
29 call disp%skip()
30
31 call disp%skip()
32 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
33 call disp%show("! Set the current RNG seed based on the input scalar `seed` and return a copy of seed vector.")
34 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
35 call disp%skip()
36
37 call disp%skip()
38 call disp%show("rng = rngf_type(seed = 12345_IK)")
39 rng = rngf_type(seed = 12345_IK)
40 call disp%show("getUnifRandState()")
41 call disp%show( getUnifRandState() )
42 call disp%skip()
43
44 call disp%skip()
45 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
46 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.")
47 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
48 call disp%skip()
49
50 call disp%skip()
51 call disp%show("rng = rngf_type(imageID = getImageID())")
52 rng = rngf_type(imageID = getImageID())
53 call disp%show("getUnifRandState()")
54 call disp%show( getUnifRandState() )
55 call disp%skip()
56
57 call disp%skip()
58 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
59 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.")
60 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
61 call disp%skip()
62
63 call disp%skip()
64 call disp%show("rng = rngf_type(seed = 12345_IK, imageID = getImageID())")
65 rng = rngf_type(seed = 12345_IK, imageID = getImageID())
66 call disp%show("getUnifRandState()")
67 call disp%show( getUnifRandState() )
68 call disp%skip()
69
70 call disp%skip()
71 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
72 call disp%show("! Reset the RNG seed deterministically on each image distinctly based on the input scalar `seed` and and the image ID and return a copy of seed vector.")
73 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
74 call disp%skip()
75
76 call disp%skip()
77 call disp%show("rng = rngf_type(seed = 12345_IK, imageID = getImageID())")
78 rng = rngf_type(seed = 12345_IK, imageID = getImageID())
79 call disp%show("getUnifRandState()")
80 call disp%show( getUnifRandState() )
81 call disp%show("getUnifRand(1, 100, 10)")
82 call disp%show( getUnifRand(1, 100, 10) )
83 call disp%show("rng = rngf_type(seed = 12345_IK, imageID = getImageID())")
84 rng = rngf_type(seed = 12345_IK, imageID = getImageID())
85 call disp%show("getUnifRandState()")
86 call disp%show( getUnifRandState() )
87 call disp%show("getUnifRand(1, 100, 10)")
88 call disp%show( getUnifRand(1, 100, 10) )
89 call disp%show("rng = rngf_type() ! random RNG seed reset.")
90 rng = rngf_type() ! random RNG seed reset.
91 call disp%show("getUnifRandState()")
92 call disp%show( getUnifRandState() )
93 call disp%show("getUnifRand(1, 100, 10)")
94 call disp%show( getUnifRand(1, 100, 10) )
95 call disp%skip()
96
97end program example
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
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
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
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! Set a random RNG seed.
4!%%%%%%%%%%%%%%%%%%%%%%%
5
6
7rng = rngf_type()
9-1068469044, +50763360, +772632643, +484198386, +1219917042, -134927478, -1658425175, +91529044
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
17rng = rngf_type(seed = 12345_IK)
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
27rng = rngf_type(imageID = getImageID())
29+498063613, -582571465, -1646054171, -1642142633, -2122310913, -396529976, -572129249, +1303191645
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
37rng = rngf_type(seed = 12345_IK, imageID = getImageID())
39-930397057, -1271367333, -839688912, +1099379960, +310102611, +1239277010, +1840798532, +493908163
40
41
42!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
43! Reset the RNG seed deterministically on each image distinctly based on the input scalar `seed` and and the image ID and return a copy of seed vector.
44!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
45
46
47rng = rngf_type(seed = 12345_IK, imageID = getImageID())
49-930397057, -1271367333, -839688912, +1099379960, +310102611, +1239277010, +1840798532, +493908163
50getUnifRand(1, 100, 10)
51+44, +74, +42, +71, +33, +82, +43, +18, +66, +10
52rng = rngf_type(seed = 12345_IK, imageID = getImageID())
54-930397057, -1271367333, -839688912, +1099379960, +310102611, +1239277010, +1840798532, +493908163
55getUnifRand(1, 100, 10)
56+44, +74, +42, +71, +33, +82, +43, +18, +66, +10
57rng = rngf_type() ! random RNG seed reset.
59-654935740, -740033451, +605112649, -1545028327, +47143573, +1897579915, -1898124548, -512840283
60getUnifRand(1, 100, 10)
61+41, +95, +43, +14, +5, +23, +11, +93, +85, +98
62
63


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


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