Generate and return the natural logarithm of the Probability Density Function (PDF) of the GenGamma distribution.
More...
Generate and return the natural logarithm of the Probability Density Function (PDF) of the GenGamma distribution.
See the documentation of pm_distGenGamma for more information on the PDF of the GenGamma 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] | invOmega | : The input scalar or array of the same shape as other array-valued arguments, containing the inverse of the second shape parameter ( \(\omega\)) of the distribution.
(optional, default = 1. ) |
[in] | invSigma | : The input scalar or array of the same shape as other array-valued arguments, containing the inverse of the scale parameter ( \(\sigma\)) of the distribution.
(optional, default = 1. ) |
- 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 GenGamma distribution.
Possible calling interfaces ⛓
logPDF
= getGenGammaLogPDF(x, kappa
= kappa, invOmega
= invOmega, invSigma
= invSigma)
Generate and return the natural logarithm of the Probability Density Function (PDF) of the GenGamma d...
This module contains classes and procedures for computing various statistical quantities related to t...
- Warning
- The conditions
kappa > 0
, invOmega > 0
, and invSigma > 0
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.
- See also
- setGenGammaLogPDF
Example usage ⛓
13 integer(IK) ,
parameter :: NP
= 1000_IK
14 real(RK) ,
allocatable :: Point(:), logPDF(:), Kappa(:), invOmega(:), invSigma(:)
16 type(display_type) :: disp
23 allocate(logPDF,
mold = Point)
26 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
28 call disp%show(
"! Compute the Probability Density Function (PDF) of GenGamma distribution at the specified values.")
29 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
30 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
34 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
35 call disp%show(
"! Compute the PDF at an input scalar real value.")
36 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
40 call disp%show(
"logPDF(1) = getGenGammaLogPDF(0.5_RK)")
49 call disp%show(
"logPDF(1) = getGenGammaLogPDF(0.5_RK, Kappa(1))")
60 call disp%show(
"logPDF(1) = getGenGammaLogPDF(0.5_RK, Kappa(1), invOmega(1))")
73 call disp%show(
"logPDF(1) = getGenGammaLogPDF(0.5_RK, Kappa(1), invOmega(1), invSigma(1))")
80 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
81 call disp%show(
"! Compute the PDF at an input vector real value with different parameter values.")
82 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
88 call disp%show(
"logPDF(1:NP:NP/5) = getGenGammaLogPDF(0.5_RK, Kappa(1:NP:NP/5))")
97 call disp%show(
"logPDF(1:NP:NP/5) = getGenGammaLogPDF(Point(1:NP:NP/5), Kappa(1:NP:NP/5))")
108 integer :: fileUnit, i
109 real(RK) :: logPDF(NP,
6)
116 open(newunit
= fileUnit, file
= "getGenGammaLogPDF.RK.txt")
117 write(fileUnit,
"(7(g0,:,' '))") (Point(i),
exp(logPDF(i,:)), i
= 1,
size(Point))
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.
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 ⛓
55-1.83669543,
-1.10657454,
-0.569179952,
-0.722454488,
-3.13467169
62+2.03589296,
+1.12029135,
+0.148524448,
-0.738731623,
-1.38531518
Postprocessing of the example output ⛓
3import matplotlib.pyplot
as plt
15xlab = {
"CK" :
"X ( real/imaginary components )"
16 ,
"IK" :
"X ( integer-valued )"
17 ,
"RK" :
"X ( real-valued )"
19legends = [
"$\kappa, 1/\omega, 1/\sigma = 1.0, 0.5, 0.5$"
20 ,
"$\kappa, 1/\omega, 1/\sigma = 2.0, 0.5, 1.0$"
21 ,
"$\kappa, 1/\omega, 1/\sigma = 0.5, 2.0, 0.5$"
22 ,
"$\kappa, 1/\omega, 1/\sigma = 0.2, 5.0, 0.2$"
23 ,
"$\kappa, 1/\omega, 1/\sigma = .14, 7.0, .14$"
24 ,
"$\kappa, 1/\omega, 1/\sigma = 2.0, 5.0, 0.3$"
27for kind
in [
"IK",
"CK",
"RK"]:
29 pattern =
"*." + kind +
".txt"
30 fileList = glob.glob(pattern)
31 if len(fileList) == 1:
33 df = pd.read_csv(fileList[0], delimiter =
" ")
35 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
39 plt.plot( df.values[:, 0]
40 , df.values[:,1:len(legends)+1]
44 plt.plot( df.values[:, 1]
45 , df.values[:,1:len(legends)+1]
50 plt.plot( df.values[:, 0]
51 , df.values[:,1:len(legends)+1]
60 plt.xticks(fontsize = fontsize - 2)
61 plt.yticks(fontsize = fontsize - 2)
62 ax.set_xlabel(xlab[kind], fontsize = 17)
63 ax.set_ylabel(
"Probability Density Function (PDF)", fontsize = 17)
67 plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
68 ax.tick_params(axis =
"y", which =
"minor")
69 ax.tick_params(axis =
"x", which =
"minor")
71 plt.savefig(fileList[0].replace(
".txt",
".png"))
73 elif len(fileList) > 1:
75 sys.exit(
"Ambiguous file list exists.")
Visualization of the example output ⛓
- Test:
- test_pm_distGenGamma
- 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 491 of file pm_distGenGamma.F90.