exppdf(3f) - [M_datapac:PROBABILITY_DENSITY] compute the exponential probability
density function
SUBROUTINE EXPPDF(X,Pdf)
REAL(kind=wp),intent(in) :: X
REAL(kind=wp),intent(out) :: Pdf
EXPPDF(3f) computes the probability density function value for the
exponential distribution with mean = 1 and standard deviation = 1.
This distribution is defined for all non-negative X, and has the
probability density function
f(X) = exp(-X)
X The value at which the probability density
function is to be evaluated. Values should be non-negative.
PDF The probability density function value.
Sample program:
program demo_exppdf
!@(#) line plotter graph of probability density function
use M_datapac, only : exppdf, plott, label
implicit none
real,allocatable :: x(:), y(:)
integer :: i
call label('exppdf')
x=[(real(i),i=0,100,1)]
if(allocated(y))deallocate(y)
allocate(y(size(x)))
do i=1,size(x)
call exppdf(x(i)/10.0,y(i))
enddo
call plott(x,y,size(x))
end program demo_exppdf
``` Results:
The following is a plot of Y(I) (vertically) versus X(I) (horizontally)
I-----------I-----------I-----------I-----------I
0.1000000E+03 - X
0.9583334E+02 I X
0.9166666E+02 I X
0.8750000E+02 I X
0.8333334E+02 I X
0.7916667E+02 I X
0.7500000E+02 - X
0.7083334E+02 I X
0.6666667E+02 I X
0.6250000E+02 I X
0.5833334E+02 I X
0.5416667E+02 I X
0.5000000E+02 - X
0.4583334E+02 I XX
0.4166667E+02 I X
0.3750000E+02 I X
0.3333334E+02 I XX
0.2916667E+02 I XX
0.2500000E+02 - XXX
0.2083334E+02 I XXX
0.1666667E+02 I XXXX
0.1250000E+02 I XXX X
0.8333336E+01 I X X X X
0.4166672E+01 I X X X X
0.0000000E+00 - X X X
I-----------I-----------I-----------I-----------I
0.4540E-04 0.2500E+00 0.5000E+00 0.7500E+00 0.1000E+01
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 EXPPDF(X,Pdf) REAL(kind=wp),intent(in) :: X REAL(kind=wp),intent(out) :: Pdf ! CHECK THE INPUT ARGUMENTS FOR ERRORS ! IF ( X<0.0_wp ) THEN WRITE (G_IO,99001) 99001 FORMAT (' ***** NON-FATAL DIAGNOSTIC--THE FIRST INPUT ARGUMENT TO EXPPDF(3f) IS NEGATIVE *****') WRITE (G_IO,99002) X 99002 FORMAT (' ***** THE VALUE OF THE ARGUMENT IS ',E15.8,' *****') Pdf = 0.0_wp RETURN ELSE Pdf = EXP(-X) ENDIF ! END SUBROUTINE EXPPDF