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

Generate and return the Cumulative Distribution Function (CDF) of the GenExpGamma distribution for an input x within the support of the distribution \(x \in (-\infty,+\infty)\). More...

Detailed Description

Generate and return the Cumulative Distribution Function (CDF) of the GenExpGamma distribution for an input x within the support of the distribution \(x \in (-\infty,+\infty)\).

See the documentation of pm_distGenExpGamma for more information on the GenExpGamma CDF.

Parameters
[in]x: The input scalar or array of the same shape as other array like arguments, of type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128), containing the values at which the CDF must be computed.
[in]kappa: The input scalar or array of the same shape as other array-like arguments, of the same type and kind as x, containing the shape parameter of the distribution.
(optional, default = 1.)
[in]invOmega: The input scalar or array of the same shape as other array-like arguments, of the same type and kind as x, containing the inverse scale parameter of the distribution.
(optional, default = 1.)
[in]logSigma: The input scalar or array of the same shape as other array-like arguments, of the same type and kind as x, containing the location parameter of the distribution.
(optional, default = 1.)
Returns
cdf : The output scalar or array of the same shape as any input array-like argument, of the same type and kind the input argument x, containing the CDF of the distribution at the specified x.


Possible calling interfaces

cdf = getGenExpGammaCDF(x, kappa = kappa, invOmega = invOmega, logSigma = logSigma)
Generate and return the Cumulative Distribution Function (CDF) of the GenExpGamma distribution for an...
This module contains classes and procedures for computing various statistical quantities related to t...
Warning
The condition kappa > 0 and invOmega > 0 must hold for the input values (kappa, invOmega).
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.
The choice of logSigma vs. sigma as the input argument to the procedures of this generic interface may seem rather odd.
This choice was made deliberately to increase the similarity and consistency of the functional interface here with the subroutine interface setGenExpGammaCDF.
Note that the primary reason for preferring logSigma vs. sigma as the choice of input argument is the enhanced computational efficiency of the routines.
Additionally, logSigma matches the alternative rate interpretation of the second parameter of the GenExpGamma distribution.
See also
setGenExpGammaCDF


Example usage

