dexpdf(3f) - [M_datapac:PROBABILITY_DENSITY] compute the double
exponential probability density function
SUBROUTINE DEXPDF(X,Pdf)
REAL(kind=wp),intent(in) :: X
REAL(kind=wp),intent(out) :: Pdf
DEXPDF(3f) computes the probability density function value for the
double exponential (Laplace) distribution with mean = 0 and standard
deviation = sqrt(2).
This distribution is defined for all X and has the probability
density function
f(X) = 0.5*exp(-abs(X))
X The value at which the probability density function is to
be evaluated.
PDF The probability density function value.
Sample program:
program demo_dexpdf
!@(#) line plotter graph
!@(#) of probability density function for Laplace distribution
use M_datapac, only : dexpdf, plott, label
implicit none
real,allocatable :: x(:), y(:)
integer :: i
call label('dexpdf')
x=[(real(i),i=-100,100,1)]
if(allocated(y))deallocate(y)
allocate(y(size(x)))
do i=1,size(x)
call dexpdf(x(i)/10.0,y(i))
enddo
call plott(x,y,size(x))
end program demo_dexpdf
Results:
The original DATAPAC library was written by James Filliben of the Statistical
Engineering Division, National Institute of Standards and Technology.
John Urban, 2022.05.31
CC0-1.0
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in) | :: | X | |||
real(kind=wp), | intent(out) | :: |
subroutine dexpdf(X,Pdf) real(kind=wp),intent(in) :: X real(kind=wp),intent(out) :: Pdf real(kind=wp) :: arg ! ! CHECK THE INPUT ARGUMENTS FOR ERRORS. -- NO INPUT ARGUMENT ERRORS POSSIBLE FOR THIS DISTRIBUTION. ! arg = X IF ( X<0.0_wp ) arg = -X Pdf = 0.5_wp*EXP(-arg) end subroutine dexpdf