Generate and return the natural logarithm of the normalization factor of the Probability Density Function (PDF) of the GenExpGamma distribution.
More...
Generate and return the natural logarithm of the normalization factor of the Probability Density Function (PDF) of the GenExpGamma distribution.
The normalization factor of the GenExpGamma PDF is defined as \(\eta = 1 / (\omega \Gamma(\kappa))\) where \((\omega, \kappa)\) are the scale and shape parameters of the GenExpGamma distribution, respectively.
For more information, see the documentation of pm_distGenExpGamma.
- Parameters
-
[in] | kappa | : The input scalar or array of the same shape as other array-like arguments, containing the shape parameter ( \(\kappa\)) of the distribution.
|
[in] | invOmega | : The input scalar or array of the same shape as other array-like arguments, containing the inverse of the scale parameter ( \(\omega\)) of the distribution.
(optional, default = 1. ) |
- Returns
logPDFNF
: The output scalar or array of the same shape as array-like input arguments, of the same type, kind, and highest rank as the input arguments containing the natural logarithm of the normalization factor of the distribution.
Possible calling interfaces ⛓
Generate and return the natural logarithm of the normalization factor of the Probability Density Func...
This module contains classes and procedures for computing various statistical quantities related to t...
- Warning
- 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.
- See also
- setGenExpGammaLogPDF
setGenExpGammaLogPDF
Example usage ⛓
13 integer(IK) ,
parameter :: NP
= 1000_IK
14 real(RK) ,
allocatable :: invOmega(:), logPDFNF(:), Kappa(:)
16 type(display_type) :: disp
21 allocate(logPDFNF,
mold = Kappa)
24 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
25 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
26 call disp%show(
"! Compute the natural logarithm of the normalization factor of the GenExpGamma distribution PDF.")
27 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
28 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
32 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
33 call disp%show(
"! Compute the PDF at an input scalar real value.")
34 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
40 call disp%show(
"logPDFNF(1) = getGenExpGammaLogPDFNF(kappa = Kappa(1))")
51 call disp%show(
"logPDFNF(1) = getGenExpGammaLogPDFNF(Kappa(1), invOmega(1))")
58 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
59 call disp%show(
"! Compute the PDF at at a mix of scalar and vector input values.")
60 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
66 call disp%show(
"logPDFNF(1:NP:NP/5) = getGenExpGammaLogPDFNF(Kappa(1:NP:NP/5))")
68 call disp%show(
"logPDFNF(1:NP:NP/5)")
69 call disp%show( logPDFNF(
1:NP:NP
/5) )
75 call disp%show(
"logPDFNF(1:NP:NP/5) = getGenExpGammaLogPDFNF(Kappa(1:NP:NP/5), 1._RK)")
77 call disp%show(
"logPDFNF(1:NP:NP/5)")
78 call disp%show( logPDFNF(
1:NP:NP
/5) )
84 call disp%show(
"invOmega(1:NP:NP/5)")
85 call disp%show( invOmega(
1:NP:NP
/5) )
86 call disp%show(
"logPDFNF(1:NP:NP/5) = getGenExpGammaLogPDFNF(1._RK, invOmega(1:NP:NP/5))")
88 call disp%show(
"logPDFNF(1:NP:NP/5)")
89 call disp%show( logPDFNF(
1:NP:NP
/5) )
95 call disp%show(
"invOmega(1:NP:NP/5)")
96 call disp%show( invOmega(
1:NP:NP
/5) )
97 call disp%show(
"logPDFNF(1:NP:NP/5) = getGenExpGammaLogPDFNF(Kappa(1:NP:NP/5), invOmega(1:NP:NP/5))")
99 call disp%show(
"logPDFNF(1:NP:NP/5)")
100 call disp%show( logPDFNF(
1:NP:NP
/5) )
108 integer :: fileUnit, i
110 open(newunit
= fileUnit, file
= "getGenExpGammaLogPDFNF.RK.txt")
111 write(fileUnit,
"(2(g0,:,' '))") (Kappa(i),
exp(logPDFNF(i)), i
= 1,
size(logPDFNF))
Generate count evenly spaced points over the interval [x1, x2] if x1 < x2, or [x2,...
Return the linSpace output argument with size(linSpace) elements of evenly-spaced values over the int...
This is a generic method of the derived type display_type with pass attribute.
This is a generic method of the derived type display_type with pass attribute.
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...
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
Generate and return an object of type display_type.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example output ⛓
36+0.999999978E-2,
+2.00999999,
+4.01000023,
+6.01000023,
+8.01000023
39-4.59948015,
-0.426001893E-2,
-1.80433512,
-4.80456209,
-8.54532528
43+0.999999978E-2,
+2.00999999,
+4.01000023,
+6.01000023,
+8.01000023
46-4.59948015,
-0.426001893E-2,
-1.80433512,
-4.80456209,
-8.54532528
50+0.999999978E-2,
+2.00999999,
+4.01000023,
+6.01000023,
+8.01000023
52+0.100000001,
+2.08198190,
+4.06396389,
+6.04594564,
+8.02792740
55-2.30258512,
+0.733320296,
+1.40215886,
+1.79938793,
+2.08292627
59+0.999999978E-2,
+2.00999999,
+4.01000023,
+6.01000023,
+8.01000023
61+0.100000001,
+2.08198190,
+4.06396389,
+6.04594564,
+8.02792740
64-6.90206528,
+0.729060292,
-0.402176261,
-3.00517416,
-6.46239901
Postprocessing of the example output ⛓
3import matplotlib.pyplot
as plt
15xlab = {
"CK" :
r"shape parameter: $\kappa$ ( real/imaginary )"
16 ,
"IK" :
r"shape parameter: $\kappa$ ( integer-valued )"
17 ,
"RK" :
r"shape parameter: $\kappa$ ( real-valued )"
20for kind
in [
"IK",
"CK",
"RK"]:
22 pattern =
"*." + kind +
".txt"
23 fileList = glob.glob(pattern)
24 if len(fileList) == 1:
26 df = pd.read_csv(fileList[0], delimiter =
" ")
28 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
32 plt.plot( df.values[:, 0]
37 plt.plot( df.values[:, 1]
43 plt.plot( df.values[:, 0]
49 plt.xticks(fontsize = fontsize - 2)
50 plt.yticks(fontsize = fontsize - 2)
51 ax.set_xlabel(xlab[kind], fontsize = fontsize)
52 ax.set_ylabel(
"Normalization Factor of the PDF", fontsize = fontsize)
54 plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
55 ax.tick_params(axis =
"y", which =
"minor")
56 ax.tick_params(axis =
"x", which =
"minor")
58 plt.savefig(fileList[0].replace(
".txt",
".png"))
60 elif len(fileList) > 1:
62 sys.exit(
"Ambiguous file list exists.")
Visualization of the example output ⛓
- Benchmarks:
Benchmark :: The runtime performance of setGenExpGammaLogPDF with and without scale
⛓
3 use iso_fortran_env,
only:
error_unit
11 integer(IK) :: fileUnit
12 integer(IK) ,
parameter :: NSIZE
= 15_IK
13 integer(IK) ,
parameter :: NBENCH
= 3_IK
14 integer(IK) :: arraySize(NSIZE)
15 real(RK) ,
allocatable :: InvScale(:)
16 real(RK) ,
allocatable ::
shape(:)
17 real(RK) :: dummy
= 0._RK
18 type(bench_type) :: Bench(NBENCH)
20 Bench(
1)
= bench_type(name
= SK_
"scalarAddition", exec
= scalarAddition, overhead
= setOverhead)
21 Bench(
2)
= bench_type(name
= SK_
"logPDFNFWithShape", exec
= logPDFNFWithShape, overhead
= setOverhead)
22 Bench(
3)
= bench_type(name
= SK_
"logPDFNFWithShapeScale", exec
= logPDFNFWithShapeScale, overhead
= setOverhead)
24 arraySize
= [(
2**isize, isize
= 1, NSIZE )]
26 write(
*,
"(*(g0,:,' '))")
27 write(
*,
"(*(g0,:,' vs. '))") (Bench(i)
%name, i
= 1, NBENCH)
28 write(
*,
"(*(g0,:,' '))")
30 open(newunit
= fileUnit, file
= "main.out", status
= "replace")
32 write(fileUnit,
"(*(g0,:,','))")
"arraySize", (Bench(i)
%name, i
= 1, NBENCH)
34 loopOverArraySize:
do isize
= 1, NSIZE
36 write(
*,
"(*(g0,:,' '))")
"Benchmarking with rank", arraySize(isize)
38 allocate(
shape(arraySize(isize)), InvScale(arraySize(isize)))
40 Bench(i)
%timing
= Bench(i)
%getTiming(minsec
= 0.05_RK)
42 deallocate(shape, InvScale)
44 write(fileUnit,
"(*(g0,:,','))") arraySize(isize), (Bench(i)
%timing
%mean, i
= 1, NBENCH)
46 end do loopOverArraySize
48 write(
*,
"(*(g0,:,' '))") dummy
49 write(
*,
"(*(g0,:,' '))")
59 subroutine setOverhead()
63 subroutine initialize()
64 call random_number(shape)
65 shape
= shape
+ epsilon(
0._RK)
66 InvScale
= shape
+ 1._RK
69 subroutine scalarAddition()
71 dummy
= dummy
+ sum(scalarAdditionCallee(shape,InvScale))
74 impure elemental function scalarAdditionCallee(shape, invScale)
result(summation)
75 real(RK),
intent(in) :: shape, invScale
77 summation
= shape
+ invScale
80 subroutine logPDFNFWithShape()
86 subroutine logPDFNFWithShapeScale()
Generate and return an object of type timing_type containing the benchmark timing information and sta...
This module contains abstract interfaces and types that facilitate benchmarking of different procedur...
This is the class for creating benchmark and performance-profiling objects.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Postprocessing of the benchmark output ⛓
3import matplotlib.pyplot
as plt
8dirname = os.path.basename(os.getcwd())
12df = pd.read_csv(
"main.out", delimiter =
",")
13colnames = list(df.columns.values)
19ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
22for colname
in colnames[1:]:
23 plt.plot( df[colnames[0]].values
28plt.xticks(fontsize = fontsize)
29plt.yticks(fontsize = fontsize)
30ax.set_xlabel(colnames[0], fontsize = fontsize)
31ax.set_ylabel(
"Runtime [ seconds ]", fontsize = fontsize)
32ax.set_title(
" vs. ".join(colnames[1:])+
"\nLower is better.", fontsize = fontsize)
36plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
37ax.tick_params(axis =
"y", which =
"minor")
38ax.tick_params(axis =
"x", which =
"minor")
39ax.legend ( colnames[1:]
46plt.savefig(
"benchmark." + dirname +
".runtime.png")
52ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
55plt.plot( df[colnames[0]].values
56 , np.ones(len(df[colnames[0]].values))
61for colname
in colnames[2:]:
62 plt.plot( df[colnames[0]].values
63 , df[colname].values / df[colnames[1]].values
67plt.xticks(fontsize = fontsize)
68plt.yticks(fontsize = fontsize)
69ax.set_xlabel(colnames[0], fontsize = fontsize)
70ax.set_ylabel(
"Runtime compared to {}".format(colnames[1]), fontsize = fontsize)
71ax.set_title(
"Runtime Ratio Comparison. Lower means faster.\nLower than 1 means faster than {}().".format(colnames[1]), fontsize = fontsize)
75plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
76ax.tick_params(axis =
"y", which =
"minor")
77ax.tick_params(axis =
"x", which =
"minor")
78ax.legend ( colnames[1:]
85plt.savefig(
"benchmark." + dirname +
".runtime.ratio.png")
Visualization of the benchmark output ⛓
Benchmark moral ⛓
- The procedures under the generic interface getGenExpGammaLogPDFNF aim to precompute the natural logarithm of the normalization factor of the GenExpGamma PDF.
Based on the results of this benchmark, passing this precomputed normalization factor to getGenExpGammaLogPDF or setGenExpGammaLogPDF appears to lead to significant performance improvement (20-30 times) when the GenExpGamma PDF procedures are to called repeatedly with the same set of input distribution parameters (kappa
, sigma
).
This speedup is entirely due to avoiding the redundant repeated computations of log_gamma()
and log()
terms in getGenExpGammaLogPDFNF.
- The initial closing in the performance gap for small array sizes (of \(\sim4\) or less elements) is likely artificial and due to extremely small cost of
simpleAddition
procedure which is likely less than the time resolution of the processor ( \(\sim10ns\) as of V2 of ParaMonte library).
- Test:
- test_pm_distGenExpGamma
- 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.
-
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.
-
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.
- Copyright
- Computational Data Science Lab
- Author:
- Amir Shahmoradi, Oct 16, 2009, 11:14 AM, Michigan
Definition at line 270 of file pm_distGenExpGamma.F90.