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

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

Detailed Description

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

logPDFNF = getGenExpGammaLogPDFNF(kappa)
logPDFNF = getGenExpGammaLogPDFNF(kappa, invOmega)
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.
Remarks
The procedures under discussion are elemental.
These procedures are particularly useful and needed in computing the PDF of the distribution.
The logic behind pre-computing the normalization factor is to speed up the calculations since the normalization factor contains log() and log_gamma() computations which are computationally costly operations.
See the benchmark below for more details.
See also
setGenExpGammaLogPDF
setGenExpGammaLogPDF


Example usage

1program example
2
3 use pm_kind, only: IK
4 use pm_kind, only: SK
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 :: invOmega(:), logPDFNF(:), Kappa(:)
15
16 type(display_type) :: disp
17 disp = display_type(file = "main.out.F90")
18
19 Kappa = getLinSpace(0.01_RK, 10._RK, count = NP)
20 invOmega = getLinSpace(0.1_RK, 10._RK, count = NP)
21 allocate(logPDFNF, mold = Kappa)
22
23 call disp%skip()
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("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
29 call disp%skip()
30
31 call disp%skip()
32 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
33 call disp%show("! Compute the PDF at an input scalar real value.")
34 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
35 call disp%skip()
36
37 call disp%skip()
38 call disp%show("Kappa(1)")
39 call disp%show( Kappa(1) )
40 call disp%show("logPDFNF(1) = getGenExpGammaLogPDFNF(kappa = Kappa(1))")
41 logPDFNF(1) = getGenExpGammaLogPDFNF(kappa = Kappa(1))
42 call disp%show("logPDFNF(1)")
43 call disp%show( logPDFNF(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("invOmega(1)")
50 call disp%show( invOmega(1) )
51 call disp%show("logPDFNF(1) = getGenExpGammaLogPDFNF(Kappa(1), invOmega(1))")
52 logPDFNF(1) = getGenExpGammaLogPDFNF(Kappa(1), invOmega(1))
53 call disp%show("logPDFNF(1)")
54 call disp%show( logPDFNF(1) )
55 call disp%skip()
56
57 call disp%skip()
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("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
61 call disp%skip()
62
63 call disp%skip()
64 call disp%show("Kappa(1:NP:NP/5)")
65 call disp%show( Kappa(1:NP:NP/5) )
66 call disp%show("logPDFNF(1:NP:NP/5) = getGenExpGammaLogPDFNF(Kappa(1:NP:NP/5))")
67 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) )
70 call disp%skip()
71
72 call disp%skip()
73 call disp%show("Kappa(1:NP:NP/5)")
74 call disp%show( Kappa(1:NP:NP/5) )
75 call disp%show("logPDFNF(1:NP:NP/5) = getGenExpGammaLogPDFNF(Kappa(1:NP:NP/5), 1._RK)")
76 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) )
79 call disp%skip()
80
81 call disp%skip()
82 call disp%show("Kappa(1:NP:NP/5)")
83 call disp%show( Kappa(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))")
87 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) )
90 call disp%skip()
91
92 call disp%skip()
93 call disp%show("Kappa(1:NP:NP/5)")
94 call disp%show( Kappa(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))")
98 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) )
101 call disp%skip()
102
103 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
104 ! Output an example array for visualization.
105 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
106
107 block
108 integer :: fileUnit, i
109 logPDFNF = getGenExpGammaLogPDFNF(Kappa)
110 open(newunit = fileUnit, file = "getGenExpGammaLogPDFNF.RK.txt")
111 write(fileUnit,"(2(g0,:,' '))") (Kappa(i), exp(logPDFNF(i)), i = 1, size(logPDFNF))
112 close(fileUnit)
113 end block
114
115end program example
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.
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 natural logarithm of the normalization factor of the GenExpGamma distribution PDF.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
10! Compute the PDF at an input scalar real value.
11!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
12
13
14Kappa(1)
15+0.999999978E-2
16logPDFNF(1) = getGenExpGammaLogPDFNF(kappa = Kappa(1))
17logPDFNF(1)
18-4.59948015
19
20
21Kappa(1)
22+0.999999978E-2
23invOmega(1)
24+0.100000001
25logPDFNF(1) = getGenExpGammaLogPDFNF(Kappa(1), invOmega(1))
26logPDFNF(1)
27-6.90206528
28
29
30!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
31! Compute the PDF at at a mix of scalar and vector input values.
32!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
33
34
35Kappa(1:NP:NP/5)
36+0.999999978E-2, +2.00999999, +4.01000023, +6.01000023, +8.01000023
37logPDFNF(1:NP:NP/5) = getGenExpGammaLogPDFNF(Kappa(1:NP:NP/5))
38logPDFNF(1:NP:NP/5)
39-4.59948015, -0.426001893E-2, -1.80433512, -4.80456209, -8.54532528
40
41
42Kappa(1:NP:NP/5)
43+0.999999978E-2, +2.00999999, +4.01000023, +6.01000023, +8.01000023
44logPDFNF(1:NP:NP/5) = getGenExpGammaLogPDFNF(Kappa(1:NP:NP/5), 1._RK)
45logPDFNF(1:NP:NP/5)
46-4.59948015, -0.426001893E-2, -1.80433512, -4.80456209, -8.54532528
47
48
49Kappa(1:NP:NP/5)
50+0.999999978E-2, +2.00999999, +4.01000023, +6.01000023, +8.01000023
51invOmega(1:NP:NP/5)
52+0.100000001, +2.08198190, +4.06396389, +6.04594564, +8.02792740
53logPDFNF(1:NP:NP/5) = getGenExpGammaLogPDFNF(1._RK, invOmega(1:NP:NP/5))
54logPDFNF(1:NP:NP/5)
55-2.30258512, +0.733320296, +1.40215886, +1.79938793, +2.08292627
56
57
58Kappa(1:NP:NP/5)
59+0.999999978E-2, +2.00999999, +4.01000023, +6.01000023, +8.01000023
60invOmega(1:NP:NP/5)
61+0.100000001, +2.08198190, +4.06396389, +6.04594564, +8.02792740
62logPDFNF(1:NP:NP/5) = getGenExpGammaLogPDFNF(Kappa(1:NP:NP/5), invOmega(1:NP:NP/5))
63logPDFNF(1:NP:NP/5)
64-6.90206528, +0.729060292, -0.402176261, -3.00517416, -6.46239901
65
66

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" : r"shape parameter: $\kappa$ ( real/imaginary )"
16 , "IK" : r"shape parameter: $\kappa$ ( integer-valued )"
17 , "RK" : r"shape parameter: $\kappa$ ( real-valued )"
18 }
19
20for kind in ["IK", "CK", "RK"]:
21
22 pattern = "*." + kind + ".txt"
23 fileList = glob.glob(pattern)
24 if len(fileList) == 1:
25
26 df = pd.read_csv(fileList[0], delimiter = " ")
27
28 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
29 ax = plt.subplot()
30
31 if kind == "CK":
32 plt.plot( df.values[:, 0]
33 , df.values[:,2]
34 , marker[kind]
35 , color = "r"
36 )
37 plt.plot( df.values[:, 1]
38 , df.values[:,3]
39 , marker[kind]
40 , color = "blue"
41 )
42 else:
43 plt.plot( df.values[:, 0]
44 , df.values[:, 1]
45 , marker[kind]
46 , color = "r"
47 )
48
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)
53
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")
57
58 plt.savefig(fileList[0].replace(".txt",".png"))
59
60 elif len(fileList) > 1:
61
62 sys.exit("Ambiguous file list exists.")

