Gamma Distribution

Table of contents


Density Function

The density function of the Gamma distribution:

\[f(x; k, \theta) = \dfrac{x^{k-1}\exp(-x/\theta)}{\theta^k \Gamma(k)} \times \mathbf{1}[ x \geq 0 ]\]

where \(\Gamma(\cdot)\) denotes the Gamma function, \(k\) is the shape parameter, and \(\theta\) is the scale 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> dgamma(const T1 x, const T2 shape_par, const T3 scale_par, const bool log_form = false) noexcept

Density function of the Gamma distribution.

Example:

stats::dgamma(2,2,3,false); 

Parameters
  • x – a real-valued input.

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

  • scale_par – the scale 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> dgamma(const std::vector<eT> &x, const T1 shape_par, const T2 scale_par, const bool log_form = false)

Density function of the Gamma distribution.

Example:

std::vector<double> x = {1.8, 0.7, 4.2};
stats::dgamma(x,3.0,2.0,false);

Parameters
  • x – a standard vector.

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

  • scale_par – the scale 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> dgamma(const ArmaMat<eT> &X, const T1 shape_par, const T2 scale_par, const bool log_form = false)

Density function of the Gamma distribution.

Example:

arma::mat X = { {1.8, 0.7, 4.2},
                {0.3, 5.3, 3.7} };
stats::dgamma(X,3.0,2.0,false);

Parameters
  • X – a matrix of input values.

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

  • scale_par – the scale 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> dgamma(const BlazeMat<eT, To> &X, const T1 shape_par, const T2 scale_par, const bool log_form = false)

Density function of the Gamma distribution.

Example:

stats::dgamma(X,3.0,2.0,false);

Parameters
  • X – a matrix of input values.

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

  • scale_par – the scale 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> dgamma(const EigenMat<eT, iTr, iTc> &X, const T1 shape_par, const T2 scale_par, const bool log_form = false)

Density function of the Gamma distribution.

Example:

stats::dgamma(X,3.0,2.0,false);

Parameters
  • X – a matrix of input values.

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

  • scale_par – the scale 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 Gamma distribution:

\[F(x; k, \theta) = \int_0^x f(z; k, \theta) dz = \frac{\gamma(k,x\theta)}{\Gamma (k)}\]

where \(\Gamma(\cdot)\) denotes the gamma function and \(\gamma(\cdot, \cdot)\) denotes the incomplete gamma function.

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> pgamma(const T1 x, const T2 shape_par, const T3 scale_par, const bool log_form = false) noexcept

Distribution function of the Gamma distribution.

Example:

stats::pgamma(2,2,3,false); 

Parameters
  • x – a real-valued input.

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

  • scale_par – the scale 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> pgamma(const std::vector<eT> &x, const T1 shape_par, const T2 scale_par, const bool log_form = false)

Distribution function of the Gamma distribution.

Example:

std::vector<double> x = {1.8, 0.7, 4.2};
stats::pgamma(x,3.0,2.0,false);

Parameters
  • x – a standard vector.

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

  • scale_par – the scale 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> pgamma(const ArmaMat<eT> &X, const T1 shape_par, const T2 scale_par, const bool log_form = false)

Distribution function of the Gamma distribution.

Example:

arma::mat X = { {1.8, 0.7, 4.2},
                {0.3, 5.3, 3.7} };
stats::pgamma(X,3.0,2.0,false);

Parameters
  • X – a matrix of input values.

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

  • scale_par – the scale 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> pgamma(const BlazeMat<eT, To> &X, const T1 shape_par, const T2 scale_par, const bool log_form = false)

Distribution function of the Gamma distribution.

Example:

stats::pgamma(X,3.0,2.0,false);

Parameters
  • X – a matrix of input values.

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

  • scale_par – the scale 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> pgamma(const EigenMat<eT, iTr, iTc> &X, const T1 shape_par, const T2 scale_par, const bool log_form = false)

Distribution function of the Gamma distribution.

Example:

stats::pgamma(X,3.0,2.0,false);

Parameters
  • X – a matrix of input values.

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

  • scale_par – the scale 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 Gamma distribution:

\[q(p; k, \theta) = \inf \left\{ x : p \leq \frac{\gamma(k,x\theta)}{\Gamma (k)} \right\}\]

