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

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...

Detailed Description

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

logPDF = getExpLogPDF(x, mu = mu, invSigma = invSigma)
Generate and return the natural logarithm of the Probability Density Function (PDF) of the Exponentia...
Definition: pm_distExp.F90:221
This module contains classes and procedures for computing various statistical quantities related to t...
Definition: pm_distExp.F90:112
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.
Remarks
The procedures under discussion are elemental.
See also
setExpLogPDF


Example usage

1program example
2
3 use pm_kind, only: SK, IK
4 use pm_io, only: display_type
5 use pm_distExp, only: getExpLogPDF
6
7 implicit none
8
9 real :: logPDF(3)
10 type(display_type) :: disp
11 disp = display_type(file = "main.out.F90")
12
13 call disp%skip()
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("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
17 call disp%skip()
18
19 call disp%skip()
20 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
21 call disp%show("! Compute the PDF at an input scalar real value.")
22 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
23 call disp%skip()
24
25 call disp%skip()
26 call disp%show("logPDF(1) = getExpLogPDF(x = +.5)")
27 logPDF(1) = getExpLogPDF(x = +.5)
28 call disp%show("logPDF(1)")
29 call disp%show( logPDF(1) )
30 call disp%skip()
31
32 call disp%skip()
33 call disp%show("logPDF(1) = getExpLogPDF(x = +.5, mu = 0.)")
34 logPDF(1) = getExpLogPDF(x = +.5, mu = 0.)
35 call disp%show("logPDF(1)")
36 call disp%show( logPDF(1) )
37 call disp%skip()
38
39 call disp%skip()
40 call disp%show("logPDF(1) = getExpLogPDF(x = +.5, invSigma = 1.)")
41 logPDF(1) = getExpLogPDF(x = +.5, invSigma = 1.)
42 call disp%show("logPDF(1)")
43 call disp%show( logPDF(1) )
44 call disp%skip()
45
46 call disp%skip()
47 call disp%show("logPDF(1) = getExpLogPDF(x = +.5, mu = 0., invSigma = 1.)")
48 logPDF(1) = getExpLogPDF(x = +.5, mu = 0., invSigma = 1.)
49 call disp%show("logPDF(1)")
50 call disp%show( logPDF(1) )
51 call disp%skip()
52
53 call disp%skip()
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("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
57 call disp%skip()
58
59 call disp%skip()
60 call disp%show("logPDF(1:3) = getExpLogPDF(x = +[.5, 1., 2.])")
61 logPDF(1:3) = getExpLogPDF(x = +[.5, 1., 2.])
62 call disp%show("logPDF(1:3)")
63 call disp%show( logPDF(1:3) )
64 call disp%skip()
65
66 call disp%skip()
67 call disp%show("logPDF(1:3) = getExpLogPDF(x = +[.5, 1., 2.], mu = 0.)")
68 logPDF(1:3) = getExpLogPDF(x = +[.5, 1., 2.], mu = 0.)
69 call disp%show("logPDF(1:3)")
70 call disp%show( logPDF(1:3) )
71 call disp%skip()
72
73 call disp%skip()
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])
76 call disp%show("logPDF(1:3)")
77 call disp%show( logPDF(1:3) )
78 call disp%skip()
79
80 call disp%skip()
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.)
83 call disp%show("logPDF(1:3)")
84 call disp%show( logPDF(1:3) )
85 call disp%skip()
86
87 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
88 ! Output an example array for visualization.
89 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
90
91 block
92 use pm_arraySpace, only: getLinSpace
93 real, parameter :: invSigma(3) = [2., 1., .5]
94 real, parameter :: mu(3) = [-2., 0., 2.]
95 real :: Point(1000)
96 integer :: fileUnit, i, j
97 Point = getLinSpace(-4., +8., size(Point,1,IK))
98 open(newunit = fileUnit, file = "getExpLogPDF.RK.txt")
99 do i = 1, size(Point)
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))
103 end do
104 write(fileUnit,"(*(g0,:,', '))") Point(i), exp(logPDF)
105 end do
106 close(fileUnit)
107 end block
108
109end program example
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.
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 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
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! Compute the Probability Density Function (PDF) of Exponential distribution at the specified values.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
8! Compute the PDF at an input scalar real value.
9!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
10
11
12logPDF(1) = getExpLogPDF(x = +.5)
13logPDF(1)
14-0.500000000
15
16
17logPDF(1) = getExpLogPDF(x = +.5, mu = 0.)
18logPDF(1)
19-0.500000000
20
21
22logPDF(1) = getExpLogPDF(x = +.5, invSigma = 1.)
23logPDF(1)
24-0.500000000
25
26
27logPDF(1) = getExpLogPDF(x = +.5, mu = 0., invSigma = 1.)
28logPDF(1)
29-0.500000000
30
31
32!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
33! Compute the PDF at an input vector real value with different parameter values.
34!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
35
36
37logPDF(1:3) = getExpLogPDF(x = +[.5, 1., 2.])
38logPDF(1:3)
39-0.500000000, -1.00000000, -2.00000000
40
41
42logPDF(1:3) = getExpLogPDF(x = +[.5, 1., 2.], mu = 0.)
43logPDF(1:3)
44-0.500000000, -1.00000000, -2.00000000
45
46
47logPDF(1:3) = getExpLogPDF(x = +[.5, 1., 2.], invSigma = [2., 1., .5])
48logPDF(1:3)
49-0.306852818, -1.00000000, -1.69314718
50
51
52logPDF(1:3) = getExpLogPDF(x = +[.5, 1., 2.], mu = [-1., 0., +1.], invSigma = 1.)
53logPDF(1:3)
54-1.50000000, -1.00000000, -1.00000000
55
56

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 = [ r"$\mu = -2., invSigma = +2.$"
20 , r"$\mu = +0., invSigma = +1.$"
21 , r"$\mu = +2., invSigma = +.5$"
22 ]
23
24for kind in ["IK", "CK", "RK"]:
25
26 pattern = "*." + kind + ".txt"
27 fileList = glob.glob(pattern)
28 if len(fileList) == 1:
29
30 df = pd.read_csv(fileList[0], delimiter = ", ")
31
32 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
33 ax = plt.subplot()
34
35 if kind == "CK":
36 plt.plot( df.values[:, 0]
37 , df.values[:,2]
38 , marker[kind]
39 , color = "r"
40 )
41 plt.plot( df.values[:,1]
42 , df.values[:,3]
43 , marker[kind]
44 , color = "blue"
45 )
46 ax.legend ( ["real", "imaginary"]
47 , fontsize = fontsize
48 )
49 else:
50 plt.plot( df.values[:, 0]
51 , df.values[:,1:4]
52 , marker[kind]
53 )
54 ax.legend ( legends
55 , fontsize = fontsize - 5
56 )
57
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)
62
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")
66
67 plt.savefig(fileList[0].replace(".txt",".png"))
68
69 elif len(fileList) > 1:
70
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.

  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 221 of file pm_distExp.F90.


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