JMSLTM Numerical Library 6.0

com.imsl.stat
Class Pdf

java.lang.Object
  extended by com.imsl.stat.Pdf

public final class Pdf
extends Object

Probability density functions.

See Also:
Example, Cumulative distribution function, Inverse Cumulative Distribution Function

Nested Class Summary
static class Pdf.AltSeriesAccuracyLossException
          The magnitude of alternating series sum is too small relative to the sum of positive terms to permit a reliable accuracy.
 
Method Summary
static double beta(double x, double pin, double qin)
          Evaluates the beta probability density function.
static double binomial(int k, int n, double pin)
          Evaluates the binomial probability density function.
static double chi(double chsq, double df)
          Evaluates the chi-squared probability density function
static double discreteUniform(int x, int n)
          Evaluates the discrete uniform probability density function.
static double exponential(double x, double scale)
          Evaluates the exponential probability density function
static double extremeValue(double x, double mu, double beta)
          Evaluates the extreme value probability density function.
static double F(double x, double dfn, double dfd)
          Evaluates the F probability density function.
static double gamma(double x, double a, double b)
          Evaluates the gamma probability density function.
static double geometric(int x, double pin)
          Evaluates the discrete geometric probability density function.
static double hypergeometric(int k, int sampleSize, int defectivesInLot, int lotSize)
          Evaluates the hypergeometric probability density function.
static double logistic(double x, double mu, double s)
          Evaluates the logistic probability density function.
static double logNormal(double x, double mu, double sigma)
          Evaluates the standard lognormal probability density function.
static double noncentralBeta(double x, double shape1, double shape2, double lambda)
          Evaluates the noncentral beta probability density function (Pdf).
static double noncentralChi(double chsq, double df, double alam)
          Evaluates the noncentral chi-squared probability density function.
static double noncentralF(double f, double df1, double df2, double lambda)
          Evaluates the noncentral F probability density function (Pdf).
static double noncentralStudentsT(double t, double df, double delta)
          Evaluates the noncentral Student's t probability density function.
static double normal(double x, double mean, double stdev)
          Evaluates the normal (Gaussian) probability density function.
static double Pareto(double x, double xm, double k)
          Evaluates the Pareto probability density function.
static double poisson(int k, double theta)
          Evaluates the Poisson probability density function.
static double Rayleigh(double x, double alpha)
          Evaluates the Rayleigh probability density function.
static double Weibull(double x, double gamma, double alpha)
          Evaluates the Weibull probability density function.
 
Methods inherited from class java.lang.Object
clone, equals, finalize, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
 

Method Detail

beta

public static double beta(double x,
                          double pin,
                          double qin)
Evaluates the beta probability density function.

Parameters:
x - a double, the argument at which the function is to be evaluated.
pin - a double, the first beta distribution parameter.
qin - a double, the second beta distribution parameter.
Returns:
a double, the value of the probability density function at x.

binomial

public static double binomial(int k,
                              int n,
                              double pin)
Evaluates the binomial probability density function.

Method binomial evaluates the probability that a binomial random variable with parameters n and p with p=pin takes on the value k. It does this by computing probabilities of the random variable taking on the values in its range less than (or the values greater than) k. These probabilities are computed by the recursive relationship

Pr left( {X = j} right) = frac{{left(
  {n + 1 - j} right)p}}{{jleft( {1 - p} right)}}Pr left( {X = j - 1}
  right)

To avoid the possibility of underflow, the probabilities are computed forward from 0, if k is not greater than n times p, and are computed backward from n, otherwise. The smallest positive machine number, varepsilon, is used as the starting value for computing the probabilities, which are rescaled by (1 - p)^n varepsilon if forward computation is performed and by p^nvarepsilon if backward computation is done.

For the special case of p = 0, binomial is set to 0 if k is greater than 0 and to 1 otherwise; and for the case p = 1, binomial is set to 0 if k is less than n and to 1 otherwise.

Plot of Binomial Probability Function

