Inverse-Gaussian Distribution

Table of contents


Density Function

The density function of the inverse Gaussian distribution:

\[f(x; \mu, \lambda) = \sqrt{ \dfrac{\lambda}{2 \pi x^3} } \exp \left[ - \dfrac{\lambda(x - \mu)^2}{2 \mu^2 x} \right] \times \mathbf{1}[ x \geq 0 ]\]

where \(\mu\) is the mean parameter and \(\lambda\) is the shape parameter.

Methods for scalar input, as well as for vector/matrix input, are listed below.

Scalar Input

template<typename T1, typename T2, typename T3>
constexpr common_return_t<T1, T2, T3> dinvgauss(const T1 x, const T2 mu_par, const T3 lambda_par, const bool log_form = false) noexcept

Density function of the inverse Gaussian distribution.

Example:

stats::dinvgauss(0.5,1.0,2.0,false); 

Parameters
  • x – a real-valued input.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • log_form – return the log-density or the true form.

Returns

the density function evaluated at x.

Vector/Matrix Input

STL Containers

template<typename eT, typename T1, typename T2, typename rT = common_return_t<eT, T1, T2>>
inline std::vector<rT> dinvgauss(const std::vector<eT> &x, const T1 mu_par, const T2 lambda_par, const bool log_form = false)

Density function of the inverse Gaussian distribution.

Example:

std::vector<double> x = {0.0, 1.0, 2.0};
stats::dinvgauss(x,1.0,2.0,false);

Parameters
  • x – a standard vector.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • log_form – return the log-density or the true form.

Returns

a vector of density function values corresponding to the elements of x.

Armadillo

template<typename eT, typename T1, typename T2, typename rT = common_return_t<eT, T1, T2>>
inline ArmaMat<rT> dinvgauss(const ArmaMat<eT> &X, const T1 mu_par, const T2 lambda_par, const bool log_form = false)

Density function of the inverse Gaussian distribution.

Example:

arma::mat X = { {0.2, -1.7,  0.1},
                {0.9,  4.0, -0.3} };
stats::dinvgauss(X,1.0,1.0,false);

Parameters
  • X – a matrix of input values.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • log_form – return the log-density or the true form.

Returns

a matrix of density function values corresponding to the elements of X.

Blaze

template<typename eT, typename T1, typename T2, typename rT = common_return_t<eT, T1, T2>, bool To = blaze::columnMajor>
inline BlazeMat<rT, To> dinvgauss(const BlazeMat<eT, To> &X, const T1 mu_par, const T2 lambda_par, const bool log_form = false)

Density function of the inverse Gaussian distribution.

Example:

stats::dinvgauss(X,1.0,1.0,false);

Parameters
  • X – a matrix of input values.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • log_form – return the log-density or the true form.

Returns

a matrix of density function values corresponding to the elements of X.

Eigen

template<typename eT, typename T1, typename T2, typename rT = common_return_t<eT, T1, T2>, int iTr = Eigen::Dynamic, int iTc = Eigen::Dynamic>
inline EigenMat<rT, iTr, iTc> dinvgauss(const EigenMat<eT, iTr, iTc> &X, const T1 mu_par, const T2 lambda_par, const bool log_form = false)

Density function of the inverse Gaussian distribution.

Example:

stats::dinvgauss(X,1.0,1.0,false);

Parameters
  • X – a matrix of input values.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • log_form – return the log-density or the true form.

Returns

a matrix of density function values corresponding to the elements of X.


Cumulative Distribution Function

The cumulative distribution function of the inverse Gaussian distribution:

\[F(x; \mu, \lambda) = \Phi \left( \sqrt{\dfrac{\lambda}{x}} \left( \frac{x}{\mu} - 1 \right) \right) + \exp \left( \dfrac{2 \lambda}{\mu} \right) \Phi \left( - \sqrt{\dfrac{\lambda}{x}} \left( \frac{x}{\mu} + 1 \right) \right)\]

where \(\Phi(\cdot)\) denotes the standard Normal CDF.

Methods for scalar input, as well as for vector/matrix input, are listed below.

Scalar Input

template<typename T1, typename T2, typename T3>
constexpr common_return_t<T1, T2, T3> pinvgauss(const T1 x, const T2 mu_par, const T3 lambda_par, const bool log_form = false) noexcept

Distribution function of the inverse Gaussian distribution.

Example:

stats::pinvgauss(2.0,1.0,2.0,false); 

Parameters
  • x – a real-valued input.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • log_form – return the log-probability or the true form.

