1#ifndef STAN_MATH_PRIM_ERR_CHECK_POS_DEFINITE_HPP
2#define STAN_MATH_PRIM_ERR_CHECK_POS_DEFINITE_HPP
33template <
typename EigMat, require_matrix_t<EigMat>* =
nullptr>
45 Eigen::LDLT<Eigen::MatrixXd> cholesky =
value_of_rec(y_ref).ldlt();
46 if (cholesky.info() != Eigen::Success || !cholesky.isPositive()
47 || (cholesky.vectorD().array() <= 0.0).any()) {
62template <
typename Derived>
64 const Eigen::LDLT<Derived>& cholesky) {
65 if (cholesky.info() != Eigen::Success || !cholesky.isPositive()
66 || !(cholesky.vectorD().array() > 0.0).all()) {
82template <
typename Derived>
84 const Eigen::LLT<Derived>& cholesky) {
85 if (cholesky.info() != Eigen::Success
86 || !(cholesky.matrixLLT().diagonal().array() > 0.0).all()) {
void check_symmetric(const char *function, const char *name, const matrix_cl< T > &y)
Check if the matrix_cl is symmetric.
double value_of_rec(const fvar< T > &v)
Return the value of the specified variable.
void check_pos_definite(const char *function, const char *name, const EigMat &y)
Check if the specified square, symmetric matrix is positive definite.
void throw_domain_error(const char *function, const char *name, const T &y, const char *msg1, const char *msg2)
Throw a domain error with a consistently formatted message.
ref_type_t< T && > to_ref(T &&a)
This evaluates expensive Eigen expressions.
void check_not_nan(const char *function, const char *name, const T_y &y)
Check if y is not NaN.
void check_positive(const char *function, const char *name, const T_y &y)
Check if y is positive.
const double CONSTRAINT_TOLERANCE
The tolerance for checking arithmetic bounds in rank and in simplexes.
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...