Visualization of the example output
Benchmarks:


Benchmark :: The runtime performance of setGenExpGammaLogPDF with and without scale

1program benchmark
2
3 use iso_fortran_env, only: error_unit
4 use pm_bench, only: bench_type
5 use pm_kind, only: IK, RK, SK
6
7 implicit none
8
9 integer(IK) :: i
10 integer(IK) :: isize
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)
19
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)
23
24 arraySize = [( 2**isize, isize = 1, NSIZE )]
25
26 write(*,"(*(g0,:,' '))")
27 write(*,"(*(g0,:,' vs. '))") (Bench(i)%name, i = 1, NBENCH)
28 write(*,"(*(g0,:,' '))")
29
30 open(newunit = fileUnit, file = "main.out", status = "replace")
31
32 write(fileUnit, "(*(g0,:,','))") "arraySize", (Bench(i)%name, i = 1, NBENCH)
33
34 loopOverArraySize: do isize = 1, NSIZE
35
36 write(*,"(*(g0,:,' '))") "Benchmarking with rank", arraySize(isize)
37
38 allocate(shape(arraySize(isize)), InvScale(arraySize(isize)))
39 do i = 1, NBENCH
40 Bench(i)%timing = Bench(i)%getTiming(minsec = 0.05_RK)
41 end do
42 deallocate(shape, InvScale)
43
44 write(fileUnit,"(*(g0,:,','))") arraySize(isize), (Bench(i)%timing%mean, i = 1, NBENCH)
45
46 end do loopOverArraySize
47
48 write(*,"(*(g0,:,' '))") dummy
49 write(*,"(*(g0,:,' '))")
50
51 close(fileUnit)
52
53contains
54
55 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
56 ! procedure wrappers.
57 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
58
59 subroutine setOverhead()
60 call initialize()
61 end subroutine
62
63 subroutine initialize()
64 call random_number(shape)
65 shape = shape + epsilon(0._RK)
66 InvScale = shape + 1._RK
67 end subroutine
68
69 subroutine scalarAddition()
70 call initialize()
71 dummy = dummy + sum(scalarAdditionCallee(shape,InvScale))
72 end subroutine
73
74 impure elemental function scalarAdditionCallee(shape, invScale) result(summation)
75 real(RK), intent(in) :: shape, invScale
76 real(RK) :: summation
77 summation = shape + invScale
78 end function
79
80 subroutine logPDFNFWithShape()
82 call initialize()
83 dummy = dummy + sum(getGenExpGammaLogPDFNF(shape))
84 end subroutine
85
86 subroutine logPDFNFWithShapeScale()
88 call initialize()
89 dummy = dummy + sum(getGenExpGammaLogPDFNF(shape, InvScale))
90 end subroutine
91
92end program benchmark
Generate and return an object of type timing_type containing the benchmark timing information and sta...
Definition: pm_bench.F90:574
This module contains abstract interfaces and types that facilitate benchmarking of different procedur...
Definition: pm_bench.F90:41
This is the class for creating benchmark and performance-profiling objects.
Definition: pm_bench.F90:386

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