Returns

the cumulative distribution function evaluated at x.

Vector/Matrix Input

STL Containers

template<typename eT, typename T1, typename T2, typename rT = common_return_t<eT, T1, T2>>
inline std::vector<rT> pinvgauss(const std::vector<eT> &x, const T1 mu_par, const T2 lambda_par, const bool log_form = false)

Distribution function of the inverse Gaussian distribution.

Example:

std::vector<double> x = {0.0, 1.0, 2.0};
stats::pinvgauss(x,1.0,2.0,false);

Parameters
  • x – a standard vector.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • log_form – return the log-probability or the true form.

Returns

a vector of CDF values corresponding to the elements of x.

Armadillo

template<typename eT, typename T1, typename T2, typename rT = common_return_t<eT, T1, T2>>
inline ArmaMat<rT> pinvgauss(const ArmaMat<eT> &X, const T1 mu_par, const T2 lambda_par, const bool log_form = false)

Distribution function of the inverse Gaussian distribution.

Example:

arma::mat X = { {0.2, -1.7,  0.1},
                {0.9,  4.0, -0.3} };
stats::pinvgauss(X,1.0,1.0,false);

Parameters
  • X – a matrix of input values.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • log_form – return the log-probability or the true form.

Returns

a matrix of CDF values corresponding to the elements of X.

Blaze

template<typename eT, typename T1, typename T2, typename rT = common_return_t<eT, T1, T2>, bool To = blaze::columnMajor>
inline BlazeMat<rT, To> pinvgauss(const BlazeMat<eT, To> &X, const T1 mu_par, const T2 lambda_par, const bool log_form = false)

Distribution function of the inverse Gaussian distribution.

Example:

stats::pinvgauss(X,1.0,1.0,false);

Parameters
  • X – a matrix of input values.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • log_form – return the log-probability or the true form.

Returns

a matrix of CDF values corresponding to the elements of X.

Eigen

template<typename eT, typename T1, typename T2, typename rT = common_return_t<eT, T1, T2>, int iTr = Eigen::Dynamic, int iTc = Eigen::Dynamic>
inline EigenMat<rT, iTr, iTc> pinvgauss(const EigenMat<eT, iTr, iTc> &X, const T1 mu_par, const T2 lambda_par, const bool log_form = false)

Distribution function of the inverse Gaussian distribution.

Example:

stats::pinvgauss(X,1.0,1.0,false);

Parameters
  • X – a matrix of input values.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • log_form – return the log-probability or the true form.

Returns

a matrix of CDF values corresponding to the elements of X.


Quantile Function

The quantile function of the inverse Gaussian distribution:

\[q(p; \mu, \lambda) = \inf \left\{ x : p \leq F(x; \mu, \lambda) \right\}\]

where \(F(\cdot)\) denotes the CDF of the inverse Gaussian distribution.

Methods for scalar input, as well as for vector/matrix input, are listed below.

Scalar Input

template<typename T1, typename T2, typename T3>
constexpr common_return_t<T1, T2, T3> qinvgauss(const T1 p, const T2 mu_par, const T3 lambda_par) noexcept

Quantile function of the inverse Gaussian distribution.

Example:

stats::qinvgauss(0.5,1.0,2.0); 

Parameters
  • p – a real-valued input.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

Returns

the quantile function evaluated at p.

Vector/Matrix Input

STL Containers

template<typename eT, typename T1, typename T2, typename rT = common_return_t<eT, T1, T2>>
inline std::vector<rT> qinvgauss(const std::vector<eT> &x, const T1 mu_par, const T2 lambda_par)

Quantile function of the inverse Gaussian distribution.

Example:

std::vector<double> x = {0.1, 0.3, 0.7};
stats::qinvgauss(x,1.0,2.0);

Parameters
  • x – a standard vector.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

Returns

a vector of quantile values corresponding to the elements of x.

Armadillo

template<typename eT, typename T1, typename T2, typename rT = common_return_t<eT, T1, T2>>
inline ArmaMat<rT> qinvgauss(const ArmaMat<eT> &X, const T1 mu_par, const T2 lambda_par)

Quantile function of the inverse Gaussian distribution.

Example:

arma::mat X = { {0.2, 0.7, 0.9},
                {0.1, 0.8, 0.3} };
stats::qinvgauss(X,1.0,1.0);

Parameters
  • X – a matrix of input values.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

Returns

a matrix of quantile values corresponding to the elements of X.

Blaze