Parameters:
k - the int argument for which the binomial distribution function is to be evaluated.
n - the int number of Bernoulli trials.
pin - a double scalar value representing the probability of success on each independent trial.
Returns:
a double scalar value representing the probability that a binomial random variable takes a value equal to k.
See Also:
Example

chi

public static double chi(double chsq,
                         double df)
Evaluates the chi-squared probability density function

Parameters:
chsq - a double scalar value representing the argument at which the function is to be evaluated.
df - a double scalar value representing the number of degrees of freedom. df must be positive.
Returns:
a double scalar value, the value of the probability density function at chsq.

discreteUniform

public static double discreteUniform(int x,
                                     int n)
Evaluates the discrete uniform probability density function.

Parameters:
x - an int argument for which the discrete uniform probability density function is to be evaluated. x should be a value between the lower limit 0 and upper limit n
n - an int scalar value representing the upper limit of the discrete uniform distribution.
Returns:
a double scalar value representing the probability that a discrete uniform random variable takes a value equal to x.
See Also:
Example

exponential

public static double exponential(double x,
                                 double scale)
Evaluates the exponential probability density function

Parameters:
x - a double scalar value representing the argument at which the function is to be evaluated.
scale - a double scalar value representing the scale parameter.
Returns:
a double scalar value, the value of the probability density function at x.

extremeValue

public static double extremeValue(double x,
                                  double mu,
                                  double beta)
Evaluates the extreme value probability density function.

Parameters:
x - a double scalar value representing the argument at which the function is to be evaluated.
mu - a double scalar value representing the location parameter.
beta - a double scalar value representing the scale parameter.
Returns:
a double scalar value representing the probability density function at x.
See Also:
Example

F

public static double F(double x,
                       double dfn,
                       double dfd)
Evaluates the F probability density function.

The probability density function of the F distribution is

{it f}(x, {it dfn}, {it dfd})=
 {frac { {Gamma}(frac {v_1 + v_2}{2})({frac {v_1}{v_2})}^{frac{v_1}{2}}
 x^{frac {v_1}{2}} } {{{Gamma}(frac {v_1}{2}) }{{Gamma}(frac {v_2}{2}) }
 {(1+frac{v_1x}{v_2})}^{frac{v_1+v_2}{2} }   }}

where v_1 and v_2 are the shape parameters dfn and dfd and Gamma is the gamma function,

{Gamma (a)} = {int}_{0}^{infty}{t^{a-1} e^{-t } it dt}

.

Parameters:
x - a double, the argument at which the function is to be evaluated.
dfn - a double, the numerator degrees of freedom. It must be positive.
dfd - a double, the denominator degrees of freedom. It must be positive.
Returns:
a double, the value of the probability density function at x.

gamma

public static double gamma(double x,
                           double a,
                           double b)
Evaluates the gamma probability density function. The probability density function of the gamma distribution is

f(x; a, b) =
     x^{a - 1} frac{1}{{b^{a} Gamma (a)}}
     e^{ -  {x}/{b}}

where a is the shape parameter and b is the scale parameter.

Parameters:
x - a double scalar value representing the argument at which the function is to be evaluated.
a - a double scalar value representing the shape parameter. This must be positive.
b - a double scalar value representing the scale parameter. This must be positive.
Returns:
a double scalar value, the probability density function at x.

geometric

public static double geometric(int x,
                               double pin)
Evaluates the discrete geometric probability density function.

Method geometric evaluates the geometric distribution for the number of trials before the first success.

Parameters:
x - the int argument for which the geometric probability function is to be evaluated
pin - a double scalar value representing the probability parameter of the geometric distribution (the probability of success for each independent trial)
Returns:
a double scalar value representing the probability that a geometric random variable takes a value equal to x.
See Also:
Example

hypergeometric

public static double hypergeometric(int k,
                                    int sampleSize,
                                    int defectivesInLot,
                                    int lotSize)
Evaluates the hypergeometric probability density function.

Method hypergeometric evaluates the probability density function of a hypergeometric random variable with parameters n, l, and m. The hypergeometric random variable X can be thought of as the number of items of a given type in a random sample of size n that is drawn without replacement from a population of size l containing m items of this type. The probability density function is:

{rm{Pr}}left( {X = k} right) =
  frac{{left( {_k^m } right)left( {_{n - k}^{l - m} } right)}}{{left(
  {_n^l } right)}}{rm{for}} ,,, k = i,;i + 1,,i + 2; ldots ,;min
  left( {n,m} right)

where i = max(0, n - l + m). hypergeometric evaluates the expression using log gamma functions.

Parameters:
k - an int, the argument at which the function is to be evaluated.
sampleSize - an int, the sample size, n.
defectivesInLot - an int, the number of defectives in the lot, m.
lotSize - an int, the lot size, l.
Returns:
a double, the probability that a hypergeometric random variable takes on a value equal to k.
See Also:
Example

logistic

public static double logistic(double x,
                              double mu,
                              double s)
Evaluates the logistic probability density function.

The probability density function of the logistic distribution is

f(x,mu,s)=frac{e^{-(x-mu)/s}}
  {sleft (1+e^{-(x-mu)/s} right )^{2}}

where mu is the location parameter and the scale parameter s gt 0.

Parameters:
x - a double scalar value representing the argument at which the function is to be evaluated.
mu - a double scalar value representing the location parameter.
s - a double scalar value representing the scale parameter.
Returns:
a double scalar value representing the probability density function at x.

logNormal

public static double logNormal(double x,
                               double mu,
                               double sigma)
Evaluates the standard lognormal probability density function.

Fleft( x right) = frac{1}{xsigmasqrt{2pi}}
  {e^{-frac{ {(ln{x}-mu)}^2 }{2{sigma}^2}} }

Parameters:
x - a double scalar value representing the argument at which the function is to be evaluated.
mu - a double scalar value representing the location parameter.
sigma - a double scalar value representing the shape parameter. sigma must be a positive.
Returns:
a double scalar value representing the probability density function at x.

noncentralBeta

public static double noncentralBeta(double x,
                                    double shape1,
                                    double shape2,
                                    double lambda)
Evaluates the noncentral beta probability density function (Pdf).

The noncentral beta distribution is a generalization of the beta distribution. If Z is a noncentral chi-square random variable with noncentrality parameter lambda and 2 alpha_1 degrees of freedom, and Y is a chi-square random variable with 2 alpha_2 degrees of freedom which is statistically independent of Z, then

X ;; = ;; frac{Z}{Z ; + ; Y} ;; = ;; frac{alpha_1 F}{alpha_1 F ; + ; alpha_2}

is a noncentral beta-distributed random variable and

F ;; = ;; frac{alpha_2 Z}{alpha_1 Y} ;; = ;; frac{alpha_2 X}{alpha_1 (1  ; -  ; X)}

is a noncentral F-distributed random variable. The Pdf for noncentral beta variable X can thus be simply defined in terms of the noncentral F Pdf:

PDF_{ncbeta}(x, ;  alpha_1, ; alpha_2, ; lambda) ;; = ;;
  PDF_{ncF}(f, ; 2 alpha_1, ; 2 alpha_2, ; lambda) ; frac{df}{dx}

where PDF_{ncbeta}(x, ;  alpha_1, ; alpha_2, ; lambda) is the noncentral beta Pdf with x = x, alpha_1 = shape1, alpha_2 = shape2, and noncentrality parameter lambda = lambda; PDF_{ncF}(f, ; 2 alpha_1, ; 2 alpha_2, ; lambda) is the noncentral F Pdf with argument f, numerator and denominator degrees of freedom 2 alpha_1 and 2 alpha_2 respectively, and noncentrality parameter lambda; and:

f ;; = ;; frac{alpha_2 x}{alpha_1 (1  ; -  ; x)}; ;;
  x ;; = ;; frac{alpha_1 f}{alpha_1 f ; + ; alpha_2}; ;;
  frac{df}{dx} ;; = ;; frac{(alpha_1 f ; + ; alpha_2)^2}{alpha_1 alpha_2} ;; = ;; frac{1}{alpha_1 (1  ; -  ; x)^2}

