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

Generate and return the natural logarithm of the Probability Density Function (PDF) of the ExpGamma distribution. More...

Detailed Description

Generate and return the natural logarithm of the Probability Density Function (PDF) of the ExpGamma distribution.

See the documentation of pm_distExpGamma for more information on the PDF of the ExpGamma distribution.

Parameters
[in]x: The input scalar or array of the same shape as other array-valued arguments, containing the values at which the log(PDF) must be computed.
[in]kappa: The input scalar or array of the same shape as other array-valued arguments, containing the shape parameter ( \(\kappa\)) of the distribution.
(optional, default = 1.)
[in]logSigma: The input scalar or array of the same shape as other array-valued arguments, containing the location parameter ( \(\log(\sigma)\)) of the distribution.
(optional, default = 0)
Returns
logPDF : The output scalar or array of the same shape as other array-valued input arguments of the same type and kind as x containing the PDF of the specified ExpGamma distribution.


Possible calling interfaces

logPDF = getExpGammaLogPDF(x, kappa = kappa, logSigma = logSigma)
Generate and return the natural logarithm of the Probability Density Function (PDF) of the ExpGamma d...
This module contains classes and procedures for computing various statistical quantities related to t...
Warning
The conditions kappa > 0 must hold for the corresponding input arguments.
This condition is 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
setExpGammaLogPDF


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(:), logPDF(:), Kappa(:), 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 LogSigma = getLinSpace(-3._RK, +3._RK, count = NP)
21 Point = getLinSpace(-10._RK, +10._RK, count = NP)
22 allocate(logPDF, mold = Point)
23
24 call disp%skip()
25 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
26 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%show("! Compute the Probability Density Function (PDF) of ExpGamma distribution at the specified values.")
28 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
29 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
30 call disp%skip()
31
32 call disp%skip()
33 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
34 call disp%show("! Compute the PDF at an input scalar real value.")
35 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
36 call disp%skip()
37
38 call disp%skip()
39 call disp%show("logPDF(1) = getExpGammaLogPDF(0.5_RK)")
40 logPDF(1) = getExpGammaLogPDF(0.5_RK)
41 call disp%show("logPDF(1)")
42 call disp%show( logPDF(1) )
43 call disp%skip()
44
45 call disp%skip()
46 call disp%show("Kappa(1)")
47 call disp%show( Kappa(1) )
48 call disp%show("logPDF(1) = getExpGammaLogPDF(0.5_RK, Kappa(1))")
49 logPDF(1) = getExpGammaLogPDF(0.5_RK, Kappa(1))
50 call disp%show("logPDF(1)")
51 call disp%show( logPDF(1) )
52 call disp%skip()
53
54 call disp%skip()
55 call disp%show("Kappa(1)")
56 call disp%show( Kappa(1) )
57 call disp%show("logPDF(1) = getExpGammaLogPDF(0.5_RK, Kappa(1), LogSigma(1))")
58 logPDF(1) = getExpGammaLogPDF(0.5_RK, Kappa(1), LogSigma(1))
59 call disp%show("logPDF(1)")
60 call disp%show( logPDF(1) )
61 call disp%skip()
62
63 call disp%skip()
64 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
65 call disp%show("! Compute the PDF at an input vector real value with different parameter values.")
66 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
67 call disp%skip()
68
69 call disp%skip()
70 call disp%show("Kappa(1)")
71 call disp%show( Kappa(1) )
72 call disp%show("logPDF(1:NP:NP/5) = getExpGammaLogPDF(0.5_RK, Kappa(1:NP:NP/5))")
73 logPDF(1:NP:NP/5) = getExpGammaLogPDF(0.5_RK, Kappa(1:NP:NP/5))
74 call disp%show("logPDF(1:NP:NP/5)")
75 call disp%show( logPDF(1:NP:NP/5) )
76 call disp%skip()
77
78 call disp%skip()
79 call disp%show("Kappa(1)")
80 call disp%show( Kappa(1) )
81 call disp%show("logPDF(1:NP:NP/5) = getExpGammaLogPDF(Point(1:NP:NP/5), Kappa(1:NP:NP/5))")
82 logPDF(1:NP:NP/5) = getExpGammaLogPDF(Point(1:NP:NP/5), Kappa(1:NP:NP/5))
83 call disp%show("logPDF(1:NP:NP/5)")
84 call disp%show( logPDF(1:NP:NP/5) )
85 call disp%skip()
86
87 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
88 ! Output an example array for visualization.
89 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
90
91 block
92 integer :: fileUnit, i
93 real :: logPDF(NP,6)
94 logPDF(:,1) = getExpGammaLogPDF(Point, .5, -.5)
95 logPDF(:,2) = getExpGammaLogPDF(Point, 1., -.5)
96 logPDF(:,3) = getExpGammaLogPDF(Point, 2., -.5)
97 logPDF(:,4) = getExpGammaLogPDF(Point, .5, 2.0)
98 logPDF(:,5) = getExpGammaLogPDF(Point, 1., 2.0)
99 logPDF(:,6) = getExpGammaLogPDF(Point, 2., 2.0)
100 open(newunit = fileUnit, file = "getExpGammaLogPDF.RK.txt")
101 write(fileUnit,"(7(g0,:,' '))") (Point(i), exp(logPDF(i,:)), i = 1, size(Point))
102 close(fileUnit)
103 end block
104
105end 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 Probability Density Function (PDF) of ExpGamma distribution at the specified values.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
10! Compute the PDF at an input scalar real value.
11!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
12
13
14logPDF(1) = getExpGammaLogPDF(0.5_RK)
15logPDF(1)
16-1.14872122
17
18
19Kappa(1)
20+0.500000000
21logPDF(1) = getExpGammaLogPDF(0.5_RK, Kappa(1))
22logPDF(1)
23-1.97108614
24
25
26Kappa(1)
27+0.500000000
28logPDF(1) = getExpGammaLogPDF(0.5_RK, Kappa(1), LogSigma(1))
29logPDF(1)
30-31.9378166
31
32
33!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
34! Compute the PDF at an input vector real value with different parameter values.
35!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
36
37
38Kappa(1)
39+0.500000000
40logPDF(1:NP:NP/5) = getExpGammaLogPDF(0.5_RK, Kappa(1:NP:NP/5))
41logPDF(1:NP:NP/5)
42-1.97108614, -1.40034103, -1.04829431, -0.828603029, -0.702564001
43
44
45Kappa(1)
46+0.500000000
47logPDF(1:NP:NP/5) = getExpGammaLogPDF(Point(1:NP:NP/5), Kappa(1:NP:NP/5))
48logPDF(1:NP:NP/5)
49-5.57241011, -4.95285606, -2.27868414, -4.54004860, -399.612122
50
51

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, \log(\sigma) = -.5$"
20 , "$\kappa = 1.0, \log(\sigma) = -.5$"
21 , "$\kappa = 2.0, \log(\sigma) = -.5$"
22 , "$\kappa = 0.5, \log(\sigma) = 2.0$"
23 , "$\kappa = 1.0, \log(\sigma) = 2.0$"
24 , "$\kappa = 2.0, \log(\sigma) = 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("Probability Density Function (PDF)", fontsize = 17)
63 ax.set_xlim([-11.0, 5])
64
65 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
66 ax.tick_params(axis = "y", which = "minor")
67 ax.tick_params(axis = "x", which = "minor")
68
69 plt.savefig(fileList[0].replace(".txt",".png"))
70
71 elif len(fileList) > 1:
72
73 sys.exit("Ambiguous file list exists.")

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

Definition at line 322 of file pm_distExpGamma.F90.


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