Generate and return the natural logarithm of the Probability Density Function (PDF) of the Exponential distribution for an input x
within the support of the distribution \([\mu, +\infty)\).
More...
Generate and return the natural logarithm of the Probability Density Function (PDF) of the Exponential distribution for an input x
within the support of the distribution \([\mu, +\infty)\).
- Parameters
-
[in] | x | : The input scalar (or array of the same rank, shape, and size 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 PDF must be computed.
|
[in] | mu | : The input scalar (or array of the same rank, shape, and size as other array-like arguments), of the same type and kind as x containing the location parameter of the distribution.
(optional, default = 0 ) |
[in] | invSigma | : The input scalar (or array of the same rank, shape, and size as other array-like arguments), of the same type and kind as x containing the inverse scale (i.e., the rate or the inverse mean) parameter of the distribution.
(optional, default = 1. ) |
- Returns
logPDF
: The output scalar (or array of the same rank, shape, and size as other array-like arguments), of the same type and kind as x
, containing the natural logarithm of the PDF of the distribution at the specified point.
Possible calling interfaces ⛓
Generate and return the natural logarithm of the Probability Density Function (PDF) of the Exponentia...
This module contains classes and procedures for computing various statistical quantities related to t...
- Warning
- The condition
mu <= x
must hold for the corresponding input arguments.
The condition 0 < invSigma
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
- setExpLogPDF
Example usage ⛓
10 type(display_type) :: disp
14 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
15 call disp%show(
"! Compute the Probability Density Function (PDF) of Exponential distribution at the specified values.")
16 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
20 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
21 call disp%show(
"! Compute the PDF at an input scalar real value.")
22 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
26 call disp%show(
"logPDF(1) = getExpLogPDF(x = +.5)")
33 call disp%show(
"logPDF(1) = getExpLogPDF(x = +.5, mu = 0.)")
40 call disp%show(
"logPDF(1) = getExpLogPDF(x = +.5, invSigma = 1.)")
47 call disp%show(
"logPDF(1) = getExpLogPDF(x = +.5, mu = 0., invSigma = 1.)")
48 logPDF(
1)
= getExpLogPDF(x
= +.
5, mu
= 0., invSigma
= 1.)
54 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
55 call disp%show(
"! Compute the PDF at an input vector real value with different parameter values.")
56 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
60 call disp%show(
"logPDF(1:3) = getExpLogPDF(x = +[.5, 1., 2.])")
67 call disp%show(
"logPDF(1:3) = getExpLogPDF(x = +[.5, 1., 2.], mu = 0.)")
74 call disp%show(
"logPDF(1:3) = getExpLogPDF(x = +[.5, 1., 2.], invSigma = [2., 1., .5])")
75 logPDF(
1:
3)
= getExpLogPDF(x
= +[.
5,
1.,
2.], invSigma
= [
2.,
1., .
5])
81 call disp%show(
"logPDF(1:3) = getExpLogPDF(x = +[.5, 1., 2.], mu = [-1., 0., +1.], invSigma = 1.)")
82 logPDF(
1:
3)
= getExpLogPDF(x
= +[.
5,
1.,
2.], mu
= [
-1.,
0.,
+1.], invSigma
= 1.)
93 real,
parameter :: invSigma(
3)
= [
2.,
1., .
5]
94 real,
parameter :: mu(
3)
= [
-2.,
0.,
2.]
96 integer :: fileUnit, i, j
98 open(newunit
= fileUnit, file
= "getExpLogPDF.RK.txt")
100 do j
= 1,
size(logPDF)
101 logPDF(j)
= -huge(
0.)
102 if (mu(j)
<= Point(i)) logPDF(j)
= getExpLogPDF(Point(i), mu(j), invSigma(j))
104 write(fileUnit,
"(*(g0,:,', '))") Point(i),
exp(logPDF)
Generate count evenly spaced points over the interval [x1, x2] if x1 < x2, or [x2,...
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 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...
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 ⛓
39-0.500000000,
-1.00000000,
-2.00000000
44-0.500000000,
-1.00000000,
-2.00000000
47logPDF(
1:
3)
= getExpLogPDF(x
= +[.
5,
1.,
2.], invSigma
= [
2.,
1., .
5])
49-0.306852818,
-1.00000000,
-1.69314718
52logPDF(
1:
3)
= getExpLogPDF(x
= +[.
5,
1.,
2.], mu
= [
-1.,
0.,
+1.], invSigma
= 1.)
54-1.50000000,
-1.00000000,
-1.00000000
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 = [
r"$\mu = -2., invSigma = +2.$"
20 ,
r"$\mu = +0., invSigma = +1.$"
21 ,
r"$\mu = +2., invSigma = +.5$"
24for kind
in [
"IK",
"CK",
"RK"]:
26 pattern =
"*." + kind +
".txt"
27 fileList = glob.glob(pattern)
28 if len(fileList) == 1:
30 df = pd.read_csv(fileList[0], delimiter =
", ")
32 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
36 plt.plot( df.values[:, 0]
41 plt.plot( df.values[:, 1]
46 ax.legend ( [
"real",
"imaginary"]
50 plt.plot( df.values[:, 0]
55 , fontsize = fontsize - 5
58 plt.xticks(fontsize = fontsize - 2)
59 plt.yticks(fontsize = fontsize - 2)
60 ax.set_xlabel(xlab[kind], fontsize = 17)
61 ax.set_ylabel(
"Probability Density Function (PDF)", fontsize = 17)
63 plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
64 ax.tick_params(axis =
"y", which =
"minor")
65 ax.tick_params(axis =
"x", which =
"minor")
67 plt.savefig(fileList[0].replace(
".txt",
".png"))
69 elif len(fileList) > 1:
71 sys.exit(
"Ambiguous file list exists.")
Visualization of the example output ⛓
- Test:
- test_pm_distExp
- 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 221 of file pm_distExp.F90.