where \(\Gamma(\cdot)\) denotes the gamma function and \(\gamma(\cdot, \cdot)\) denotes the incomplete gamma function.

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> qgamma(const T1 p, const T2 shape_par, const T3 scale_par) noexcept

Quantile function of the Gamma distribution.

Example:

stats::qgamma(0.4,2,3); 

Parameters
  • p – a real-valued input.

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

  • scale_par – the scale 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> qgamma(const std::vector<eT> &x, const T1 shape_par, const T2 scale_par)

Quantile function of the Gamma distribution.

Example:

std::vector<double> x = {0.3, 0.5, 0.9};
stats::qgamma(x,3.0,2.0);

Parameters
  • x – a standard vector.

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

  • scale_par – the scale 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> qgamma(const ArmaMat<eT> &X, const T1 shape_par, const T2 scale_par)

Quantile function of the Gamma distribution.

Example:

arma::mat X = { {0.2,  0.7,  0.1},
                {0.9,  0.3,  0.87} };
stats::qgamma(X,3.0,2.0);

Parameters
  • X – a matrix of input values.

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

  • scale_par – the scale 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> qgamma(const BlazeMat<eT, To> &X, const T1 shape_par, const T2 scale_par)

Quantile function of the Gamma distribution.

Example:

stats::qgamma(X,3.0,2.0);

Parameters
  • X – a matrix of input values.

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

  • scale_par – the scale 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> qgamma(const EigenMat<eT, iTr, iTc> &X, const T1 shape_par, const T2 scale_par)

Quantile function of the Gamma distribution.

Example:

stats::qgamma(X,3.0,2.0);

Parameters
  • X – a matrix of input values.

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

  • scale_par – the scale parameter, a real-valued input.

Returns

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


Random Sampling

Random sampling for the Gamma distribution is achieved via the Ziggurat method of Marsaglia and Tsang (2000).

Scalar Output

  1. Random number engines

template<typename T1, typename T2>
inline common_return_t<T1, T2> rgamma(const T1 shape_par, const T2 scale_par, rand_engine_t &engine)

Random sampling function for the Gamma distribution.

Example:

stats::rand_engine_t engine(1776);
stats::rgamma(3.0,2.0,engine);

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

  • scale_par – the scale parameter, a real-valued input.

  • engine – a random engine, passed by reference.

Returns

a pseudo-random draw from the Gamma distribution.

  1. Seed values

template<typename T1, typename T2>
inline common_return_t<T1, T2> rgamma(const T1 shape_par, const T2 scale_par, const ullint_t seed_val = std::random_device{}())

Random sampling function for the Gamma distribution.

Example:

stats::rgamma(3.0,2.0,1776);

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

  • scale_par – the scale 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 Gamma distribution.

Vector/Matrix Output

  1. Random number engines

template<typename mT, typename T1, typename T2>
inline mT rgamma(const ullint_t n, const ullint_t k, const T1 shape_par, const T2 scale_par, rand_engine_t &engine)

Random matrix sampling function for the Gamma distribution.

Example:

stats::rand_engine_t engine(1776);
// std::vector
stats::rgamma<std::vector<double>>(5,4,3.0,2.0,engine);
// Armadillo matrix
stats::rgamma<arma::mat>(5,4,3.0,2.0,engine);
// Blaze dynamic matrix
stats::rgamma<blaze::DynamicMatrix<double,blaze::columnMajor>>(5,4,3.0,2.0,engine);
// Eigen dynamic matrix
stats::rgamma<Eigen::MatrixXd>(5,4,3.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

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

  • scale_par – the scale parameter, a real-valued input.

  • engine – a random engine, passed by reference.

Returns

a matrix of pseudo-random draws from the Gamma distribution.

  1. Seed values

template<typename mT, typename T1, typename T2>
inline mT rgamma(const ullint_t n, const ullint_t k, const T1 shape_par, const T2 scale_par, const ullint_t seed_val = std::random_device{}())

Random matrix sampling function for the Gamma distribution.

Example:

// std::vector
stats::rgamma<std::vector<double>>(5,4,3.0,2.0);
// Armadillo matrix
stats::rgamma<arma::mat>(5,4,3.0,2.0);
// Blaze dynamic matrix
stats::rgamma<blaze::DynamicMatrix<double,blaze::columnMajor>>(5,4,3.0,2.0);
// Eigen dynamic matrix
stats::rgamma<Eigen::MatrixXd>(5,4,3.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

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

  • scale_par – the scale 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 Gamma distribution.