Beta Distribution

Table of contents


Density Function

The density function of the Beta distribution:

\[f(x; a,b) = \frac{1}{\mathcal{B}(a,b)} x^{a-1} (1-x)^{b-1} \times \mathbf{1}[0 \leq x \leq 1]\]

where \(\mathcal{B}(a,b)\) denotes the Beta 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> dbeta(const T1 x, const T2 a_par, const T3 b_par, const bool log_form = false) noexcept

Density function of the Beta distribution.

Example:

stats::dbeta(0.5,3.0,2.0,false); 

Parameters
  • x – a real-valued input.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

  • 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> dbeta(const std::vector<eT> &x, const T1 a_par, const T2 b_par, const bool log_form = false)

Density function of the Beta distribution.

Example:

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

Parameters
  • x – a standard vector.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

  • 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> dbeta(const ArmaMat<eT> &X, const T1 a_par, const T2 b_par, const bool log_form = false)

Density function of the Beta distribution.

Example:

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

Parameters
  • X – a matrix of input values.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

  • 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> dbeta(const BlazeMat<eT, To> &X, const T1 a_par, const T2 b_par, const bool log_form = false)

Density function of the Beta distribution.

Example:

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

Parameters
  • X – a matrix of input values.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

  • 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> dbeta(const EigenMat<eT, iTr, iTc> &X, const T1 a_par, const T2 b_par, const bool log_form = false)

Density function of the Beta distribution.

Example:

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

Parameters
  • X – a matrix of input values.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

  • 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 Beta distribution:

\[F(x; a, b) = \int_0^x f(z; a,b) dz = I_x (a,b)\]

where \(I_x (a,b)\) denotes the regularized incomplete Beta 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> pbeta(const T1 x, const T2 a_par, const T3 b_par, const bool log_form = false) noexcept

Distribution function of the Beta distribution.

Example:

stats::pbeta(0.5,3.0,2.0,false); 

Parameters
  • x – a real-valued input.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

  • 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> pbeta(const std::vector<eT> &x, const T1 a_par, const T2 b_par, const bool log_form = false)

Distribution function of the Beta distribution.

Example:

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

Parameters
  • x – a standard vector.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

  • 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> pbeta(const ArmaMat<eT> &X, const T1 a_par, const T2 b_par, const bool log_form = false)

Distribution function of the Beta distribution.

Example:

arma::mat X = { {0.2,  0.7,  0.1},
                {0.9, -0.3,  1.3} };
stats::pbeta(X,3.0,2.0,false);

Parameters
  • X – a matrix of input values.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

  • 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> pbeta(const BlazeMat<eT, To> &X, const T1 a_par, const T2 b_par, const bool log_form = false)

Distribution function of the Beta distribution.

Example:

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

Parameters
  • X – a matrix of input values.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

  • 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> pbeta(const EigenMat<eT, iTr, iTc> &X, const T1 a_par, const T2 b_par, const bool log_form = false)

Distribution function of the Beta distribution.

Example:

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

Parameters
  • X – a matrix of input values.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

  • 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 Beta distribution:

\[q(p; a,b) = \inf \left\{ x : p \leq I_x (a,b) \right\}\]

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> qbeta(const T1 p, const T2 a_par, const T3 b_par) noexcept

Quantile function of the Beta distribution.

Example:

stats::qbeta(0.5,3.0,2.0); 

Parameters
  • p – a real-valued input.

  • a_par – shape parameter, a real-valued input.

  • b_par – 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> qbeta(const std::vector<eT> &x, const T1 a_par, const T2 b_par)

Quantile function of the Beta distribution.

Example:

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

Parameters
  • x – a standard vector.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

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> qbeta(const ArmaMat<eT> &X, const T1 a_par, const T2 b_par)

Quantile function of the Beta distribution.

Example:

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

Parameters
  • X – a matrix of input values.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

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> qbeta(const BlazeMat<eT, To> &X, const T1 a_par, const T2 b_par)

Quantile function of the Beta distribution.

Example:

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

Parameters
  • X – a matrix of input values.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

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> qbeta(const EigenMat<eT, iTr, iTc> &X, const T1 a_par, const T2 b_par)

Quantile function of the Beta distribution.

Example:

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

Parameters
  • X – a matrix of input values.

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

Returns

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


Random Sampling

Random sampling for the Beta distribution is achieved by simulating two independent gamma-distributed random variables, \(X \sim G(a,1), Y \sim G(a,1)\), then returning:

\[Z = \frac{X}{X+Y} \sim B(a,b)\]

Scalar Output

  1. Random number engines

template<typename T1, typename T2>
inline common_return_t<T1, T2> rbeta(const T1 a_par, const T2 b_par, rand_engine_t &engine)

Random sampling function for the Beta distribution.

Example:

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

Parameters
  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

  • engine – a random engine, passed by reference.

Returns

a pseudo-random draw from the Beta distribution.

  1. Seed values

template<typename T1, typename T2>
inline common_return_t<T1, T2> rbeta(const T1 a_par, const T2 b_par, const ullint_t seed_val = std::random_device{}())

Random sampling function for the Beta distribution.

Example:

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

Parameters
  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

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

Returns

a pseudo-random draw from the Beta distribution.

Vector/Matrix Output

  1. Random number engines

template<typename mT, typename T1, typename T2>
inline mT rbeta(const ullint_t n, const ullint_t k, const T1 a_par, const T2 b_par, rand_engine_t &engine)

Random matrix sampling function for the Beta distribution.

Example:

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

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

  • engine – a random engine, passed by reference.

Returns

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

  1. Seed values

template<typename mT, typename T1, typename T2>
inline mT rbeta(const ullint_t n, const ullint_t k, const T1 a_par, const T2 b_par, const ullint_t seed_val = std::random_device{}())

Random matrix sampling function for the Beta distribution.

Example:

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

  • a_par – a real-valued shape parameter.

  • b_par – a real-valued shape parameter.

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

Returns

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