1 #ifndef __STAN__PROB__TRANSFORM_HPP__
2 #define __STAN__PROB__TRANSFORM_HPP__
9 #include <boost/multi_array.hpp>
10 #include <boost/throw_exception.hpp>
11 #include <boost/math/tools/promotion.hpp>
43 Eigen::Array<T,Eigen::Dynamic,1>& sds,
44 const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& Sigma) {
46 size_t K = sds.rows();
48 sds = Sigma.diagonal().array();
49 if( (sds <= 0.0).any() )
return false;
52 Eigen::DiagonalMatrix<T,Eigen::Dynamic> D(K);
53 D.diagonal() = sds.inverse();
56 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> R = D * Sigma * D;
58 R.diagonal().setOnes();
59 Eigen::LDLT<Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> > ldlt;
61 if( !ldlt.isPositive() )
return false;
62 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> U = ldlt.matrixU();
67 Eigen::Array<T,Eigen::Dynamic,1> temp = U.row(0).tail(pull);
68 CPCs.head(pull) = temp;
70 Eigen::Array<T,Eigen::Dynamic,1> acc(K);
72 acc.tail(pull) = 1.0 - temp.square();
73 for(
size_t i = 1; i < (K - 1); i++) {
76 temp = U.row(i).tail(pull);
77 temp /=
sqrt(acc.tail(pull) / acc(i));
78 CPCs.segment(position, pull) = temp;
79 acc.tail(pull) *= 1.0 - temp.square();
81 CPCs = 0.5 * ( (1.0 + CPCs) / (1.0 - CPCs) ).
log();
108 template <
typename T>
109 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
112 Eigen::Array<T,Eigen::Dynamic,1> temp;
113 Eigen::Array<T,Eigen::Dynamic,1> acc(K-1);
116 Eigen::Array<T,Eigen::Dynamic,Eigen::Dynamic> L(K,K);
123 L.col(0).tail(pull) = temp = CPCs.head(pull);
124 acc.tail(pull) = 1.0 - temp.square();
125 for(
size_t i = 1; i < (K - 1); i++) {
128 temp = CPCs.segment(position, pull);
129 L(i,i) =
sqrt(acc(i-1));
130 L.col(i).tail(pull) = temp * acc.tail(pull).sqrt();
131 acc.tail(pull) *= 1.0 - temp.square();
133 L(K-1,K-1) =
sqrt(acc(K-2));
150 template <
typename T>
151 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
154 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L
186 template <
typename T>
187 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
195 double lead = K - 2.0;
201 j < (CPCs.rows() - 1);
207 log_prob += lead / 2.0 * log_1cpc2;
236 template <
typename T>
237 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
242 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L
258 template <
typename T>
259 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
261 const Eigen::Array<T,Eigen::Dynamic,1>& sds,
263 size_t K = sds.rows();
266 return sds.matrix().asDiagonal() *
read_corr_L(CPCs, K, log_prob);
278 template <
typename T>
279 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
281 const Eigen::Array<T,Eigen::Dynamic,1>& sds,
284 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L
298 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
300 const Eigen::Array<T,Eigen::Dynamic,1>& sds) {
302 size_t K = sds.rows();
303 Eigen::DiagonalMatrix<T,Eigen::Dynamic> D(K);
305 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L
322 const Eigen::Array<T,Eigen::Dynamic,1>
325 Eigen::Array<T,Eigen::Dynamic,1> nu(K * (K - 1) / 2);
327 T alpha = eta + (K - 2.0) / 2.0;
331 T alpha2 = 2.0 * alpha;
334 for (
size_type j = 0; j < (K - 1); j++) {
337 size_t counter = K - 1;
338 for (
size_type i = 1; i < (K - 1); i++) {
340 alpha2 = 2.0 * alpha;
342 nu(counter) = alpha2;
366 template <
typename T>
385 template <
typename T>
402 template <
typename T>
421 template <
typename T>
443 template <
typename T>
466 template <
typename T>
494 template <
typename T,
typename TL>
497 if (lb == -std::numeric_limits<double>::infinity())
518 template <
typename T,
typename TL>
520 typename boost::math::tools::promote_args<T,TL>::type
522 if (lb == -std::numeric_limits<double>::infinity())
543 template <
typename T,
typename TL>
545 typename boost::math::tools::promote_args<T,TL>::type
547 if (lb == -std::numeric_limits<double>::infinity())
550 y, lb,
"Lower bounded variable");
576 template <
typename T,
typename TU>
578 typename boost::math::tools::promote_args<T,TU>::type
580 if (ub == std::numeric_limits<double>::infinity())
608 template <
typename T,
typename TU>
610 typename boost::math::tools::promote_args<T,TU>::type
612 if (ub == std::numeric_limits<double>::infinity())
640 template <
typename T,
typename TU>
642 typename boost::math::tools::promote_args<T,TU>::type
644 if (ub == std::numeric_limits<double>::infinity())
647 y, ub,
"Upper bounded variable");
681 template <
typename T,
typename TL,
typename TU>
683 typename boost::math::tools::promote_args<T,TL,TU>::type
687 if (lb == -std::numeric_limits<double>::infinity())
689 if (ub == std::numeric_limits<double>::infinity())
694 T exp_minus_x =
exp(-x);
695 inv_logit_x = 1.0 / (1.0 + exp_minus_x);
697 if ((x < std::numeric_limits<double>::infinity())
698 && (inv_logit_x == 1))
699 inv_logit_x = 1 - 1
e-15;
702 inv_logit_x = 1.0 - 1.0 / (1.0 + exp_x);
704 if ((x > -std::numeric_limits<double>::infinity())
705 && (inv_logit_x== 0))
708 return lb + (ub - lb) * inv_logit_x;
752 template <
typename T,
typename TL,
typename TU>
753 typename boost::math::tools::promote_args<T,TL,TU>::type
757 s <<
"domain error in lub_constrain; lower bound = " << lb
758 <<
" must be strictly less than upper bound = " << ub;
759 throw std::domain_error(s.str());
761 if (lb == -std::numeric_limits<double>::infinity())
763 if (ub == std::numeric_limits<double>::infinity())
767 T exp_minus_x =
exp(-x);
768 inv_logit_x = 1.0 / (1.0 + exp_minus_x);
769 lp +=
log(ub - lb) - x - 2 *
log1p(exp_minus_x);
771 if ((x < std::numeric_limits<double>::infinity())
772 && (inv_logit_x == 1))
773 inv_logit_x = 1 - 1
e-15;
776 inv_logit_x = 1.0 - 1.0 / (1.0 + exp_x);
777 lp +=
log(ub - lb) + x - 2 *
log1p(exp_x);
779 if ((x > -std::numeric_limits<double>::infinity())
780 && (inv_logit_x== 0))
783 return lb + (ub - lb) * inv_logit_x;
816 template <
typename T,
typename TL,
typename TU>
818 typename boost::math::tools::promote_args<T,TL,TU>::type
822 y, lb, ub,
"Bounded variable");
823 if (lb == -std::numeric_limits<double>::infinity())
825 if (ub == std::numeric_limits<double>::infinity())
827 return logit((y - lb) / (ub - lb));
846 template <
typename T>
874 template <
typename T>
880 lp +=
log(inv_logit_x) +
log1m(inv_logit_x);
898 template <
typename T>
903 y, 0, 1,
"Probability variable");
922 template <
typename T>
940 template <
typename T>
945 lp +=
log1m(tanh_x * tanh_x);
965 template <
typename T>
969 y, -1, 1,
"Correlation variable");
989 template <
typename T>
990 Eigen::Matrix<T,Eigen::Dynamic,1>
998 Eigen::Matrix<T,Eigen::Dynamic,1> x(Km1 + 1);
1002 x(k) = stick_len * z_k;
1022 template <
typename T>
1023 Eigen::Matrix<T,Eigen::Dynamic,1>
1032 Eigen::Matrix<T,Eigen::Dynamic,1> x(Km1 + 1);
1035 double eq_share = -
log(Km1 - k);
1036 T adj_y_k(y(k) + eq_share);
1038 x(k) = stick_len * z_k;
1039 lp +=
log(stick_len);
1062 template <
typename T>
1063 Eigen::Matrix<T,Eigen::Dynamic,1>
1068 int Km1 = x.size() - 1;
1069 Eigen::Matrix<T,Eigen::Dynamic,1> y(Km1);
1070 T stick_len(x(Km1));
1073 T z_k(x(k) / stick_len);
1092 template <
typename T>
1093 Eigen::Matrix<T,Eigen::Dynamic,1>
1097 Eigen::Matrix<T,Eigen::Dynamic,1> y(k);
1104 y[i] = y[i-1] +
exp(x[i]);
1120 template <
typename T>
1122 Eigen::Matrix<T,Eigen::Dynamic,1>
1125 for (
size_type i = 1; i < x.size(); ++i)
1145 template <
typename T>
1146 Eigen::Matrix<T,Eigen::Dynamic,1>
1149 y,
"Ordered variable");
1152 Eigen::Matrix<T,Eigen::Dynamic,1> x(k);
1157 x[i] =
log(y[i] - y[i-1]);
1173 template <
typename T>
1174 Eigen::Matrix<T,Eigen::Dynamic,1>
1178 Eigen::Matrix<T,Eigen::Dynamic,1> y(k);
1185 y[i] = y[i-1] +
exp(x[i]);
1201 template <
typename T>
1203 Eigen::Matrix<T,Eigen::Dynamic,1>
1206 for (
size_type i = 0; i < x.size(); ++i)
1226 template <
typename T>
1227 Eigen::Matrix<T,Eigen::Dynamic,1>
1230 y,
"Positive ordered variable");
1233 Eigen::Matrix<T,Eigen::Dynamic,1> x(k);
1238 x[i] =
log(y[i] - y[i-1]);
1268 template <
typename T>
1269 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1273 size_type k_choose_2 = (k * (k - 1)) / 2;
1274 if (k_choose_2 != x.size())
1275 throw std::invalid_argument (
"x is not a valid correlation matrix");
1276 Eigen::Array<T,Eigen::Dynamic,1> cpcs(k_choose_2);
1277 for (
size_type i = 0; i < k_choose_2; ++i)
1301 template <
typename T>
1302 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1307 size_type k_choose_2 = (k * (k - 1)) / 2;
1308 if (k_choose_2 != x.size())
1309 throw std::invalid_argument (
"x is not a valid correlation matrix");
1311 Eigen::Array<T,Eigen::Dynamic,1> cpcs(k_choose_2);
1312 for (
size_type i = 0; i < k_choose_2; ++i)
1337 template <
typename T>
1338 Eigen::Matrix<T,Eigen::Dynamic,1>
1344 throw std::domain_error(
"y is not a square matrix");
1346 throw std::domain_error(
"y has no elements");
1348 Eigen::Array<T,Eigen::Dynamic,1> x(k_choose_2);
1349 Eigen::Array<T,Eigen::Dynamic,1> sds(k);
1352 throw std::runtime_error(
"factor_cov_matrix failed on y");
1356 std::stringstream s;
1357 s <<
"all standard deviations must be zero."
1358 <<
" found log(sd[" << i <<
"])=" << sds[i] << std::endl;
1359 BOOST_THROW_EXCEPTION(std::runtime_error(s.str()));
1380 template <
typename T>
1381 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1386 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L(K,K);
1387 if (x.size() != (K * (K + 1)) / 2)
1388 throw std::domain_error(
"x.size() != K + (K choose 2)");
1391 for (
int n = 0; n < m; ++n)
1393 L(m,m) =
exp(x(i++));
1397 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> M = L * L.transpose();
1398 return L * L.transpose();
1414 template <
typename T>
1415 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1421 if (x.size() != (K * (K + 1)) / 2)
1422 throw std::domain_error(
"x.size() != K + (K choose 2)");
1423 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L(K,K);
1428 L(m,m) =
exp(x(i++));
1434 for (
int k = 0; k < K; ++k)
1435 lp += (K - k + 1) *
log(L(k,k));
1436 return L * L.transpose();
1462 template <
typename T>
1463 Eigen::Matrix<T,Eigen::Dynamic,1>
1468 throw std::domain_error(
"y is not a square matrix");
1470 throw std::domain_error(
"y has no elements");
1471 for (
int k = 0; k < K; ++k)
1472 if (!(y(k,k) > 0.0))
1473 throw std::domain_error(
"y has non-positive diagonal");
1474 Eigen::Matrix<T,Eigen::Dynamic,1> x((K * (K + 1)) / 2);
1477 Eigen::LLT<Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> >
1480 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L = llt.matrixL();
1482 for (
int m = 0; m < K; ++m) {
1483 for (
int n = 0; n < m; ++n)
1485 x(i++) =
log(L(m,m));
1509 template <
typename T>
1510 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1513 size_t k_choose_2 = (k * (k - 1)) / 2;
1514 Eigen::Array<T,Eigen::Dynamic,1> cpcs(k_choose_2);
1516 for (
size_t i = 0; i < k_choose_2; ++i)
1518 Eigen::Array<T,Eigen::Dynamic,1> sds(k);
1519 for (
size_t i = 0; i < k; ++i)
1548 template <
typename T>
1549 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1553 size_t k_choose_2 = (k * (k - 1)) / 2;
1554 Eigen::Array<T,Eigen::Dynamic,1> cpcs(k_choose_2);
1556 for (
size_t i = 0; i < k_choose_2; ++i)
1558 Eigen::Array<T,Eigen::Dynamic,1> sds(k);
1559 for (
size_t i = 0; i < k; ++i)
1582 template <
typename T>
1583 Eigen::Matrix<T,Eigen::Dynamic,1>
1585 const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& y) {
1590 throw std::domain_error(
"y is not a square matrix");
1592 throw std::domain_error(
"y has no elements");
1594 Eigen::Array<T,Eigen::Dynamic,1> cpcs(k_choose_2);
1595 Eigen::Array<T,Eigen::Dynamic,1> sds(k);
1598 throw std::runtime_error (
"factor_cov_matrix failed on y");
1599 Eigen::Matrix<T,Eigen::Dynamic,1> x(k_choose_2 + k);
1601 for (
size_type i = 0; i < k_choose_2; ++i)
internal::traits< Derived >::Index size_type
var inv_logit(const stan::agrad::var &a)
The inverse logit function for variables (stan).
matrix_v multiply_lower_tri_self_transpose(const matrix_v &L)
var log1p_exp(const stan::agrad::var &a)
Return the log of 1 plus the exponential of the specified variable.
var sqrt(const var &a)
Return the square root of the specified variable (cmath).
var fabs(const var &a)
Return the absolute value of the variable (cmath).
var square(const var &x)
Return the square of the input variable.
var log1m(const stan::agrad::var &a)
The log (1 - x) function for variables.
var log1p(const stan::agrad::var &a)
The log (1 + x) function for variables (C99).
var tanh(const var &a)
Return the hyperbolic tangent of the specified variable (cmath).
var log(const var &a)
Return the natural log of the specified variable (cmath).
var atanh(const stan::agrad::var &a)
The inverse hyperbolic tangent function for variables (C99).
var exp(const var &a)
Return the exponentiation of the specified variable (cmath).
bool check_simplex(const char *function, const Eigen::Matrix< T_prob, Eigen::Dynamic, 1 > &theta, const char *name, T_result *result, const Policy &)
Return true if the specified vector is simplex.
boost::math::tools::promote_args< T >::type logit(T a)
Returns the logit function applied to the argument.
bool check_ordered(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, 1 > &y, const char *name, T_result *result, const Policy &)
Return true if the specified vector is sorted into increasing order.
double e()
Return the base of the natural logarithm.
matrix_d multiply_lower_tri_self_transpose(const matrix_d &L)
Returns the result of multiplying the lower triangular portion of the input matrix by its own transpo...
const double LOG_2
The natural logarithm of 2, .
boost::math::tools::promote_args< T >::type log1m(T x)
Return the natural logarithm of one minus the specified value.
bool check_greater_or_equal(const char *function, const T_y &y, const T_low &low, const char *name, T_result *result, const Policy &)
bool check_positive_ordered(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, 1 > &y, const char *name, T_result *result, const Policy &)
Return true if the specified vector contains only non-negative values and is sorted into increasing o...
const double E
The base of the natural logarithm, .
boost::math::tools::promote_args< T >::type inv_logit(T a)
Returns the inverse logit function applied to the argument.
void validate_less(const T1 &x, const T2 &y, const char *x_name, const char *y_name, const char *fun_name)
T square(T x)
Return the square of the specified argument.
double log1p_exp(const double &a)
Calculates the log of 1 plus the exponential of the specified value without overflow.
bool check_less_or_equal(const char *function, const T_y &y, const T_high &high, const char *name, T_result *result, const Policy &)
bool check_bounded(const char *function, const T_y &y, const T_low &low, const T_high &high, const char *name, T_result *result, const Policy &)
bool check_positive(const char *function, const T_y &y, const char *name, T_result *result, const Policy &)
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > read_corr_matrix(const Eigen::Array< T, Eigen::Dynamic, 1 > &CPCs, const size_t K)
Return the correlation matrix of the specified dimensionality corresponding to the specified canonica...
bool factor_cov_matrix(Eigen::Array< T, Eigen::Dynamic, 1 > &CPCs, Eigen::Array< T, Eigen::Dynamic, 1 > &sds, const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &Sigma)
This function is intended to make starting values, given a covariance matrix Sigma.
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > read_corr_L(const Eigen::Array< T, Eigen::Dynamic, 1 > &CPCs, const size_t K)
Return the Cholesky factor of the correlation matrix of the specified dimensionality corresponding to...
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > read_cov_L(const Eigen::Array< T, Eigen::Dynamic, 1 > &CPCs, const Eigen::Array< T, Eigen::Dynamic, 1 > &sds, T &log_prob)
This is the function that should be called prior to evaluating the density of any elliptical distribu...
T corr_constrain(const T x)
Return the result of transforming the specified scalar to have a valid correlation value between -1 a...
Eigen::Matrix< T, Eigen::Dynamic, 1 > positive_ordered_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x)
Return an increasing positive ordered vector derived from the specified free vector.
Eigen::Matrix< T, Eigen::Dynamic, 1 > simplex_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &y)
Return the simplex corresponding to the specified free vector.
boost::math::tools::promote_args< T, TU >::type ub_constrain(const T x, const TU ub)
Return the upper-bounded value for the specified unconstrained scalar and upper bound.
Eigen::Matrix< T, Eigen::Dynamic, 1 > cov_matrix_free(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &y)
The covariance matrix derived from the symmetric view of the lower-triangular view of the K by K spec...
const double CONSTRAINT_TOLERANCE
Eigen::Matrix< T, Eigen::Dynamic, 1 > ordered_free(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &y)
Return the vector of unconstrained scalars that transform to the specified positive ordered vector.
const Eigen::Array< T, Eigen::Dynamic, 1 > make_nu(const T eta, const size_t K)
This function calculates the degrees of freedom for the t distribution that corresponds to the shape ...
boost::math::tools::promote_args< T, TL, TU >::type lub_free(const T y, TL lb, TU ub)
Return the unconstrained scalar that transforms to the specified lower- and upper-bounded scalar give...
Eigen::Matrix< T, Eigen::Dynamic, 1 > ordered_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x)
Return an increasing ordered vector derived from the specified free vector.
boost::math::tools::promote_args< T, TL >::type lb_free(const T y, const TL lb)
Return the unconstrained value that produces the specified lower-bound constrained value.
Eigen::Matrix< T, Eigen::Dynamic, 1 > simplex_free(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x)
Return an unconstrained vector that when transformed produces the specified simplex.
T prob_constrain(const T x)
Return a probability value constrained to fall between 0 and 1 (inclusive) for the specified free sca...
T lb_constrain(const T x, const TL lb)
Return the lower-bounded value for the specified unconstrained input and specified lower bound.
T prob_free(const T y)
Return the free scalar that when transformed to a probability produces the specified scalar.
boost::math::tools::promote_args< T, TL, TU >::type lub_constrain(const T x, TL lb, TU ub)
Return the lower- and upper-bounded scalar derived by transforming the specified free scalar given th...
Eigen::Matrix< T, Eigen::Dynamic, 1 > cov_matrix_free_lkj(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &y)
Return the vector of unconstrained partial correlations and deviations that transform to the specifie...
T positive_constrain(const T x)
Return the positive value for the specified unconstrained input.
boost::math::tools::promote_args< T, TU >::type ub_free(const T y, const TU ub)
Return the free scalar that corresponds to the specified upper-bounded value with respect to the spec...
Eigen::Matrix< T, Eigen::Dynamic, 1 > corr_matrix_free(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &y)
Return the vector of unconstrained partial correlations that define the specified correlation matrix ...
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > read_cov_matrix(const Eigen::Array< T, Eigen::Dynamic, 1 > &CPCs, const Eigen::Array< T, Eigen::Dynamic, 1 > &sds, T &log_prob)
A generally worse alternative to call prior to evaluating the density of an elliptical distribution.
T identity_constrain(T x)
Returns the result of applying the identity constraint transform to the input.
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > corr_matrix_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x, typename Eigen::Matrix< T, Eigen::Dynamic, 1 >::size_type k)
Return the correlation matrix of the specified dimensionality derived from the specified vector of un...
T corr_free(const T y)
Return the unconstrained scalar that when transformed to a valid correlation produces the specified v...
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > cov_matrix_constrain_lkj(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x, size_t k)
Return the covariance matrix of the specified dimensionality derived from constraining the specified ...
T identity_free(const T y)
Returns the result of applying the inverse of the identity constraint transform to the input.
Eigen::Matrix< T, Eigen::Dynamic, 1 > positive_ordered_free(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &y)
Return the vector of unconstrained scalars that transform to the specified positive ordered vector.
T positive_free(const T y)
Return the unconstrained value corresponding to the specified positive-constrained value.
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > cov_matrix_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x, typename Eigen::Matrix< T, Eigen::Dynamic, 1 >::size_type K)
Return the symmetric, positive-definite matrix of dimensions K by K resulting from transforming the s...
Probability, optimization and sampling library.