unicdf(3f) - [M_datapac:CUMULATIVE_DISTRIBUTION] trivially compute the Uniform
cumulative distribution function
subroutine unicdf(X,Cdf)
real(kind=wp), intent(in) :: X
real(kind=wp), intent(out) :: Cdf
UNICDF(3f) computes the cumulative distribution function value for
the uniform (rectangular) distribution on the unit interval (0,1).
This distribution has mean = 0.5 and standard deviation = sqrt(1/12)
= 0.28867513.
This distribution has the probability density function f(x) = x.
That is, this is a trivial function as the output equals the input.
X The value at which the cumulative distribution function is to
be evaluated. X should be between 0 and 1, inclusively.
CDF the REAL cumulative distribution function value.
Sample program:
program demo_unicdf
!@(#) line plotter graph of function
use M_datapac, only : unicdf, plott, label
implicit none
integer,parameter :: n=40
real :: x(0:n), y(0:n)
integer :: i
call label('unicdf')
x=[(real(i)/real(n),i=0,n)]
do i=0,n
call unicdf(x(i),y(i))
enddo
call plott(x,y,n+1)
end program demo_unicdf
Result:
The following is a plot of Y(I) (vertically) versus X(I) (horizontally) I-----------I-----------I-----------I-----------I 0.1000000E+01 - X 0.9583333E+00 I XX 0.9166667E+00 I XX 0.8750000E+00 I X 0.8333333E+00 I XX 0.7916666E+00 I XX 0.7500000E+00 - X 0.7083333E+00 I XX 0.6666666E+00 I XX 0.6250000E+00 I X 0.5833333E+00 I XX 0.5416666E+00 I XX 0.5000000E+00 - X 0.4583333E+00 I XX 0.4166666E+00 I XX 0.3750000E+00 I X 0.3333333E+00 I XX 0.2916666E+00 I XX 0.2500000E+00 - X 0.2083333E+00 I XX 0.1666666E+00 I XX 0.1250000E+00 I X 0.8333331E-01 I XX 0.4166663E-01 I XX 0.0000000E+00 - X I-----------I-----------I-----------I-----------I 0.0000E+00 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 unicdf(X,Cdf) real(kind=wp), intent(in) :: X real(kind=wp), intent(out) :: Cdf if ( X < 0.0_wp .or. X > 1.0_wp ) then ! CHECK THE INPUT ARGUMENTS FOR ERRORS write (g_io,99001) 99001 format (' ***** NON-FATAL DIAGNOSTIC--The first input argument to UNICDF(3f) is outside the usual (0,1) interval *****') write (g_io,99002) X 99002 format (' ***** The value of the argument is ',E15.8,' *****') if ( X < 0.0_wp ) then Cdf = 0.0_wp else if ( X > 1.0_wp ) then Cdf = 1.0_wp else stop '**unicdf** should not get here' endif else Cdf = X endif end subroutine unicdf