1 #ifndef STAN__PROB__TRANSFORM_HPP
2 #define STAN__PROB__TRANSFORM_HPP
11 #include <boost/throw_exception.hpp>
12 #include <boost/math/tools/promotion.hpp>
43 factor_U(
const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& U,
44 Eigen::Array<T,Eigen::Dynamic,1>& CPCs) {
51 CPCs(0) =
atanh(U(0,1));
55 Eigen::Array<T,1,Eigen::Dynamic> temp = U.row(0).tail(pull);
57 CPCs.head(pull) = temp;
59 Eigen::Array<T,Eigen::Dynamic,1> acc(K);
61 acc.tail(pull) = 1.0 - temp.square();
62 for(
size_t i = 1; i < (K - 1); i++) {
65 temp = U.row(i).tail(pull);
66 temp /=
sqrt(acc.tail(pull) / acc(i));
67 CPCs.segment(position, pull) = temp;
68 acc.tail(pull) *= 1.0 - temp.square();
70 CPCs = 0.5 * ( (1.0 + CPCs) / (1.0 - CPCs) ).
log();
89 Eigen::Array<T,Eigen::Dynamic,1>& CPCs,
90 Eigen::Array<T,Eigen::Dynamic,1>& sds) {
92 size_t K = sds.rows();
94 sds = Sigma.diagonal().array();
95 if( (sds <= 0.0).any() )
return false;
98 Eigen::DiagonalMatrix<T,Eigen::Dynamic> D(K);
99 D.diagonal() = sds.inverse();
102 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> R = D * Sigma * D;
104 R.diagonal().setOnes();
105 Eigen::LDLT<Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> > ldlt;
107 if (!ldlt.isPositive())
109 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> U = ldlt.matrixU();
137 template <
typename T>
138 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
141 Eigen::Array<T,Eigen::Dynamic,1> temp;
142 Eigen::Array<T,Eigen::Dynamic,1> acc(K-1);
145 Eigen::Array<T,Eigen::Dynamic,Eigen::Dynamic> L(K,K);
152 L.col(0).tail(pull) = temp = CPCs.head(pull);
153 acc.tail(pull) = 1.0 - temp.square();
154 for(
size_t i = 1; i < (K - 1); i++) {
157 temp = CPCs.segment(position, pull);
158 L(i,i) =
sqrt(acc(i-1));
159 L.col(i).tail(pull) = temp * acc.tail(pull).sqrt();
160 acc.tail(pull) *= 1.0 - temp.square();
162 L(K-1,K-1) =
sqrt(acc(K-2));
179 template <
typename T>
180 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
183 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L
213 template <
typename T>
214 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
223 Eigen::Matrix<T,Eigen::Dynamic,1> values(CPCs.rows() - 1);
228 for (
size_t k = 1; k <= (K - 2); k++)
229 for (
size_t i = k + 1; i <= K; i++) {
230 values(pos) = (K - k - 1) *
log1m(
square(CPCs(pos)));
234 log_prob += 0.5 *
sum(values);
256 template <
typename T>
257 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
262 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L
278 template <
typename T>
279 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
281 const Eigen::Array<T,Eigen::Dynamic,1>& sds,
283 size_t K = sds.rows();
286 return sds.matrix().asDiagonal() *
read_corr_L(CPCs, K, log_prob);
298 template <
typename T>
299 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
301 const Eigen::Array<T,Eigen::Dynamic,1>& sds,
304 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L
318 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
320 const Eigen::Array<T,Eigen::Dynamic,1>& sds) {
322 size_t K = sds.rows();
323 Eigen::DiagonalMatrix<T,Eigen::Dynamic> D(K);
325 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L
342 const Eigen::Array<T,Eigen::Dynamic,1>
345 using Eigen::Dynamic;
348 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
350 Array<T,Dynamic,1> nu(K * (K - 1) / 2);
352 T alpha = eta + (K - 2.0) / 2.0;
356 T alpha2 = 2.0 * alpha;
357 for (size_type j = 0; j < (K - 1); j++) {
360 size_t counter = K - 1;
361 for (size_type i = 1; i < (K - 1); i++) {
363 alpha2 = 2.0 * alpha;
364 for (size_type j = i + 1; j < K; j++) {
365 nu(counter) = alpha2;
389 template <
typename T>
408 template <
typename T>
425 template <
typename T>
444 template <
typename T>
466 template <
typename T>
489 template <
typename T>
493 "Positive variable", (
double*)0);
518 template <
typename T,
typename TL>
521 if (lb == -std::numeric_limits<double>::infinity())
542 template <
typename T,
typename TL>
544 typename boost::math::tools::promote_args<T,TL>::type
546 if (lb == -std::numeric_limits<double>::infinity())
567 template <
typename T,
typename TL>
569 typename boost::math::tools::promote_args<T,TL>::type
571 if (lb == -std::numeric_limits<double>::infinity())
574 y, lb,
"Lower bounded variable",
601 template <
typename T,
typename TU>
603 typename boost::math::tools::promote_args<T,TU>::type
605 if (ub == std::numeric_limits<double>::infinity())
633 template <
typename T,
typename TU>
635 typename boost::math::tools::promote_args<T,TU>::type
637 if (ub == std::numeric_limits<double>::infinity())
665 template <
typename T,
typename TU>
667 typename boost::math::tools::promote_args<T,TU>::type
669 if (ub == std::numeric_limits<double>::infinity())
672 y, ub,
"Upper bounded variable",
707 template <
typename T,
typename TL,
typename TU>
709 typename boost::math::tools::promote_args<T,TL,TU>::type
713 if (lb == -std::numeric_limits<double>::infinity())
715 if (ub == std::numeric_limits<double>::infinity())
720 T exp_minus_x =
exp(-x);
721 inv_logit_x = 1.0 / (1.0 + exp_minus_x);
723 if ((x < std::numeric_limits<double>::infinity())
724 && (inv_logit_x == 1))
725 inv_logit_x = 1 - 1
e-15;
728 inv_logit_x = 1.0 - 1.0 / (1.0 + exp_x);
730 if ((x > -std::numeric_limits<double>::infinity())
731 && (inv_logit_x== 0))
734 return lb + (ub - lb) * inv_logit_x;
778 template <
typename T,
typename TL,
typename TU>
779 typename boost::math::tools::promote_args<T,TL,TU>::type
783 s <<
"domain error in lub_constrain; lower bound = " << lb
784 <<
" must be strictly less than upper bound = " << ub;
785 throw std::domain_error(s.str());
787 if (lb == -std::numeric_limits<double>::infinity())
789 if (ub == std::numeric_limits<double>::infinity())
793 T exp_minus_x =
exp(-x);
794 inv_logit_x = 1.0 / (1.0 + exp_minus_x);
795 lp +=
log(ub - lb) - x - 2 *
log1p(exp_minus_x);
797 if ((x < std::numeric_limits<double>::infinity())
798 && (inv_logit_x == 1))
799 inv_logit_x = 1 - 1
e-15;
802 inv_logit_x = 1.0 - 1.0 / (1.0 + exp_x);
803 lp +=
log(ub - lb) + x - 2 *
log1p(exp_x);
805 if ((x > -std::numeric_limits<double>::infinity())
806 && (inv_logit_x== 0))
809 return lb + (ub - lb) * inv_logit_x;
842 template <
typename T,
typename TL,
typename TU>
844 typename boost::math::tools::promote_args<T,TL,TU>::type
847 stan::math::check_bounded<T, TL, TU, typename scalar_type<T>::type>
848 (
"stan::prob::lub_free(%1%)",
849 y, lb, ub,
"Bounded variable", (
double*)0);
850 if (lb == -std::numeric_limits<double>::infinity())
852 if (ub == std::numeric_limits<double>::infinity())
854 return logit((y - lb) / (ub - lb));
873 template <
typename T>
901 template <
typename T>
907 lp +=
log(inv_logit_x) +
log1m(inv_logit_x);
925 template <
typename T>
929 stan::math::check_bounded<T,double,double,T>
930 (
"stan::prob::prob_free(%1%)",
931 y, 0, 1,
"Probability variable",(
double*)0);
950 template <
typename T>
968 template <
typename T>
973 lp +=
log1m(tanh_x * tanh_x);
993 template <
typename T>
996 stan::math::check_bounded<T,double,double,double>
997 (
"stan::prob::lub_free(%1%)",
998 y, -1, 1,
"Correlation variable", (
double*)0);
1013 template <
typename T>
1014 Eigen::Matrix<T,Eigen::Dynamic,1>
1016 using Eigen::Matrix;
1017 using Eigen::Dynamic;
1019 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1021 Matrix<T,Dynamic,1> x(Km1 + 1);
1023 const T half_pi = T(M_PI/2.0);
1024 for (size_type k = 1; k <= Km1; ++k) {
1025 T yk_1 = y(k-1) + half_pi;
1026 T sin_yk_1 =
sin(yk_1);
1027 x(k) = x(k-1)*sin_yk_1;
1028 x(k-1) *=
cos(yk_1);
1042 template <
typename T>
1043 Eigen::Matrix<T,Eigen::Dynamic,1>
1045 using Eigen::Matrix;
1046 using Eigen::Dynamic;
1048 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1051 Matrix<T,Dynamic,1> x(Km1 + 1);
1053 const T half_pi = T(0.5 * M_PI);
1054 for (size_type k = 1; k <= Km1; ++k) {
1055 T yk_1 = y(k-1) + half_pi;
1056 T sin_yk_1 =
sin(yk_1);
1057 x(k) = x(k-1) * sin_yk_1;
1058 x(k-1) *=
cos(yk_1);
1060 lp += (Km1 - k) *
log(
fabs(sin_yk_1));
1065 template <
typename T>
1066 Eigen::Matrix<T,Eigen::Dynamic,1>
1068 using Eigen::Matrix;
1069 using Eigen::Dynamic;
1071 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1074 x,
"Unit vector variable", (
double*)0);
1075 int Km1 = x.size() - 1;
1076 Matrix<T,Dynamic,1> y(Km1);
1077 T sumSq = x(Km1)*x(Km1);
1078 const T half_pi = T(M_PI/2.0);
1079 for (size_type k = Km1; --k >= 0; ) {
1080 y(k) =
atan2(
sqrt(sumSq),x(k)) - half_pi;
1101 template <
typename T>
1102 Eigen::Matrix<T,Eigen::Dynamic,1>
1106 using Eigen::Matrix;
1107 using Eigen::Dynamic;
1112 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1116 Matrix<T,Dynamic,1> x(Km1 + 1);
1118 for (size_type k = 0; k < Km1; ++k) {
1120 x(k) = stick_len * z_k;
1140 template <
typename T>
1141 Eigen::Matrix<T,Eigen::Dynamic,1>
1145 using Eigen::Dynamic;
1146 using Eigen::Matrix;
1153 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1156 Matrix<T,Dynamic,1> x(Km1 + 1);
1158 for (size_type k = 0; k < Km1; ++k) {
1159 double eq_share = -
log(Km1 - k);
1160 T adj_y_k(y(k) + eq_share);
1162 x(k) = stick_len * z_k;
1163 lp +=
log(stick_len);
1186 template <
typename T>
1187 Eigen::Matrix<T,Eigen::Dynamic,1>
1189 using Eigen::Dynamic;
1190 using Eigen::Matrix;
1194 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1198 int Km1 = x.size() - 1;
1199 Eigen::Matrix<T,Eigen::Dynamic,1> y(Km1);
1200 T stick_len(x(Km1));
1201 for (size_type k = Km1; --k >= 0; ) {
1203 T z_k(x(k) / stick_len);
1222 template <
typename T>
1223 Eigen::Matrix<T,Eigen::Dynamic,1>
1225 using Eigen::Matrix;
1226 using Eigen::Dynamic;
1229 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1231 size_type k = x.size();
1232 Matrix<T,Dynamic,1> y(k);
1236 for (size_type i = 1; i < k; ++i)
1237 y[i] = y[i-1] +
exp(x[i]);
1253 template <
typename T>
1255 Eigen::Matrix<T,Eigen::Dynamic,1>
1257 using Eigen::Matrix;
1258 using Eigen::Dynamic;
1261 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1263 for (size_type i = 1; i < x.size(); ++i)
1283 template <
typename T>
1284 Eigen::Matrix<T,Eigen::Dynamic,1>
1287 y,
"Ordered variable");
1288 using Eigen::Matrix;
1289 using Eigen::Dynamic;
1291 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1293 size_type k = y.size();
1294 Matrix<T,Dynamic,1> x(k);
1298 for (size_type i = 1; i < k; ++i)
1299 x[i] =
log(y[i] - y[i-1]);
1315 template <
typename T>
1316 Eigen::Matrix<T,Eigen::Dynamic,1>
1318 using Eigen::Matrix;
1319 using Eigen::Dynamic;
1321 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1323 size_type k = x.size();
1324 Matrix<T,Dynamic,1> y(k);
1328 for (size_type i = 1; i < k; ++i)
1329 y[i] = y[i-1] +
exp(x[i]);
1345 template <
typename T>
1347 Eigen::Matrix<T,Eigen::Dynamic,1>
1349 using Eigen::Matrix;
1350 using Eigen::Dynamic;
1352 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1354 for (size_type i = 0; i < x.size(); ++i)
1374 template <
typename T>
1375 Eigen::Matrix<T,Eigen::Dynamic,1>
1377 using Eigen::Matrix;
1378 using Eigen::Dynamic;
1381 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1384 y,
"Positive ordered variable", (
double*)0);
1386 size_type k = y.size();
1387 Matrix<T,Dynamic,1> x(k);
1391 for (size_type i = 1; i < k; ++i)
1392 x[i] =
log(y[i] - y[i-1]);
1410 template <
typename T>
1411 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1417 throw std::domain_error(
"cholesky_factor_constrain: num rows must be >= num cols");
1418 if (x.size() != ((N * (N + 1)) / 2 + (M - N) * N))
1419 throw std::domain_error(
"cholesky_factor_constrain: x.size() must"
1420 " be (N * (N + 1)) / 2 + (M - N) * N");
1421 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> y(M,N);
1425 for (
int m = 0; m < N; ++m) {
1426 for (
int n = 0; n < m; ++n)
1428 y(m,m) =
exp(x(pos++));
1429 for (
int n = m + 1; n < N; ++n)
1433 for (
int m = N; m < M; ++m)
1434 for (
int n = 0; n < N; ++n)
1453 template <
typename T>
1454 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1462 if (x.size() != ((N * (N + 1)) / 2 + (M - N) * N))
1463 throw std::domain_error(
"cholesky_factor_constrain: x.size() must be (k choose 2) + k");
1465 std::vector<T> log_jacobians(N);
1466 for (
int n = 0; n < N; ++n) {
1468 log_jacobians[n] = x(pos++);
1470 lp +=
sum(log_jacobians);
1483 template <
typename T>
1484 Eigen::Matrix<T,Eigen::Dynamic,1>
1489 throw std::domain_error(
"cholesky_factor_free: y is not a Cholesky factor");
1492 Eigen::Matrix<T,Eigen::Dynamic,1> x((N * (N + 1)) / 2 + (M - N) * N);
1495 for (
int m = 0; m < N; ++m) {
1496 for (
int n = 0; n < m; ++n)
1499 x(pos++) =
log(y(m,m));
1502 for (
int m = N; m < M; ++m)
1503 for (
int n = 0; n < N; ++n)
1510 template <
typename T>
1511 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1515 using Eigen::Matrix;
1516 using Eigen::Dynamic;
1518 int k_choose_2 = (K * (K - 1)) / 2;
1519 if (k_choose_2 != y.size()) {
1520 std::cout <<
"k_choose_2 = " << k_choose_2 <<
" y.size()=" << y.size() << std::endl;
1521 throw std::domain_error(
"y is not a valid unconstrained cholesky correlation matrix."
1522 "Require (K choose 2) elements in y.");
1524 Matrix<T,Dynamic,1> z(k_choose_2);
1525 for (
int i = 0; i < k_choose_2; ++i)
1527 Matrix<T,Dynamic,Dynamic> x(K,K);
1528 if (K == 0)
return x;
1530 for (
int j = 1; j < K; ++j)
1531 for (
int i = 0; i < j; ++i)
1535 for (
int i = 1; i < K; ++i) {
1537 T sum_sqs(
square(x(i,0)));
1538 for (
int j = 1; j < i; ++j) {
1539 x(i,j) = z(k++) *
sqrt(1.0 - sum_sqs);
1540 sum_sqs +=
square(x(i,j));
1542 x(i,i) =
sqrt(1.0 - sum_sqs);
1548 template <
typename T>
1549 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1554 using Eigen::Matrix;
1555 using Eigen::Dynamic;
1558 int k_choose_2 = (K * (K - 1)) / 2;
1559 if (k_choose_2 != y.size()) {
1560 std::cout <<
"k_choose_2 = " << k_choose_2 <<
" y.size()=" << y.size() << std::endl;
1561 throw std::domain_error(
"y is not a valid unconstrained cholesky correlation matrix."
1562 " Require (K choose 2) elements in y.");
1564 Matrix<T,Dynamic,1> z(k_choose_2);
1565 for (
int i = 0; i < k_choose_2; ++i)
1567 Matrix<T,Dynamic,Dynamic> x(K,K);
1568 if (K == 0)
return x;
1570 for (
int j = 1; j < K; ++j)
1571 for (
int i = 0; i < j; ++i)
1575 for (
int i = 1; i < K; ++i) {
1577 T sum_sqs =
square(x(i,0));
1578 for (
int j = 1; j < i; ++j) {
1579 lp += 0.5 *
log1m(sum_sqs);
1580 x(i,j) = z(k++) *
sqrt(1.0 - sum_sqs);
1581 sum_sqs +=
square(x(i,j));
1583 x(i,i) =
sqrt(1.0 - sum_sqs);
1590 template <
typename T>
1591 Eigen::Matrix<T,Eigen::Dynamic,1>
1594 using Eigen::Matrix;
1595 using Eigen::Dynamic;
1602 int K = (x.rows() * (x.rows() - 1)) / 2;
1603 Matrix<T,Dynamic,1> z(K);
1605 for (
int i = 1; i < x.rows(); ++i) {
1607 double sum_sqs =
square(x(i,0));
1608 for (
int j = 1; j < i; ++j) {
1610 sum_sqs +=
square(x(i,j));
1642 template <
typename T>
1643 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1647 using Eigen::Dynamic;
1648 using Eigen::Matrix;
1650 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1652 size_type k_choose_2 = (k * (k - 1)) / 2;
1653 if (k_choose_2 != x.size())
1654 throw std::invalid_argument (
"x is not a valid correlation matrix");
1655 Eigen::Array<T,Eigen::Dynamic,1> cpcs(k_choose_2);
1656 for (size_type i = 0; i < k_choose_2; ++i)
1680 template <
typename T>
1681 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1686 using Eigen::Dynamic;
1687 using Eigen::Matrix;
1689 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1691 size_type k_choose_2 = (k * (k - 1)) / 2;
1692 if (k_choose_2 != x.size())
1693 throw std::invalid_argument (
"x is not a valid correlation matrix");
1694 Array<T,Dynamic,1> cpcs(k_choose_2);
1695 for (size_type i = 0; i < k_choose_2; ++i)
1720 template <
typename T>
1721 Eigen::Matrix<T,Eigen::Dynamic,1>
1724 using Eigen::Dynamic;
1725 using Eigen::Matrix;
1727 typedef typename index_type<Matrix<T,Dynamic,1> >::type
size_type;
1729 size_type k = y.rows();
1731 throw std::domain_error(
"y is not a square matrix");
1733 throw std::domain_error(
"y has no elements");
1734 size_type k_choose_2 = (k * (k-1)) / 2;
1735 Array<T,Dynamic,1> x(k_choose_2);
1736 Array<T,Dynamic,1> sds(k);
1739 throw std::runtime_error(
"factor_cov_matrix failed on y");
1740 for (size_type i = 0; i < k; ++i) {
1743 std::stringstream s;
1744 s <<
"all standard deviations must be zero."
1745 <<
" found log(sd[" << i <<
"])=" << sds[i] << std::endl;
1746 BOOST_THROW_EXCEPTION(std::runtime_error(s.str()));
1767 template <
typename T>
1768 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1773 using Eigen::Dynamic;
1774 using Eigen::Matrix;
1777 typedef typename index_type<Matrix<T,Dynamic,Dynamic> >::type
size_type;
1779 Matrix<T,Dynamic,Dynamic> L(K,K);
1780 if (x.size() != (K * (K + 1)) / 2)
1781 throw std::domain_error(
"x.size() != K + (K choose 2)");
1783 for (size_type m = 0; m < K; ++m) {
1784 for (
int n = 0; n < m; ++n)
1786 L(m,m) =
exp(x(i++));
1787 for (size_type n = m + 1; n < K; ++n)
1806 template <
typename T>
1807 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1811 Eigen::Dynamic> >::type K,
1815 using Eigen::Dynamic;
1816 using Eigen::Matrix;
1818 typedef typename index_type<Matrix<T,Dynamic,Dynamic> >::type
size_type;
1820 if (x.size() != (K * (K + 1)) / 2)
1821 throw std::domain_error(
"x.size() != K + (K choose 2)");
1822 Matrix<T,Dynamic,Dynamic> L(K,K);
1824 for (size_type m = 0; m < K; ++m) {
1825 for (size_type n = 0; n < m; ++n)
1827 L(m,m) =
exp(x(i++));
1828 for (size_type n = m + 1; n < K; ++n)
1833 for (
int k = 0; k < K; ++k)
1834 lp += (K - k + 1) *
log(L(k,k));
1835 return L * L.transpose();
1861 template <
typename T>
1862 Eigen::Matrix<T,Eigen::Dynamic,1>
1867 throw std::domain_error(
"y is not a square matrix");
1869 throw std::domain_error(
"y has no elements");
1870 for (
int k = 0; k < K; ++k)
1871 if (!(y(k,k) > 0.0))
1872 throw std::domain_error(
"y has non-positive diagonal");
1873 Eigen::Matrix<T,Eigen::Dynamic,1> x((K * (K + 1)) / 2);
1876 Eigen::LLT<Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> >
1879 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L = llt.matrixL();
1881 for (
int m = 0; m < K; ++m) {
1882 for (
int n = 0; n < m; ++n)
1884 x(i++) =
log(L(m,m));
1908 template <
typename T>
1909 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1912 size_t k_choose_2 = (k * (k - 1)) / 2;
1913 Eigen::Array<T,Eigen::Dynamic,1> cpcs(k_choose_2);
1915 for (
size_t i = 0; i < k_choose_2; ++i)
1917 Eigen::Array<T,Eigen::Dynamic,1> sds(k);
1918 for (
size_t i = 0; i < k; ++i)
1947 template <
typename T>
1948 Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1952 size_t k_choose_2 = (k * (k - 1)) / 2;
1953 Eigen::Array<T,Eigen::Dynamic,1> cpcs(k_choose_2);
1955 for (
size_t i = 0; i < k_choose_2; ++i)
1957 Eigen::Array<T,Eigen::Dynamic,1> sds(k);
1958 for (
size_t i = 0; i < k; ++i)
1981 template <
typename T>
1982 Eigen::Matrix<T,Eigen::Dynamic,1>
1984 const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& y) {
1987 using Eigen::Dynamic;
1988 using Eigen::Matrix;
1990 typedef typename index_type<Matrix<T,Dynamic,Dynamic> >::type
size_type;
1992 size_type k = y.rows();
1994 throw std::domain_error(
"y is not a square matrix");
1996 throw std::domain_error(
"y has no elements");
1997 size_type k_choose_2 = (k * (k-1)) / 2;
1998 Array<T,Dynamic,1> cpcs(k_choose_2);
1999 Array<T,Dynamic,1> sds(k);
2002 throw std::runtime_error (
"factor_cov_matrix failed on y");
2003 Matrix<T,Dynamic,1> x(k_choose_2 + k);
2005 for (size_type i = 0; i < k_choose_2; ++i)
2007 for (size_type i = 0; i < k; ++i)
fvar< T > log1m(const fvar< T > &x)
T square(const T x)
Return the square of the specified argument.
const double LOG_2
The natural logarithm of 2, .
fvar< T > fabs(const fvar< T > &x)
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...
bool check_ordered(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, 1 > &y, const char *name, T_result *result)
Return true if the specified vector is sorted into increasing order.
fvar< T > log1p_exp(const fvar< T > &x)
boost::math::tools::promote_args< T >::type log1p_exp(const T a)
Calculates the log of 1 plus the exponential of the specified value without overflow.
T corr_constrain(const T x)
Return the result of transforming the specified scalar to have a valid correlation value between -1 a...
fvar< T > sin(const fvar< T > &x)
void factor_U(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &U, Eigen::Array< T, Eigen::Dynamic, 1 > &CPCs)
This function is intended to make starting values, given a unit upper-triangular matrix U such that U...
T identity_free(const T y)
Returns the result of applying the inverse of the identity constraint transform to the input...
boost::math::tools::promote_args< T >::type logit(const T a)
Returns the logit function applied to the argument.
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.
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > cholesky_factor_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x, int M, int N)
Return the Cholesky factor of the specified size read from the specified vector.
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 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, 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 > ordered_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x)
Return an increasing ordered vector derived from the specified free vector.
bool check_simplex(const char *function, const Eigen::Matrix< T_prob, Eigen::Dynamic, 1 > &theta, const char *name, T_result *result)
Return true if the specified vector is simplex.
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, 1 > simplex_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &y)
Return the simplex corresponding to the specified free vector.
fvar< T > square(const fvar< T > &x)
boost::math::tools::promote_args< T >::type inv_logit(const T a)
Returns the inverse logit function applied to the argument.
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 prob_constrain(const T x)
Return a probability value constrained to fall between 0 and 1 (inclusive) for the specified free sca...
fvar< T > sum(const Eigen::Matrix< fvar< T >, R, C > &m)
bool check_unit_vector(const char *function, const Eigen::Matrix< T_prob, Eigen::Dynamic, 1 > &theta, const char *name, T_result *result)
Return true if the specified vector is unit vector.
fvar< T > inv_logit(const fvar< T > &x)
fvar< T > tanh(const fvar< T > &x)
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...
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, TU >::type ub_constrain(const T x, const TU ub)
Return the upper-bounded value for the specified unconstrained scalar and upper bound.
Primary template class for the metaprogram to compute the index type of a container.
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic >::Index size_type
T corr_free(const T y)
Return the unconstrained scalar that when transformed to a valid correlation produces the specified v...
T lb_constrain(const T x, const TL lb)
Return the lower-bounded value for the specified unconstrained input and specified lower bound...
fvar< T > sqrt(const fvar< T > &x)
bool check_cholesky_factor(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const char *name, T_result *result)
Return true if the specified matrix is a valid Cholesky factor.
Eigen::Matrix< T, Eigen::Dynamic, 1 > unit_vector_free(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x)
Eigen::Matrix< T, Eigen::Dynamic, 1 > cholesky_corr_free(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &x)
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 > unit_vector_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &y)
Return the unit length vector corresponding to the free vector y.
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 > 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...
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...
fvar< T > atan2(const fvar< T > &x1, const fvar< T > &x2)
fvar< T > atanh(const fvar< T > &x)
bool check_positive(const char *function, const T_y &y, const char *name, T_result *result)
T identity_constrain(T x)
Returns the result of applying the identity constraint transform to the input.
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > cov_matrix_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x, typename math::index_type< Eigen::Matrix< T, Eigen::Dynamic, 1 > >::type K)
Return the symmetric, positive-definite matrix of dimensions K by K resulting from transforming the s...
fvar< T > cos(const fvar< T > &x)
fvar< T > logit(const fvar< T > &x)
Eigen::Matrix< T, Eigen::Dynamic, 1 > cholesky_factor_free(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &y)
Return the unconstrained vector of parameters correspdonding to the specified Cholesky factor...
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > corr_matrix_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x, typename math::index_type< Eigen::Matrix< T, Eigen::Dynamic, 1 > >::type k)
Return the correlation matrix of the specified dimensionality derived from the specified vector of un...
double sum(std::vector< double > &x)
bool check_less_or_equal(const char *function, const T_y &y, const T_high &high, const char *name, T_result *result)
double e()
Return the base of the natural logarithm.
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...
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...
const double E
The base of the natural logarithm, .
Eigen::Matrix< fvar< T >, R, R > multiply_lower_tri_self_transpose(const Eigen::Matrix< fvar< T >, R, C > &m)
bool check_positive_ordered(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, 1 > &y, const char *name, T_result *result)
Return true if the specified vector contains only non-negative values and is sorted into increasing o...
bool factor_cov_matrix(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &Sigma, Eigen::Array< T, Eigen::Dynamic, 1 > &CPCs, Eigen::Array< T, Eigen::Dynamic, 1 > &sds)
This function is intended to make starting values, given a covariance matrix Sigma.
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...
T positive_free(const T y)
Return the unconstrained value corresponding to the specified positive-constrained value...
T positive_constrain(const T x)
Return the positive value for the specified unconstrained input.
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...
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...
bool check_less(const char *function, const T_y &y, const T_high &high, const char *name, T_result *result)
const double CONSTRAINT_TOLERANCE
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > cholesky_corr_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &y, int K)
fvar< T > log(const fvar< T > &x)
fvar< T > exp(const fvar< T > &x)
bool check_square(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const char *name, T_result *result)
Return true if the specified matrix is square.
bool check_greater_or_equal(const char *function, const T_y &y, const T_low &low, const char *name, T_result *result)
fvar< T > log1p(const fvar< T > &x)
boost::math::tools::promote_args< T >::type log1m(T x)
Return the natural logarithm of one minus the specified value.
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...