Student’s t-Distribution¶
Table of contents
Density Function¶
The density function of the Student’s t distribution:
where \(\Gamma(\cdot)\) denotes the gamma function.
Methods for scalar input, as well as for vector/matrix input, are listed below.
Scalar Input¶
-
template<typename T1, typename T2>
constexpr common_return_t<T1, T2> dt(const T1 x, const T2 dof_par, const bool log_form = false) noexcept¶ Density function of the t-distribution.
Example:
stats::dt(0.37,11,false);
- Parameters
x – a real-valued input.
dof_par – the degrees of freedom 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 rT = common_return_t<eT, T1>>
inline std::vector<rT> dt(const std::vector<eT> &x, const T1 dof_par, const bool log_form = false)¶ Density function of the t-distribution.
Example:
std::vector<double> x = {1.8, 0.7, 4.2}; stats::dt(x,4,false);
- Parameters
x – a standard vector.
dof_par – the degrees of freedom 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 rT = common_return_t<eT, T1>>
inline ArmaMat<rT> dt(const ArmaMat<eT> &X, const T1 dof_par, const bool log_form = false)¶ Density function of the t-distribution.
Example:
arma::mat X = { {1.8, 0.7, 4.2}, {0.3, 5.3, 3.7} }; stats::dt(X,4,false);
- Parameters
X – a matrix of input values.
dof_par – the degrees of freedom 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 rT = common_return_t<eT, T1>, bool To = blaze::columnMajor>
inline BlazeMat<rT, To> dt(const BlazeMat<eT, To> &X, const T1 dof_par, const bool log_form = false)¶ Density function of the t-distribution.
Example:
stats::dt(X,4,false);
- Parameters
X – a matrix of input values.
dof_par – the degrees of freedom 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 rT = common_return_t<eT, T1>, int iTr = Eigen::Dynamic, int iTc = Eigen::Dynamic>
inline EigenMat<rT, iTr, iTc> dt(const EigenMat<eT, iTr, iTc> &X, const T1 dof_par, const bool log_form = false)¶ Density function of the t-distribution.
Example:
stats::dt(X,4,false);
- Parameters
X – a matrix of input values.
dof_par – the degrees of freedom 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 Student’s t distribution:
where \(\Gamma(\cdot)\) denotes the gamma function and \({}_2 F_1\) denotes the hypergeometric function.
Methods for scalar input, as well as for vector/matrix input, are listed below.
Scalar Input¶
-
template<typename T1, typename T2>
constexpr common_return_t<T1, T2> pt(const T1 x, const T2 dof_par, const bool log_form = false) noexcept¶ Distribution function of the t-distribution.
Example:
stats::pt(0.37,11,false);
- Parameters
x – a real-valued input.
dof_par – the degrees of freedom 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 rT = common_return_t<eT, T1>>
inline std::vector<rT> pt(const std::vector<eT> &x, const T1 dof_par, const bool log_form = false)¶ Distribution function of the t-distribution.
Example:
std::vector<double> x = {0.0, 1.0, 2.0}; stats::pt(x,4,false);
- Parameters
x – a standard vector.
dof_par – the degrees of freedom 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 rT = common_return_t<eT, T1>>
inline ArmaMat<rT> pt(const ArmaMat<eT> &X, const T1 dof_par, const bool log_form = false)¶ Distribution function of the t-distribution.
Example:
arma::mat X = { {0.2, -1.7, 0.1}, {0.9, 4.0, -0.3} }; stats::pt(X,4,false);
- Parameters
X – a matrix of input values.
dof_par – the degrees of freedom 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 rT = common_return_t<eT, T1>, bool To = blaze::columnMajor>
inline BlazeMat<rT, To> pt(const BlazeMat<eT, To> &X, const T1 dof_par, const bool log_form = false)¶ Distribution function of the t-distribution.
Example:
stats::pt(X,4,false);
- Parameters
X – a matrix of input values.
dof_par – the degrees of freedom 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 rT = common_return_t<eT, T1>, int iTr = Eigen::Dynamic, int iTc = Eigen::Dynamic>
inline EigenMat<rT, iTr, iTc> pt(const EigenMat<eT, iTr, iTc> &X, const T1 dof_par, const bool log_form = false)¶ Distribution function of the t-distribution.
Example:
stats::pt(X,4,false);
- Parameters
X – a matrix of input values.
dof_par – the degrees of freedom 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 Student’s t distribution:
Methods for scalar input, as well as for vector/matrix input, are listed below.
Scalar Input¶
-
template<typename T1, typename T2>
constexpr common_return_t<T1, T2> qt(const T1 p, const T2 dof_par) noexcept¶ Quantile function of the t-distribution.
Example:
stats::qt(0.5,11);
- Parameters
p – a real-valued input.
dof_par – the degrees of freedom parameter, a real-valued input.
- Returns
the quantile function evaluated at
p
.
Vector/Matrix Input¶
STL Containers¶
-
template<typename eT, typename T1, typename rT = common_return_t<eT, T1>>
inline std::vector<rT> qt(const std::vector<eT> &x, const T1 dof_par)¶ Quantile function of the t-distribution.
Example:
std::vector<double> x = {0.3, 0.5, 0.8}; stats::qt(x,4);
- Parameters
x – a standard vector.
dof_par – the degrees of freedom parameter, a real-valued input.
- Returns
a vector of quantile values corresponding to the elements of
x
.
Armadillo¶
-
template<typename eT, typename T1, typename rT = common_return_t<eT, T1>>
inline ArmaMat<rT> qt(const ArmaMat<eT> &X, const T1 dof_par)¶ Quantile function of the t-distribution.
Example:
arma::mat X = { {0.2, 0.7, 0.9}, {0.1, 0.8, 0.3} }; stats::qt(X,4);
- Parameters
X – a matrix of input values.
dof_par – the degrees of freedom parameter, a real-valued input.
- Returns
a matrix of quantile values corresponding to the elements of
X
.
Blaze¶
-
template<typename eT, typename T1, typename rT = common_return_t<eT, T1>, bool To = blaze::columnMajor>
inline BlazeMat<rT, To> qt(const BlazeMat<eT, To> &X, const T1 dof_par)¶ Quantile function of the t-distribution.
Example:
stats::qt(X,4);
- Parameters
X – a matrix of input values.
dof_par – the degrees of freedom parameter, a real-valued input.
- Returns
a matrix of quantile values corresponding to the elements of
X
.
Eigen¶
-
template<typename eT, typename T1, typename rT = common_return_t<eT, T1>, int iTr = Eigen::Dynamic, int iTc = Eigen::Dynamic>
inline EigenMat<rT, iTr, iTc> qt(const EigenMat<eT, iTr, iTc> &X, const T1 dof_par)¶ Quantile function of the t-distribution.
Example:
stats::qt(X,4);
- Parameters
X – a matrix of input values.
dof_par – the degrees of freedom parameter, a real-valued input.
- Returns
a matrix of quantile values corresponding to the elements of
X
.
Random Sampling¶
Scalar Output¶
Random number engines
-
template<typename T>
inline return_t<T> rt(const T dof_par, rand_engine_t &engine)¶ Random sampling function for Student’s t-distribution.
Example:
stats::rand_engine_t engine(1776); stats::rt(4,engine);
- Parameters
dof_par – the degrees of freedom parameter, a real-valued input.
engine – a random engine, passed by reference.
- Returns
a pseudo-random draw from Student’s t-distribution.
Seed values
-
template<typename T>
inline return_t<T> rt(const T dof_par, const ullint_t seed_val = std::random_device{}())¶ Random sampling function for Student’s t-distribution.
Example:
stats::rt(4,1776);
- Parameters
dof_par – the degrees of freedom parameter, a real-valued input.
seed_val – initialize the random engine with a non-negative integral-valued seed.
- Returns
a pseudo-random draw from Student’s t-distribution.
Vector/Matrix Output¶
Random number engines
-
template<typename mT, typename T1>
inline mT rt(const ullint_t n, const ullint_t k, const T1 dof_par, rand_engine_t &engine)¶ Random matrix sampling function for Student’s t-distribution.
Example:
stats::rand_engine_t engine(1776); // std::vector stats::rt<std::vector<double>>(5,4,12,engine); // Armadillo matrix stats::rt<arma::mat>(5,4,12,engine); // Blaze dynamic matrix stats::rt<blaze::DynamicMatrix<double,blaze::columnMajor>>(5,4,12,engine); // Eigen dynamic matrix stats::rt<Eigen::MatrixXd>(5,4,12,engine);
Note
This function requires template instantiation; acceptable output types include:
std::vector
, with element typefloat
,double
, etc., as well as Armadillo, Blaze, and Eigen dense matrices.- Parameters
n – the number of output rows
k – the number of output columns
dof_par – the degrees of freedom parameter, a real-valued input.
engine – a random engine, passed by reference.
- Returns
a matrix of pseudo-random draws from the t-distribution.
Seed values
-
template<typename mT, typename T1>
inline mT rt(const ullint_t n, const ullint_t k, const T1 dof_par, const ullint_t seed_val = std::random_device{}())¶ Random matrix sampling function for Student’s t-distribution.
Example:
// std::vector stats::rt<std::vector<double>>(5,4,12); // Armadillo matrix stats::rt<arma::mat>(5,4,12); // Blaze dynamic matrix stats::rt<blaze::DynamicMatrix<double,blaze::columnMajor>>(5,4,12); // Eigen dynamic matrix stats::rt<Eigen::MatrixXd>(5,4,12);
Note
This function requires template instantiation; acceptable output types include:
std::vector
, with element typefloat
,double
, etc., as well as Armadillo, Blaze, and Eigen dense matrices.- Parameters
n – the number of output rows
k – the number of output columns
dof_par – the degrees of freedom 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 Student’s t-distribution.