NORPDF Subroutine

public subroutine NORPDF(X, Pdf)

NAME

norpdf(3f) - [M_datapac:PROBABILITY_DENSITY] compute the normal
probability density function

SYNOPSIS

   SUBROUTINE NORPDF(X,Pdf)

    REAL(kind=wp),intent(in)  :: X
    REAL(kind=wp),intent(out) :: Pdf

DESCRIPTION

NORPDF(3f) computes the probability density function value for the
normal (Gaussian) distribution with mean = 0 and standard deviation
= 1.

This distribution is defined for all X and has the probability
density function

    f(X) = (1/sqrt(2*pi))*exp(-X*X/2)

INPUT ARGUMENTS

X      The value at which the probability density function is to
       be evaluated.

OUTPUT ARGUMENTS

PDF    The probability density function value.

EXAMPLES

Sample program:

program demo_norpdf
!@(#) line plotter graph of probability density function
use M_datapac, only : norpdf, plott, label
implicit none
real,allocatable  :: x(:), y(:)
integer           :: i
   call label('norpdf')
   x=[(real(i),i=-100,100,1)]
   if(allocated(y))deallocate(y)
   allocate(y(size(x)))
   do i=1,size(x)
      call norpdf(x(i)/10.0,y(i))
   enddo
   call plott(x,y,size(x))
end program demo_norpdf

Results:

AUTHOR

The original DATAPAC library was written by James Filliben of the
Statistical Engineering Division, National Institute of Standards
and Technology.

MAINTAINER

John Urban, 2022.05.31

LICENSE

CC0-1.0

REFERENCES

  • Johnson and Kotz, Continuous Univariate Distributions–1, 1970, pages 40-111.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: X
real(kind=wp), intent(out) :: Pdf

Source Code

SUBROUTINE NORPDF(X,Pdf)
REAL(kind=wp),intent(in)  :: X
REAL(kind=wp),intent(out) :: Pdf

REAL(kind=wp),parameter   :: c=0.3989422804_wp
!
!     CHECK THE INPUT ARGUMENTS FOR ERRORS. -- NO INPUT ARGUMENT ERRORS POSSIBLE FOR THIS DISTRIBUTION.
!
      Pdf = c*EXP(-(X*X)/2.0_wp)
!
END SUBROUTINE NORPDF