template<typename eT, typename T1, typename T2, typename rT = common_return_t<eT, T1, T2>, bool To = blaze::columnMajor>
inline BlazeMat<rT, To> qinvgauss(const BlazeMat<eT, To> &X, const T1 mu_par, const T2 lambda_par)

Quantile function of the inverse Gaussian distribution.

Example:

stats::qinvgauss(X,1.0,1.0);

Parameters
  • X – a matrix of input values.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

Returns

a matrix of quantile values corresponding to the elements of X.

Eigen

template<typename eT, typename T1, typename T2, typename rT = common_return_t<eT, T1, T2>, int iTr = Eigen::Dynamic, int iTc = Eigen::Dynamic>
inline EigenMat<rT, iTr, iTc> qinvgauss(const EigenMat<eT, iTr, iTc> &X, const T1 mu_par, const T2 lambda_par)

Quantile function of the inverse Gaussian distribution.

Example:

stats::qinvgauss(X,1.0,1.0);

Parameters
  • X – a matrix of input values.

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

Returns

a matrix of quantile values corresponding to the elements of X.


Random Sampling

Scalar Output

  1. Random number engines

template<typename T1, typename T2>
inline common_return_t<T1, T2> rinvgauss(const T1 mu_par, const T2 lambda_par, rand_engine_t &engine)

Random sampling function for the inverse Gaussian distribution.

Example:

stats::rand_engine_t engine(1776);
stats::rinvgauss(1.0,2.0,engine);

Parameters
  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • engine – a random engine, passed by reference.

Returns

a pseudo-random draw from the inverse Gaussian distribution.

  1. Seed values

template<typename T1, typename T2>
inline common_return_t<T1, T2> rinvgauss(const T1 mu_par, const T2 lambda_par, const ullint_t seed_val = std::random_device{}())

Random sampling function for the inverse Gaussian distribution.

Example:

stats::rinvgauss(1.0,2.0,1776);

Parameters
  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • seed_val – initialize the random engine with a non-negative integral-valued seed.

Returns

a pseudo-random draw from the inverse Gaussian distribution.

Vector/Matrix Output

  1. Random number engines

template<typename mT, typename T1 = double, typename T2 = double>
inline mT rinvgauss(const ullint_t n, const ullint_t k, const T1 mu_par, const T2 lambda_par, rand_engine_t &engine)

Random matrix sampling function for the inverse Gaussian distribution.

Example:

stats::rand_engine_t engine(1776);
// std::vector
stats::rinvgauss<std::vector<double>>(5,4,1.0,2.0,engine);
// Armadillo matrix
stats::rinvgauss<arma::mat>(5,4,1.0,2.0,engine);
// Blaze dynamic matrix
stats::rinvgauss<blaze::DynamicMatrix<double,blaze::columnMajor>>(5,4,1.0,2.0,engine);
// Eigen dynamic matrix
stats::rinvgauss<Eigen::MatrixXd>(5,4,1.0,2.0,engine);

Note

This function requires template instantiation; acceptable output types include: std::vector, with element type float, double, etc., as well as Armadillo, Blaze, and Eigen dense matrices.

Parameters
  • n – the number of output rows

  • k – the number of output columns

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • engine – a random engine, passed by reference.

Returns

a matrix of pseudo-random draws from the inverse Gaussian distribution.

  1. Seed values

template<typename mT, typename T1 = double, typename T2 = double>
inline mT rinvgauss(const ullint_t n, const ullint_t k, const T1 mu_par, const T2 lambda_par, const ullint_t seed_val = std::random_device{}())

Random matrix sampling function for the inverse Gaussian distribution.

Example:

// std::vector
stats::rinvgauss<std::vector<double>>(5,4,1.0,2.0);
// Armadillo matrix
stats::rinvgauss<arma::mat>(5,4,1.0,2.0);
// Blaze dynamic matrix
stats::rinvgauss<blaze::DynamicMatrix<double,blaze::columnMajor>>(5,4,1.0,2.0);
// Eigen dynamic matrix
stats::rinvgauss<Eigen::MatrixXd>(5,4,1.0,2.0);

Note

This function requires template instantiation; acceptable output types include: std::vector, with element type float, double, etc., as well as Armadillo, Blaze, and Eigen dense matrices.

Parameters
  • n – the number of output rows

  • k – the number of output columns

  • mu_par – the mean parameter, a real-valued input.

  • lambda_par – the shape parameter, a real-valued input.

  • seed_val – initialize the random engine with a non-negative integral-valued seed.

Returns

a matrix of pseudo-random draws from the inverse Gaussian distribution.