(See documentation for class Cdf method noncentralF for a discussion of how the noncentral F Pdf is defined and calculated.)

With a noncentrality parameter of zero, the noncentral beta distribution is the same as the beta distribution.

Parameters:
x - a double scalar value representing the argument at which the function is to be evaluated. x must be nonnegative and less than or equal to 1.
shape1 - a double scalar value representing the first shape parameter. shape1 must be positive.
shape2 - a double scalar value representing the second shape parameter. shape2 must be positive.
lambda - a double scalar value representing the noncentrality parameter. lambda must nonnegative.
Returns:
a double scalar value representing the probability density associated with a noncentral beta random variable with value x.

noncentralChi

public static double noncentralChi(double chsq,
                                   double df,
                                   double alam)
Evaluates the noncentral chi-squared probability density function.

The noncentral chi-squared distribution is a generalization of the chi-squared distribution. If {X_i} are k independent, normally distributed random variables with means mu_i and variances sigma^2_i, then the random variable

X ;; = ;; sum_{i = 1}^k left(frac{X_i}{sigma_i}right)^2

is distributed according to the noncentral chi-squared distribution. The noncentral chi-squared distribution has two parameters, k which specifies the number of degrees of freedom (i.e. the number of X_i), and lambda which is related to the mean of the random variables X_i by

lambda ;; = ;; sum_{i = 1}^k left(frac{mu_i}{sigma_i}right)^2

The noncentral chi-squared distribution is equivalent to a (central) chi-squared distribution with k + 2i degrees of freedom, where i is the value of a Poisson distributed random variable with parameter lambda/2. Thus, the probability density function is given by:

F(x,k,lambda) ;; = ;; sum_{i = 0}^infty {frac{e^{-lambda/2} (lambda/2)^i}{i!}} f(x,k+2i)

where the (central) chi-squared Pdf f(x, k) is given by:

f(x, k) ;; = ;;  frac{(x/2)^{k/2} ; e^{-x/2}}{x ; Gamma(k/2)} quad for ;; x ; > ; 0, ;; else ;; 0

where Gamma (cdot) is the gamma function. The above representation of F(x,k,lambda) can be shown to be equivalent to the representation:

F(x,k,lambda) ;; = ;; frac{e^{-(lambda+x)/2} ; (x/2)^{k/2}}{x} ; sum_{i = 0}^infty {phi_i}

phi_i ;; = ;; frac{(lambda x / 4)^i}{i! ;  Gamma(k/2 ;; + ;; i)}

Method noncentralChi evaluates the probability density function, F(x,k,lambda), of a noncentral chi-squared random variable with df degrees of freedom and noncentrality parameter alam, corresponding to k = df, lambda = alam, and x = chsq.

Method noncentralChi evaluates the cumulative distribution function incorporating the above probability density function.

With a noncentrality parameter of zero, the noncentral chi-squared distribution is the same as the central chi-squared distribution.

Parameters:
chsq - a double scalar value at which the function is to be evaluated. chsq must be nonnegative.
df - a double scalar value representing the number of degrees of freedom. df must be positive.
alam - a double scalar value representing the noncentrality parameter. alam must be nonnegative.
Returns:
a double scalar value representing the probability density associated with a noncentral chi-squared random variable with value chsq.

noncentralF

public static double noncentralF(double f,
                                 double df1,
                                 double df2,
                                 double lambda)
Evaluates the noncentral F probability density function (Pdf).

The noncentral F distribution is a generalization of the F distribution. If x is a noncentral chi-square random variable with noncentrality parameter lambda and nu_1 degrees of freedom, and y is a chi-square random variable with nu_2 degrees of freedom which is statistically independent of X, then

F ;; = ;; (x/nu_1)/(y/nu_2)

is a noncentral F-distributed random variable whose Pdf is given by:

Pdf(f, nu_1, nu_2, lambda) ;; = ;; Psi ; sum_{k = 0}^infty {Phi_k}

where

