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

Return the natural logarithm of the Probability Mass Function (PMF) of the Cyclic Geometric distribution for an input stepSuccess within the discrete integer support of the distribution \([0, period]\). More...

Detailed Description

Return the natural logarithm of the Probability Mass Function (PMF) of the Cyclic Geometric distribution for an input stepSuccess within the discrete integer support of the distribution \([0, period]\).

Parameters
[out]logPMF: The output scalar (or array of the same rank, shape, and size as other array-like arguments), of the same type and kind as logProbSuccess, containing the natural logarithm of the PMF of the distribution at the specified stepSuccess.
[in]stepSuccess: The input positive scalar (or array of the same rank, shape, and size as other array-like arguments), of
  1. type real of the same kind as that of logProbSuccess,
  2. type integer of default kind IK,
containing the values at which the PMF must be computed.
Note that stepSuccess represents the Bernoulli trial step at which the first success occurs.
[in]logProbSuccess: 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),
containing the natural logarithm of the probability of success at each Bernoulli step.
[in]logProbFailure: The input scalar (or array of the same rank, shape, and size as other array-like arguments), of the same type and kind as logProbSuccess containing the natural logarithm of the probability of failure at each Bernoulli step.
(optional, default = log(get1mexp(logProbSuccess)))


Possible calling interfaces