Postprocessing of the benchmark output
1#!/usr/bin/env python
2
3import matplotlib.pyplot as plt
4import pandas as pd
5import numpy as np
6
7import os
8dirname = os.path.basename(os.getcwd())
9
10fontsize = 14
11
12df = pd.read_csv("main.out", delimiter = ",")
13colnames = list(df.columns.values)
14
15
18
19ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
20ax = plt.subplot()
21
22for colname in colnames[1:]:
23 plt.plot( df[colnames[0]].values
24 , df[colname].values
25 , linewidth = 2
26 )
27
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)
33ax.set_xscale("log")
34ax.set_yscale("log")
35plt.minorticks_on()
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:]
40 #, loc='center left'
41 #, bbox_to_anchor=(1, 0.5)
42 , fontsize = fontsize
43 )
44
45plt.tight_layout()
46plt.savefig("benchmark." + dirname + ".runtime.png")
47
48
51
52ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
53ax = plt.subplot()
54
55plt.plot( df[colnames[0]].values
56 , np.ones(len(df[colnames[0]].values))
57 , linestyle = "--"
58 #, color = "black"
59 , linewidth = 2
60 )
61for colname in colnames[2:]:
62 plt.plot( df[colnames[0]].values
63 , df[colname].values / df[colnames[1]].values
64 , linewidth = 2
65 )
66
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)
72ax.set_xscale("log")
73#ax.set_yscale("log")
74plt.minorticks_on()
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:]
79 #, bbox_to_anchor = (1, 0.5)
80 #, loc = "center left"
81 , fontsize = fontsize
82 )
83
84plt.tight_layout()
85plt.savefig("benchmark." + dirname + ".runtime.ratio.png")

Visualization of the benchmark output

Benchmark moral
  1. 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.
  2. 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.

  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 270 of file pm_distGenExpGamma.F90.


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