Psi ;; = ;; frac{ e^{-lambda/2}(nu_1 f)^{nu_1/2}(nu_2)^{nu_2/2} }
  { f ; (nu_1 f ; + ; nu_2)^{(nu_1 + nu_2)/2} ; Gamma(nu_2/2) }

Phi_k ;; = ;; frac{ R^k ; Gamma(frac{nu_1 + nu_2}{2} ; + ; k) }
  { k! ; Gamma(frac{nu_1}{2} ; + ; k) }

R ;; = ;; frac{ lambda nu_1 f }{ 2 (nu_1 f ; + ; nu_2)}

where Gamma (cdot) is the gamma function, nu_1 = df1, nu_2 = df2, lambda = lambda, and f = f.

With a noncentrality parameter of zero, the noncentral F distribution is the same as the F distribution.

The efficiency of the calculation of the above series is enhanced by:

  1. calculating each term Phi_k in the series recursively in terms of either the term Phi_{k-1} preceding it or the term Phi_{k+1} following it, and
  2. initializing the sum with the largest series term and adding the subsequent terms in order of decreasing magnitude.

Special cases:

For R ;; = ;; lambda f ;; = ;; 0:

Pdf(f, nu_1, nu_2, lambda) ;; = ;; Psi ; Phi_0 ;; = ;; Psi ;
  frac{ Gamma([nu_1 + nu_2]/2) }{ Gamma(nu_1/2) }

For lambda ;; = ;; 0:

Pdf(f, nu_1, nu_2, lambda) ;; = ;;
  frac{ (nu_1 f)^{nu_1/2} ; (nu_2)^{nu_2/2} ; Gamma([nu_1 + nu_2]/2) }
  { f ; (nu_1 f ; + ; nu_2)^{(nu_1 + nu_2)/2} ; Gamma(nu_1/2) ; Gamma(nu_2/2) }

For f ;; = ;; 0:

Pdf(f, nu_1, nu_2, lambda) ;; = ;; frac{{e^{ - lambda /2} ;f^{nu _1 /2;; - ;;1} ;(nu _1 /nu _2 )^{nu _1 /2} ;Gamma ([nu _1 ; + ;nu _2 ]/2)}} {{;Gamma (nu _1 /2);Gamma (nu _2 /2)}}