call setGeomCyclicLogPMF(logPMF, stepSuccess, logProbSuccess, period) ! elemental
call setGeomCyclicLogPMF(logPMF, stepSuccess, logProbSuccess, logProbFailure, period) ! elemental
call setGeomCyclicLogPMF(logPMF(:), stepSuccess(:), logProbSuccess, period)
call setGeomCyclicLogPMF(logPMF(:), stepSuccess(:), logProbSuccess, logProbFailure, period)
Return the natural logarithm of the Probability Mass Function (PMF) of the Cyclic Geometric distribut...
This module contains classes and procedures for computing various statistical quantities related to t...
Warning
The condition all(0 < [stepSuccess]) must hold for the corresponding input arguments.
The condition all([stepSuccess] <= period) must hold for the corresponding input arguments.
The condition logProbSuccess <= 0. must hold for the corresponding input arguments.
The condition logProbFailure < 0. must hold for the corresponding input arguments.
The condition exp(logProbFailure) + exp(logProbSuccess) == 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.
The pure procedure(s) documented herein become impure when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1.
By default, these procedures are pure in release build and impure in debug and testing builds.
Remarks
The procedures under discussion are elemental.
See also
getGeomCyclicLogCDF
setGeomCyclicLogCDF
getGeomCyclicLogPMF
setGeomCyclicLogPMF


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_kind, only: RKG => RK ! all processor kinds are supported.
5 use pm_io, only: display_type
7
8 implicit none
9
10 real(RKG) :: logPMF(3)
11 type(display_type) :: disp
12 disp = display_type(file = "main.out.F90")
13
14 call disp%skip()
15 call disp%show("call setGeomCyclicLogPMF(logPMF(1), 1_IK, logProbSuccess = log(.2_RKG), period = 1_IK)")
16 call setGeomCyclicLogPMF(logPMF(1), 1_IK, logProbSuccess = log(.2_RKG), period = 1_IK)
17 call disp%show("logPMF(1)")
18 call disp%show( logPMF(1) )
19 call disp%skip()
20
21 call disp%skip()
22 call disp%show("call setGeomCyclicLogPMF(logPMF(1), 1_IK, logProbSuccess = log(.2_RKG), period = 2_IK)")
23 call setGeomCyclicLogPMF(logPMF(1), 1_IK, logProbSuccess = log(.2_RKG), period = 2_IK)
24 call disp%show("logPMF(1)")
25 call disp%show( logPMF(1) )
26 call disp%skip()
27
28 call disp%skip()
29 call disp%show("call setGeomCyclicLogPMF(logPMF(1), 1_IK, logProbSuccess = log(.2_RKG), period = 20_IK)")
30 call setGeomCyclicLogPMF(logPMF(1), 1_IK, logProbSuccess = log(.2_RKG), period = 20_IK)
31 call disp%show("logPMF(1)")
32 call disp%show( logPMF(1) )
33 call disp%skip()
34
35 call disp%skip()
36 call disp%show("call setGeomCyclicLogPMF(logPMF(1), 2_IK, logProbSuccess = log(.2_RKG), period = 2_IK)")
37 call setGeomCyclicLogPMF(logPMF(1), 2_IK, logProbSuccess = log(.2_RKG), period = 2_IK)
38 call disp%show("logPMF(1)")
39 call disp%show( logPMF(1) )
40 call disp%skip()
41
42 call disp%skip()
43 call disp%show("call setGeomCyclicLogPMF(logPMF(1:3), [integer(IK) :: 1, 2, 3], logProbSuccess = log(.2_RKG), logProbFailure = log(1._RKG - .2_RKG), period = 4_IK)")
44 call setGeomCyclicLogPMF(logPMF(1:3), [integer(IK) :: 1, 2, 3], logProbSuccess = log(.2_RKG), logProbFailure = log(1._RKG - .2_RKG), period = 4_IK)
45 call disp%show("logPMF(1:3)")
46 call disp%show( logPMF(1:3) )
47 call disp%skip()
48
49 call disp%skip()
50 call disp%show("call setGeomCyclicLogPMF(logPMF(1:3), [integer(IK) :: 1, 2, 3], logProbSuccess = log([0.1_RKG, .2_RKG, 1._RKG]), period = 5_IK)")
51 call setGeomCyclicLogPMF(logPMF(1:3), [integer(IK) :: 1, 2, 3], logProbSuccess = log([0.1_RKG, .2_RKG, 1._RKG]), period = 5_IK)
52 call disp%show("logPMF(1:3)")
53 call disp%show( logPMF(1:3) )
54 call disp%skip()
55
56 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
57 ! Output an example array for visualization.
58 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
59
60 block
61
62 use pm_arrayRange, only: getRange
63 integer(IK), allocatable :: stepSuccess(:)
64 real(RKG) :: logPMF(4)
65 integer :: fileUnit, i
66
67 stepSuccess = getRange(1_IK, 10_IK)
68 open(newunit = fileUnit, file = "setGeomCyclicLogPMF.IK.txt")
69 do i = 1, size(stepSuccess)
70 call setGeomCyclicLogPMF(logPMF, stepSuccess(i), log([.05_RKG, .25_RKG, .05_RKG, .25_RKG]), period = [10_IK, 10_IK, 10000_IK, 10000_IK])
71 write(fileUnit, "(*(g0,:,','))") stepSuccess(i), exp(logPMF)
72 end do
73 close(fileUnit)
74
75 end block
76
77end program example
Generate minimally-spaced character, integer, real sequences or sequences at fixed intervals of size ...
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 ranges of discrete character,...
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 RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
Definition: pm_kind.F90:543
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
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
2call setGeomCyclicLogPMF(logPMF(1), 1_IK, logProbSuccess = log(.2_RKG), period = 1_IK)
3logPMF(1)
4+0.22204460492503131E-15
5
6
7call setGeomCyclicLogPMF(logPMF(1), 1_IK, logProbSuccess = log(.2_RKG), period = 2_IK)
8logPMF(1)
9-0.58778666490211884
10
11
12call setGeomCyclicLogPMF(logPMF(1), 1_IK, logProbSuccess = log(.2_RKG), period = 20_IK)
13logPMF(1)
14-1.5978417206981419
15
16
17call setGeomCyclicLogPMF(logPMF(1), 2_IK, logProbSuccess = log(.2_RKG), period = 2_IK)
18logPMF(1)
19-0.81093021621632855
20
21
22call setGeomCyclicLogPMF(logPMF(1:3), [integer(IK) :: 1, 2, 3], logProbSuccess = log(.2_RKG), logProbFailure = log(1._RKG - .2_RKG), period = 4_IK)
23logPMF(1:3)
24-1.0824829067382260, -1.3056264580524357, -1.5287700093666454
25
26
27call setGeomCyclicLogPMF(logPMF(1:3), [integer(IK) :: 1, 2, 3], logProbSuccess = log([0.1_RKG, .2_RKG, 1._RKG]), period = 5_IK)
28logPMF(1:3)
29-1.4097911370312888, -1.4355606024228085, -709.78271289338397
30
31

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
9fontsize = 17
10
11kind = "IK"
12label = [ r"probSuccess = .05, period = 10"
13 , r"probSuccess = .25, period = 10"
14 , r"probSuccess = .05, period = 10000"
15 , r"probSuccess = .25, period = 10000"
16 ]
17
18pattern = "*." + kind + ".txt"
19fileList = glob.glob(pattern)
20if len(fileList) == 1:
21
22 df = pd.read_csv(fileList[0], delimiter = ",", header = None)
23
24 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
25 ax = plt.subplot()
26
27 for i in range(1,len(df.values[0,:]+1)):
28
29 plt.plot( df.values[:, 0]
30 , df.values[:,i]
31 , marker = "."
32 )
33
34 plt.xticks(fontsize = fontsize - 2)
35 plt.yticks(fontsize = fontsize - 2)
36 ax.set_xlim([0, None])
37 ax.set_xlabel("stepSuccess", fontsize = fontsize)
38 ax.set_ylabel("PMF", fontsize = fontsize)
39
40 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
41 ax.tick_params(axis = "y", which = "minor")
42 ax.tick_params(axis = "x", which = "minor")
43
44 ax.legend ( label
45 , fontsize = fontsize
46 #, loc = "center left"
47 #, bbox_to_anchor = (1, 0.5)
48 )
49
50 plt.savefig(fileList[0].replace(".txt",".png"))
51
52else:
53
54 sys.exit("Ambiguous file list exists.")

Visualization of the example output
Test:
test_pm_distGeomCyclic
Todo:
Normal Priority: This generic interface can be extended to complex arguments.


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, Monday March 6, 2017, 3:22 pm, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin.

Definition at line 413 of file pm_distGeomCyclic.F90.


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