hfncdf(3f) - [M_datapac:CUMULATIVE_DISTRIBUTION] compute the half-normal cumulative
distribution function
SUBROUTINE HFNCDF(X,Cdf)
REAL(kind=wp),intent(in) :: X
REAL(kind=wp),intent(out) :: Cdf
HFNCDF(3f) computes the cumulative distribution function value for
the halfnormal distribution.
The halfnormal distribution used herein has mean = sqrt(2/pi) =
0.79788456 and standard deviation = 1.
This distribution is defined for all non-negative x and has the
probability density function
f(X) = (2/sqrt(2*pi)) * exp(-X*X/2).
The halfnormal distribution used herein is the distribution of the
variate X = abs(z) where the variate z is normally distributed with
mean = 0 and standard deviation = 1.
X The value at which the cumulative distribution function is
to be evaluated. X should be non-negative.
CDF The cumulative distribution function value.
for the halfnormal distribution
Sample program:
program demo_hfncdf
!@(#) line plotter graph of cumulative distribution function
!@(#) for the halfnormal distribution
use M_datapac, only : hfncdf, plott, label
implicit none
real,allocatable :: x(:), y(:)
integer :: i
call label('hfncdf')
x=[(real(i),i=0,100,1)]
if(allocated(y))deallocate(y)
allocate(y(size(x)))
do i=1,size(x)
call hfncdf(x(i)/10.0,y(i))
enddo
call plott(x,y,size(x))
end program demo_hfncdf
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 X
0.4166667E+02 I X
0.3750000E+02 I X
0.3333334E+02 I X
0.2916667E+02 I X
0.2500000E+02 - XX
0.2083334E+02 I XXX
0.1666667E+02 I XXXX
0.1250000E+02 I X X XX
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.1192E-06 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) | :: | Cdf |
subroutine hfncdf(X,Cdf) real(kind=wp),intent(in) :: X real(kind=wp),intent(out) :: Cdf ! ! 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 HFNCDF(3f) is negative *****') write (G_io,99002) X 99002 format (' ***** The value of the argument is ',E15.8,' *****') Cdf = 0.0_wp return else call norcdf(X,Cdf) Cdf = 2.0_wp*Cdf - 1.0_wp endif end subroutine hfncdf