Pdf(f, nu_1, nu_2, lambda) ;; = left{ begin{array}{ll}
         0 ,,,,,,,,,,  mbox{if} ;; nu_1 ;; > ;; 2; \
          e^{-lambda/2} ;;    mbox{if} ;; nu_1 ;; = ;; 2; \
    infty ,,,,,,,,,,  mbox{if} ;; nu_1 ;; lt ;; 2 
      end{array}  right.

Parameters:
f - a double value representing the argument at which the function is to be evaluated. f must be nonnegative.
df1 - a double value representing the number of numerator degrees of freedom. df1 must be positive.
df2 - a double value representing the number of denominator degrees of freedom. df2 must be positive.
lambda - a double value representing the noncentrality parameter. lambda must be nonnegative.
Returns:
a double value representing the probability density associated with a noncentral F random variable with value f.

noncentralStudentsT

public static double noncentralStudentsT(double t,
                                         double df,
                                         double delta)
                                  throws Pdf.AltSeriesAccuracyLossException
Evaluates the noncentral Student's t probability density function.

The noncentral Student's t-distribution is a generalization of the Student's t-distribution. If w is a normally distributed random variable with unit variance and mean delta and u is a chi-square random variable with nu degrees of freedom that is statistically independent of w, then

T ;; = ;; w/sqrt{u/nu}

is a noncentral t-distributed random variable with nu degrees of freedom and noncentrality parameter delta, that is, with nu = df, and delta = delta. The probability density function for the noncentral t-distribution is:

f(t,nu,delta) ;; = ;;
  frac{nu^{nu/2} ; e^{-delta^2/2}}{sqrt{pi} ; Gamma(nu/2) ; ( nu + t^2 ) ^ {(nu + 1)/2}}
  ; sum_{i = 0}^infty {Phi_i}

where

Phi_i ;; = ;;
  frac{Gamma((nu + i + 1)/2)}{i!} ; [delta t]^i ; left(frac{2}{nu + t^2}right)^{i/2}

and t = t.

For noncentrality parameter delta = 0, the Pdf reduces to the (central) Student's t Pdf:

f(t,nu,0) ;; = ;;
  frac{Gamma((nu+1)/2) ; left( 1 ; + ; (t^2/nu) right)^{-(nu+1)/2}}{sqrt{nu pi} ; Gamma(nu/2)}

and, for t = 0, the Pdf becomes:

f(0,nu,delta) ;; = ;;
  frac{Gamma((nu+1)/2) ; e^{-delta^2/2}}{sqrt{nu pi} ; Gamma(nu/2)}

Method noncentralStudentsT evaluates the cumulative distribution function incorporating the above probability density function.

Parameters:
t - a double value representing the argument at which the function is to be evaluated.
df - a double value representing the number of degrees of freedom. df must be positive.
delta - a double value representing the noncentrality parameter.
Returns:
a double value representing the probability density associated with a noncentral Student's t random variable with value t.
Throws:
Pdf.AltSeriesAccuracyLossException - is thrown when the magnitude of alternating series sum is too small relative to the sum of positive terms to permit a reliable accuracy.

normal

public static double normal(double x,
                            double mean,
                            double stdev)
Evaluates the normal (Gaussian) probability density function.

The probability density function for a normal distribution is given by

frac{1}{sigma sqrt{2pi}} {e}^{
  frac{{-(x - mu)}^2}{{2 {sigma}^2}} }

where mu and sigma are the conditional mean and standard deviation.

Parameters:
x - a double scalar value representing the argument at which the function is to be evaluated.
mean - a double scalar value containing the mean.
stdev - a double scalar value containing the standard deviation.
Returns:
a double containing the value of the probability density function at x

Pareto

public static double Pareto(double x,
                            double xm,
                            double k)
Evaluates the Pareto probability density function.

The probability density function of the Pareto distribution is

f(x,x_m,k)=1-frac{kx_m^{k}}{x^{k+1}}

where the scale parameter x_m>0 and the shape parameter k>0. The function is only defined for x geq x_m

Parameters:
x - a double scalar value representing the argument at which the function is to be evaluated.
xm - a double scalar value representing the scale parameter, x_m.
k - a double scalar value representing the shape parameter.
Returns:
a double scalar value representing the probability density function at x.

poisson

public static double poisson(int k,
                             double theta)
Evaluates the Poisson probability density function.

Method poisson evaluates the probability density function of a Poisson random variable with parameter theta. theta, which is the mean of the Poisson random variable, must be positive. The probability function (with theta = theta) is

f(x) = e^{- theta} ,,theta ^k /k!,,,,,,
  for,k = 0,,,1,,,2,, ldots

poisson evaluates this function directly, taking logarithms and using the log gamma function.

Poisson ProbabilityFunction

Parameters:
k - the int argument for which the Poisson probability function is to be evaluated.
theta - a double scalar value representing the mean of the Poisson distribution.
Returns:
a double scalar value representing the probability that a Poisson random variable takes a value equal to k.
See Also:
Example

Rayleigh

public static double Rayleigh(double x,
                              double alpha)
Evaluates the Rayleigh probability density function.

Parameters:
x - a double scalar value representing the argument at which the function is to be evaluated. It must be non-negative.
alpha - a double scalar value representing the scale parameter.
Returns:
a double scalar value representing the probability density function at x.

Weibull

public static double Weibull(double x,
                             double gamma,
                             double alpha)
Evaluates the Weibull probability density function.

Parameters:
x - a double scalar value representing the argument at which the function is to be evaluated. It must be non-negative.
gamma - a double scalar value representing the shape parameter.
alpha - a double scalar value representing the scale parameter.
Returns:
a double scalar value, the probability density function at x.

JMSLTM Numerical Library 6.0

Copyright © 1970-2009 Visual Numerics, Inc.
Built September 1 2009.