caupdf(3f) - [M_datapac:PROBABILITY_DENSITY] compute the Cauchy probability
density function
subroutine caupdf(X,Pdf)
real(kind=wp),intent(in) :: X
real(kind=wp),intent(out):: Pdf
CAUPDF(3f) computes the probability density function value for the
Cauchy distribution with median = 0 and 75% point = 1.
This distribution is defined for all X and has the probability
density function
f(x) = (1/pi)*(1/(1+x*x))
X The value at which the probability density function is to be
evaluated.
PDF The probability density function value.
Sample program:
program demo_caupdf
!@(#) line plotter graph of probability density function
use M_datapac, only : caupdf, plott, label
implicit none
real,allocatable :: x(:), y(:)
integer :: i
call label('caupdf')
x=[(real(i),i=-100,100,1)]
if(allocated(y))deallocate(y)
allocate(y(size(x)))
do i=1,size(x)
call caupdf(x(i)/10.0,y(i))
enddo
call plott(x,y,size(x))
end program demo_caupdf
Results:
The following is a plot of Y(i) (vertically) versus X(i) (horizontally)
I-----------I-----------I-----------I-----------I
0.1000000E+03 - X
0.9166666E+02 I X
0.8333334E+02 I X
0.7500000E+02 I X
0.6666667E+02 I XX
0.5833334E+02 I X
0.5000000E+02 - XX
0.4166667E+02 I XX
0.3333334E+02 I XX
0.2500000E+02 I XXXX
0.1666667E+02 I XXXXXX X X
0.8333336E+01 I X X X X X X X X
0.0000000E+00 - X X X X
-0.8333328E+01 I X X X X X X X X
-0.1666666E+02 I XXXXXX X X
-0.2499999E+02 I XXXX
-0.3333333E+02 I XX
-0.4166666E+02 I XX
-0.5000000E+02 - XX
-0.5833333E+02 I X
-0.6666666E+02 I XX
-0.7500000E+02 I X
-0.8333333E+02 I X
-0.9166666E+02 I X
-0.1000000E+03 - X
I-----------I-----------I-----------I-----------I
0.3152E-02 0.8194E-01 0.1607E+00 0.2395E+00 0.3183E+00
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 caupdf(X,Pdf) real(kind=wp),intent(in) :: X real(kind=wp),intent(out):: Pdf real(kind=wp),parameter :: c = 0.31830988618379_wp ! ! CHECK THE INPUT ARGUMENTS FOR ERRORS -- NO INPUT ARGUMENT ERRORS POSSIBLE FOR THIS DISTRIBUTION. ! Pdf = c*(1.0_wp/(1.0_wp+X*X)) end subroutine caupdf