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

Generate and return a scalar (or array of arbitrary rank of) random value(s) from the Cyclic Geometric distribution.
More...

Detailed Description

Generate and return a scalar (or array of arbitrary rank of) random value(s) from the Cyclic Geometric distribution.

See the documentation of pm_distGeomCyclic for more details.

Parameters
[in]probSuccess: The input scalar (or array of the same rank, shape, and size as other array-like arguments), of
  1. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
representing the probability of success parameter of the distribution.
[in]period: The input positive scalar (or array of the same rank, shape, and size as other array-like arguments), of
type integer of default kind IK, containing the period of the Cyclic Geometric distribution.
The period represents the maximum number of Bernoulli trials in the experiment.
Returns
rand : The output positive scalar (or array of the same rank, shape, and size as other array-like arguments), of the same type and kind as probSuccess, containing random value(s) from the distribution.


Possible calling interfaces

rand = getGeomCyclicRand(probSuccess, period)
Generate and return a scalar (or array of arbitrary rank of) random value(s) from the Cyclic Geometri...
This module contains classes and procedures for computing various statistical quantities related to t...
Warning
The condition 0 < period must hold for the corresponding input arguments.
The condition 0 < probSuccess .and. probSuccess <= 1 must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Remarks
The procedures under discussion are impure.
The procedures under discussion are elemental.
See also
setGeomCyclicRand


Example usage

1program example
2
3 use pm_kind, only: SK, IK
4 use pm_kind, only: RKG => RKS ! all real kinds are supported.
8 use pm_io, only: display_type
9
10 implicit none
11
12 integer(IK), parameter :: NP = 1000_IK
13 integer(IK) :: rand(NP)
14 real(RKG) :: probSuccess(NP)
15
16 type(display_type) :: disp
17 disp = display_type(file = "main.out.F90")
18
19 call setLinSpace(probSuccess, x1 = 0.001_RKG, x2 = 1._RKG)
20
21 call disp%skip()
22 call disp%show("probSuccess(1)")
23 call disp%show( probSuccess(1) )
24 call disp%show("rand(1) = getGeomCyclicRand(probSuccess(1), period = 10_IK)")
25 rand(1) = getGeomCyclicRand(probSuccess(1), period = 10_IK)
26 call disp%show("rand(1)")
27 call disp%show( rand(1) )
28 call disp%skip()
29
30 call disp%skip()
31 call disp%show("probSuccess(1)")
32 call disp%show( probSuccess(1) )
33 call disp%show("rand(1:2) = getGeomCyclicRand(probSuccess(1), period = 10_IK)")
34 rand(1:2) = getGeomCyclicRand(probSuccess(1), period = 10_IK)
35 call disp%show("rand(1:2)")
36 call disp%show( rand(1:2) )
37 call disp%skip()
38
39 call disp%skip()
40 call disp%show("probSuccess(1)")
41 call disp%show( probSuccess(1) )
42 call disp%show("rand(1:2) = getGeomCyclicRand(probSuccess(1), period = 100_IK)")
43 rand(1:2) = getGeomCyclicRand(probSuccess(1), period = 100_IK)
44 call disp%show("rand(1:2)")
45 call disp%show( rand(1:2) )
46 call disp%skip()
47
48 call disp%skip()
49 call disp%show("probSuccess(NP-2:NP)")
50 call disp%show( probSuccess(NP-2:NP) )
51 call disp%show("rand(NP-2:NP) = getGeomCyclicRand(probSuccess(NP-2:NP), period = 100_IK)")
52 rand(NP-2:NP) = getGeomCyclicRand(probSuccess(NP-2:NP), period = 100_IK)
53 call disp%show("rand(NP-2:NP)")
54 call disp%show( rand(NP-2:NP) )
55 call disp%skip()
56
57 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
58 ! Output an example rand array for visualization.
59 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
60
61 block
62 integer(IK) :: rand(4)
63 integer :: fileUnit, i
64 open(newunit = fileUnit, file = "getGeomCyclicRand.IK.txt")
65 do i = 1, 5000
66 rand = getGeomCyclicRand([.05_RKG, .25_RKG, .05_RKG, .25_RKG], period = [10_IK, 10_IK, 10000_IK, 10000_IK])
67 write(fileUnit, "(*(g0,:,','))") rand
68 end do
69 close(fileUnit)
70 end block
71
72end program example
Return the linSpace output argument with size(linSpace) elements of evenly-spaced values over the int...
Return the logSpace output argument with size(logSpace) elements of logarithmically-evenly-spaced val...
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 procedures and generic interfaces for generating arrays with linear or logarithm...
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 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
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
Definition: pm_kind.F90:567
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
2probSuccess(1)
3+0.100000005E-2
4rand(1) = getGeomCyclicRand(probSuccess(1), period = 10_IK)
5rand(1)
6+6
7
8
9probSuccess(1)
10+0.100000005E-2
11rand(1:2) = getGeomCyclicRand(probSuccess(1), period = 10_IK)
12rand(1:2)
13+6, +6
14
15
16probSuccess(1)
17+0.100000005E-2
18rand(1:2) = getGeomCyclicRand(probSuccess(1), period = 100_IK)
19rand(1:2)
20+21, +21
21
22
23probSuccess(NP-2:NP)
24+0.998000026, +0.999000013, +1.00000000
25rand(NP-2:NP) = getGeomCyclicRand(probSuccess(NP-2:NP), period = 100_IK)
26rand(NP-2:NP)
27+1, +1, +1
28
29

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" : "Cyclic Geometric Random Value ( real/imaginary components )"
17 , "IK" : "Cyclic Geometric Random Value ( integer-valued )"
18 , "RK" : "Cyclic Geometric Random Value ( real-valued )"
19 }
20label = [ r"probSuccess = .05, period = 10"
21 , r"probSuccess = .25, period = 10"
22 , r"probSuccess = .05, period = 10000"
23 , r"probSuccess = .25, period = 10000"
24 ]
25
26for kind in ["IK", "CK", "RK"]:
27
28 pattern = "*." + kind + ".txt"
29 fileList = glob.glob(pattern)
30 if len(fileList) == 1:
31
32 df = pd.read_csv(fileList[0], delimiter = ",", header = None)
33
34 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
35 ax = plt.subplot()
36
37 for j in range(len(df.values[0,:])):
38 if kind == "CK":
39 plt.hist( df.values[:,j]
40 , histtype = "stepfilled"
41 , density = False
42 , alpha = 0.5
43 , bins = 75
44 )
45 else:
46 plt.hist( df.values[:,j]
47 , histtype = "stepfilled"
48 , density = False
49 , alpha = 0.5
50 , bins = 75
51 )
52 ax.legend ( label
53 , fontsize = fontsize
54 )
55 plt.xticks(fontsize = fontsize - 2)
56 plt.yticks(fontsize = fontsize - 2)
57 ax.set_xlim(0, 10)
58 ax.set_xlabel(xlab[kind], fontsize = 17)
59 ax.set_ylabel("Count", fontsize = 17)
60 ax.set_title("Histograms of {} Cyclic Geometric random values".format(len(df.values[:, 0])), fontsize = 17)
61
62 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
63 ax.tick_params(axis = "y", which = "minor")
64 ax.tick_params(axis = "x", which = "minor")
65
66 plt.savefig(fileList[0].replace(".txt",".png"))
67
68 elif len(fileList) > 1:
69
70 sys.exit("Ambiguous file list exists.")

Visualization of the example output
Test:
test_pm_distGeomCyclic


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, Oct 16, 2009, 11:14 AM, Michigan

Definition at line 1908 of file pm_distGeomCyclic.F90.


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