1program example
2
3 use pm_kind, only: SK
4 use pm_kind, only: IK
5 use pm_kind, only: RK => RKS ! all other real kinds are also supported.
6 use pm_io, only: display_type
10
11 implicit none
12
13 integer(IK) , parameter :: NP = 1000_IK
14 real(RK) , allocatable :: point(:), cdf(:), kappa(:), invOmega(:), logSigma(:)
15
16 type(display_type) :: disp
17 disp = display_type(file = "main.out.F90")
18
19 kappa = getLinSpace(+0.5_RK, +2._RK, count = NP)
20 invOmega = getLogSpace(logx1 = log(0.1_RK), logx2 = log(10._RK), count = NP)
21 logSigma = getLinSpace(-3._RK, +3._RK, count = NP)
22 point = getLinSpace(-10._RK, +10._RK, count = NP)
23 allocate(cdf, mold = point)
24
25 call disp%skip()
26 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
28 call disp%show("! Compute the Cumulative Distribution Function (cdf) of GenExpGamma distribution at the specified values.")
29 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
30 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
31 call disp%skip()
32
33 call disp%skip()
34 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
35 call disp%show("! Compute the cdf at an input scalar real value.")
36 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
37 call disp%skip()
38
39 call disp%skip()
40 call disp%show("cdf(1) = getGenExpGammaCDF(0.5_RK)")
41 cdf(1) = getGenExpGammaCDF(0.5_RK)
42 call disp%show("cdf(1)")
43 call disp%show( cdf(1) )
44 call disp%skip()
45
46 call disp%skip()
47 call disp%show("kappa(1)")
48 call disp%show( kappa(1) )
49 call disp%show("cdf(1) = getGenExpGammaCDF(0.5_RK, kappa(1))")
50 cdf(1) = getGenExpGammaCDF(0.5_RK, kappa(1))
51 call disp%show("cdf(1)")
52 call disp%show( cdf(1) )
53 call disp%skip()
54
55 call disp%skip()
56 call disp%show("kappa(1)")
57 call disp%show( kappa(1) )
58 call disp%show("invOmega(1)")
59 call disp%show( invOmega(1) )
60 call disp%show("cdf(1) = getGenExpGammaCDF(0.5_RK, kappa(1), invOmega(1))")
61 cdf(1) = getGenExpGammaCDF(0.5_RK, kappa(1), invOmega(1))
62 call disp%show("cdf(1)")
63 call disp%show( cdf(1) )
64 call disp%skip()
65
66 call disp%skip()
67 call disp%show("kappa(1)")
68 call disp%show( kappa(1) )
69 call disp%show("invOmega(1)")
70 call disp%show( invOmega(1) )
71 call disp%show("cdf(1) = getGenExpGammaCDF(0.5_RK, kappa(1), invOmega(1), logSigma(1))")
72 cdf(1) = getGenExpGammaCDF(0.5_RK, kappa(1), invOmega(1), logSigma(1))
73 call disp%show("cdf(1)")
74 call disp%show( cdf(1) )
75 call disp%skip()
76
77 call disp%skip()
78 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
79 call disp%show("! Compute the cdf at an input vector real value with different parameter values.")
80 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
81 call disp%skip()
82
83 call disp%skip()
84 call disp%show("kappa(1)")
85 call disp%show( kappa(1) )
86 call disp%show("cdf(1:NP:NP/5) = getGenExpGammaCDF(0.5_RK, kappa(1:NP:NP/5))")
87 cdf(1:NP:NP/5) = getGenExpGammaCDF(0.5_RK, kappa(1:NP:NP/5))
88 call disp%show("cdf(1:NP:NP/5)")
89 call disp%show( cdf(1:NP:NP/5) )
90 call disp%skip()
91
92 call disp%skip()
93 call disp%show("kappa(1)")
94 call disp%show( kappa(1) )
95 call disp%show("cdf(1:NP:NP/5) = getGenExpGammaCDF(point(1:NP:NP/5), kappa(1:NP:NP/5))")
96 cdf(1:NP:NP/5) = getGenExpGammaCDF(point(1:NP:NP/5), kappa(1:NP:NP/5))
97 call disp%show("cdf(1:NP:NP/5)")
98 call disp%show( cdf(1:NP:NP/5) )
99 call disp%skip()
100
101 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
102 ! Output an example array for visualization.
103 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
104
105 block
106 real(RK) :: cdf(NP, 6)
107 integer(IK) :: i, fileUnit
108 cdf(:,1) = getGenExpGammaCDF(point, .5_RK, 2._RK)
109 cdf(:,2) = getGenExpGammaCDF(point, 1._RK, 2._RK)
110 cdf(:,3) = getGenExpGammaCDF(point, 2._RK, 2._RK)
111 cdf(:,4) = getGenExpGammaCDF(point, .5_RK, .5_RK)
112 cdf(:,5) = getGenExpGammaCDF(point, 1._RK, .5_RK)
113 cdf(:,6) = getGenExpGammaCDF(point, 2._RK, .5_RK)
114 open(newunit = fileUnit, file = "getGenExpGammaCDF.RK.txt")
115 write(fileUnit,"(7(g0,:,' '))") (point(i), cdf(i,:), i = 1, size(point))
116 close(fileUnit)
117 end block
118
119end program example
Generate count evenly spaced points over the interval [x1, x2] if x1 < x2, or [x2,...
Generate count evenly-logarithmically-spaced points over the interval [base**logx1,...
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 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 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
2!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
4! Compute the Cumulative Distribution Function (cdf) of GenExpGamma distribution at the specified values.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
10! Compute the cdf at an input scalar real value.
11!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
12
13
14cdf(1) = getGenExpGammaCDF(0.5_RK)
15cdf(1)
16+0.807704329
17
18
19kappa(1)
20+0.500000000
21cdf(1) = getGenExpGammaCDF(0.5_RK, kappa(1))
22cdf(1)
23+0.930612206
24
25
26kappa(1)
27+0.500000000
28invOmega(1)
29+0.999999940E-1
30cdf(1) = getGenExpGammaCDF(0.5_RK, kappa(1), invOmega(1))
31cdf(1)
32+0.852945745
33
34
35kappa(1)
36+0.500000000
37invOmega(1)
38+0.999999940E-1
39cdf(1) = getGenExpGammaCDF(0.5_RK, kappa(1), invOmega(1), logSigma(1))
40cdf(1)
41+0.907949090
42
43
44!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
45! Compute the cdf at an input vector real value with different parameter values.
46!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
47
48
49kappa(1)
50+0.500000000
51cdf(1:NP:NP/5) = getGenExpGammaCDF(0.5_RK, kappa(1:NP:NP/5))
52cdf(1:NP:NP/5)
53+0.930612206, +0.862185597, +0.778141320, +0.684153199, +0.586314142
54
55
56kappa(1)
57+0.500000000
58cdf(1:NP:NP/5) = getGenExpGammaCDF(point(1:NP:NP/5), kappa(1:NP:NP/5))
59cdf(1:NP:NP/5)
60+0.760284485E-2, +0.883790851E-2, +0.993758515E-1, +0.998501480, +1.00000000
61
62

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
11marker ={ "CK" : "-"
12 , "IK" : "."
13 , "RK" : "-"
14 }
15xlab = { "CK" : "X ( real/imaginary components )"
16 , "IK" : "X ( integer-valued )"
17 , "RK" : "X ( real-valued )"
18 }
19legends = [ "$\kappa = 0.5, \omega = 0.5$"
20 , "$\kappa = 1.0, \omega = 0.5$"
21 , "$\kappa = 2.0, \omega = 0.5$"
22 , "$\kappa = 0.5, \omega = 2.0$"
23 , "$\kappa = 1.0, \omega = 2.0$"
24 , "$\kappa = 2.0, \omega = 2.0$"
25 ]
26
27for kind in ["IK", "CK", "RK"]:
28
29 pattern = "*." + kind + ".txt"
30 fileList = glob.glob(pattern)
31 if len(fileList) == 1:
32
33 df = pd.read_csv(fileList[0], delimiter = " ")
34
35 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
36 ax = plt.subplot()
37
38 if kind == "CK":
39 plt.plot( df.values[:, 0]
40 , df.values[:,1:len(legends)+1]
41 , marker[kind]
42 #, color = "r"
43 )
44 plt.plot( df.values[:, 1]
45 , df.values[:,1:len(legends)+1]
46 , marker[kind]
47 #, color = "blue"
48 )
49 else:
50 plt.plot( df.values[:, 0]
51 , df.values[:,1:len(legends)+1]
52 , marker[kind]
53 #, color = "r"
54 )
55 ax.legend ( legends
56 , fontsize = fontsize
57 )
58
59 plt.xticks(fontsize = fontsize - 2)
60 plt.yticks(fontsize = fontsize - 2)
61 ax.set_xlabel(xlab[kind], fontsize = 17)
62 ax.set_ylabel("Cumulative Distribution Function (CDF)", fontsize = 17)
63
64 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
65 ax.tick_params(axis = "y", which = "minor")
66 ax.tick_params(axis = "x", which = "minor")
67
68 plt.savefig(fileList[0].replace(".txt",".png"))
69
70 elif len(fileList) > 1:
71
72 sys.exit("Ambiguous file list exists.")

Visualization of the example output
Test:
test_pm_distGenExpGamma
Todo:
Low 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, Oct 16, 2009, 11:14 AM, Michigan

Definition at line 933 of file pm_distGenExpGamma.F90.


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