Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
transform.hpp
Go to the documentation of this file.
1 #ifndef STAN__PROB__TRANSFORM_HPP
2 #define STAN__PROB__TRANSFORM_HPP
3 
4 #include <cmath>
5 #include <cstddef>
6 #include <limits>
7 #include <stdexcept>
8 #include <sstream>
9 #include <vector>
10 
11 #include <boost/throw_exception.hpp>
12 #include <boost/math/tools/promotion.hpp>
13 
15 
16 #include <stan/math.hpp>
19 #include <stan/math/matrix.hpp>
20 
21 #include <stan/math/matrix/sum.hpp>
25 
26 
27 namespace stan {
28 
29  namespace prob {
30 
31 
32  const double CONSTRAINT_TOLERANCE = 1E-8;
33 
41  template<typename T>
42  void
43  factor_U(const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& U,
44  Eigen::Array<T,Eigen::Dynamic,1>& CPCs) {
45 
46  size_t K = U.rows();
47  size_t position = 0;
48  size_t pull = K - 1;
49 
50  if (K == 2) {
51  CPCs(0) = atanh(U(0,1));
52  return;
53  }
54 
55  Eigen::Array<T,1,Eigen::Dynamic> temp = U.row(0).tail(pull);
56 
57  CPCs.head(pull) = temp;
58 
59  Eigen::Array<T,Eigen::Dynamic,1> acc(K);
60  acc(0) = -0.0;
61  acc.tail(pull) = 1.0 - temp.square();
62  for(size_t i = 1; i < (K - 1); i++) {
63  position += pull;
64  pull--;
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();
69  }
70  CPCs = 0.5 * ( (1.0 + CPCs) / (1.0 - CPCs) ).log(); // now unbounded
71  }
72 
73 
86  template<typename T>
87  bool
88  factor_cov_matrix(const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& Sigma,
89  Eigen::Array<T,Eigen::Dynamic,1>& CPCs,
90  Eigen::Array<T,Eigen::Dynamic,1>& sds) {
91 
92  size_t K = sds.rows();
93 
94  sds = Sigma.diagonal().array();
95  if( (sds <= 0.0).any() ) return false;
96  sds = sds.sqrt();
97 
98  Eigen::DiagonalMatrix<T,Eigen::Dynamic> D(K);
99  D.diagonal() = sds.inverse();
100  sds = sds.log(); // now unbounded
101 
102  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> R = D * Sigma * D;
103  // to hopefully prevent pivoting due to floating point error
104  R.diagonal().setOnes();
105  Eigen::LDLT<Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> > ldlt;
106  ldlt = R.ldlt();
107  if (!ldlt.isPositive())
108  return false;
109  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> U = ldlt.matrixU();
110  factor_U(U, CPCs);
111  return true;
112  }
113 
114  // MATRIX TRANSFORMS +/- JACOBIANS
115 
137  template <typename T>
138  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
139  read_corr_L(const Eigen::Array<T,Eigen::Dynamic,1>& CPCs, // on (-1,1)
140  const size_t K) {
141  Eigen::Array<T,Eigen::Dynamic,1> temp;
142  Eigen::Array<T,Eigen::Dynamic,1> acc(K-1);
143  acc.setOnes();
144  // Cholesky factor of correlation matrix
145  Eigen::Array<T,Eigen::Dynamic,Eigen::Dynamic> L(K,K);
146  L.setZero();
147 
148  size_t position = 0;
149  size_t pull = K - 1;
150 
151  L(0,0) = 1.0;
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++) {
155  position += pull;
156  pull--;
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();
161  }
162  L(K-1,K-1) = sqrt(acc(K-2));
163  return L.matrix();
164  }
165 
179  template <typename T>
180  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
181  read_corr_matrix(const Eigen::Array<T,Eigen::Dynamic,1>& CPCs,
182  const size_t K) {
183  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L
184  = read_corr_L(CPCs, K);
187  }
188 
213  template <typename T>
214  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
215  read_corr_L(const Eigen::Array<T,Eigen::Dynamic,1>& CPCs,
216  const size_t K,
217  T& log_prob) {
218 
219  using stan::math::log1m;
220  using stan::math::square;
221  using stan::math::sum;
222 
223  Eigen::Matrix<T,Eigen::Dynamic,1> values(CPCs.rows() - 1);
224  size_t pos = 0;
225  // no need to abs() because this Jacobian determinant
226  // is strictly positive (and triangular)
227  // see inverse of Jacobian in equation 11 of LKJ paper
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)));
231  pos++;
232  }
233 
234  log_prob += 0.5 * sum(values);
235  return read_corr_L(CPCs,K);
236  }
237 
256  template <typename T>
257  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
258  read_corr_matrix(const Eigen::Array<T,Eigen::Dynamic,1>& CPCs,
259  const size_t K,
260  T& log_prob) {
261 
262  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L
263  = read_corr_L(CPCs, K, log_prob);
266  }
267 
278  template <typename T>
279  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
280  read_cov_L(const Eigen::Array<T,Eigen::Dynamic,1>& CPCs,
281  const Eigen::Array<T,Eigen::Dynamic,1>& sds,
282  T& log_prob) {
283  size_t K = sds.rows();
284  // adjust due to transformation from correlations to covariances
285  log_prob += (sds.log().sum() + stan::math::LOG_2) * K;
286  return sds.matrix().asDiagonal() * read_corr_L(CPCs, K, log_prob);
287  }
288 
298  template <typename T>
299  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
300  read_cov_matrix(const Eigen::Array<T,Eigen::Dynamic,1>& CPCs,
301  const Eigen::Array<T,Eigen::Dynamic,1>& sds,
302  T& log_prob) {
303 
304  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L
305  = read_cov_L(CPCs, sds, log_prob);
308  }
309 
317  template<typename T>
318  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
319  read_cov_matrix(const Eigen::Array<T,Eigen::Dynamic,1>& CPCs,
320  const Eigen::Array<T,Eigen::Dynamic,1>& sds) {
321 
322  size_t K = sds.rows();
323  Eigen::DiagonalMatrix<T,Eigen::Dynamic> D(K);
324  D.diagonal() = sds;
325  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L
326  = D * read_corr_L(CPCs, K);
329  }
330 
331 
341  template<typename T>
342  const Eigen::Array<T,Eigen::Dynamic,1>
343  make_nu(const T eta, const size_t K) {
344  using Eigen::Array;
345  using Eigen::Dynamic;
346  using Eigen::Matrix;
348  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
349 
350  Array<T,Dynamic,1> nu(K * (K - 1) / 2);
351 
352  T alpha = eta + (K - 2.0) / 2.0; // from Lewandowski et. al.
353 
354  // Best (1978) implies nu = 2 * alpha for the dof in a t
355  // distribution that generates a beta variate on (-1,1)
356  T alpha2 = 2.0 * alpha;
357  for (size_type j = 0; j < (K - 1); j++) {
358  nu(j) = alpha2;
359  }
360  size_t counter = K - 1;
361  for (size_type i = 1; i < (K - 1); i++) {
362  alpha -= 0.5;
363  alpha2 = 2.0 * alpha;
364  for (size_type j = i + 1; j < K; j++) {
365  nu(counter) = alpha2;
366  counter++;
367  }
368  }
369  return nu;
370  }
371 
372 
373 
374 
375  // IDENTITY
376 
389  template <typename T>
390  inline
392  return x;
393  }
394 
408  template <typename T>
409  inline
410  T identity_constrain(const T x, T& /*lp*/) {
411  return x;
412  }
413 
425  template <typename T>
426  inline
427  T identity_free(const T y) {
428  return y;
429  }
430 
431 
432  // POSITIVE
433 
444  template <typename T>
445  inline
446  T positive_constrain(const T x) {
447  return exp(x);
448  }
449 
466  template <typename T>
467  inline
468  T positive_constrain(const T x, T& lp) {
469  lp += x;
470  return exp(x);
471  }
472 
489  template <typename T>
490  inline
491  T positive_free(const T y) {
492  stan::math::check_positive("stan::prob::positive_free(%1%)", y,
493  "Positive variable", (double*)0);
494  return log(y);
495  }
496 
497  // LOWER BOUND
498 
518  template <typename T, typename TL>
519  inline
520  T lb_constrain(const T x, const TL lb) {
521  if (lb == -std::numeric_limits<double>::infinity())
522  return identity_constrain(x);
523  return exp(x) + lb;
524  }
525 
542  template <typename T, typename TL>
543  inline
544  typename boost::math::tools::promote_args<T,TL>::type
545  lb_constrain(const T x, const TL lb, T& lp) {
546  if (lb == -std::numeric_limits<double>::infinity())
547  return identity_constrain(x,lp);
548  lp += x;
549  return exp(x) + lb;
550  }
551 
567  template <typename T, typename TL>
568  inline
569  typename boost::math::tools::promote_args<T,TL>::type
570  lb_free(const T y, const TL lb) {
571  if (lb == -std::numeric_limits<double>::infinity())
572  return identity_free(y);
573  stan::math::check_greater_or_equal("stan::prob::lb_free(%1%)",
574  y, lb, "Lower bounded variable",
575  (double*)0);
576  return log(y - lb);
577  }
578 
579 
580  // UPPER BOUND
581 
601  template <typename T, typename TU>
602  inline
603  typename boost::math::tools::promote_args<T,TU>::type
604  ub_constrain(const T x, const TU ub) {
605  if (ub == std::numeric_limits<double>::infinity())
606  return identity_constrain(x);
607  return ub - exp(x);
608  }
609 
633  template <typename T, typename TU>
634  inline
635  typename boost::math::tools::promote_args<T,TU>::type
636  ub_constrain(const T x, const TU ub, T& lp) {
637  if (ub == std::numeric_limits<double>::infinity())
638  return identity_constrain(x,lp);
639  lp += x;
640  return ub - exp(x);
641  }
642 
665  template <typename T, typename TU>
666  inline
667  typename boost::math::tools::promote_args<T,TU>::type
668  ub_free(const T y, const TU ub) {
669  if (ub == std::numeric_limits<double>::infinity())
670  return identity_free(y);
671  stan::math::check_less_or_equal("stan::prob::ub_free(%1%)",
672  y, ub, "Upper bounded variable",
673  (double*)0);
674  return log(ub - y);
675  }
676 
677 
678  // LOWER & UPPER BOUNDS
679 
707  template <typename T, typename TL, typename TU>
708  inline
709  typename boost::math::tools::promote_args<T,TL,TU>::type
710  lub_constrain(const T x, TL lb, TU ub) {
711  stan::math::check_less("lub_constrain(%1%)",lb,ub,"lb",(double*)0);
712 
713  if (lb == -std::numeric_limits<double>::infinity())
714  return ub_constrain(x,ub);
715  if (ub == std::numeric_limits<double>::infinity())
716  return lb_constrain(x,lb);
717 
718  T inv_logit_x;
719  if (x > 0) {
720  T exp_minus_x = exp(-x);
721  inv_logit_x = 1.0 / (1.0 + exp_minus_x);
722  // Prevent x from reaching one unless it really really should.
723  if ((x < std::numeric_limits<double>::infinity())
724  && (inv_logit_x == 1))
725  inv_logit_x = 1 - 1e-15;
726  } else {
727  T exp_x = exp(x);
728  inv_logit_x = 1.0 - 1.0 / (1.0 + exp_x);
729  // Prevent x from reaching zero unless it really really should.
730  if ((x > -std::numeric_limits<double>::infinity())
731  && (inv_logit_x== 0))
732  inv_logit_x = 1e-15;
733  }
734  return lb + (ub - lb) * inv_logit_x;
735  }
736 
778  template <typename T, typename TL, typename TU>
779  typename boost::math::tools::promote_args<T,TL,TU>::type
780  lub_constrain(const T x, const TL lb, const TU ub, T& lp) {
781  if (!(lb < ub)) {
782  std::stringstream s;
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());
786  }
787  if (lb == -std::numeric_limits<double>::infinity())
788  return ub_constrain(x,ub,lp);
789  if (ub == std::numeric_limits<double>::infinity())
790  return lb_constrain(x,lb,lp);
791  T inv_logit_x;
792  if (x > 0) {
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);
796  // Prevent x from reaching one unless it really really should.
797  if ((x < std::numeric_limits<double>::infinity())
798  && (inv_logit_x == 1))
799  inv_logit_x = 1 - 1e-15;
800  } else {
801  T exp_x = exp(x);
802  inv_logit_x = 1.0 - 1.0 / (1.0 + exp_x);
803  lp += log(ub - lb) + x - 2 * log1p(exp_x);
804  // Prevent x from reaching zero unless it really really should.
805  if ((x > -std::numeric_limits<double>::infinity())
806  && (inv_logit_x== 0))
807  inv_logit_x = 1e-15;
808  }
809  return lb + (ub - lb) * inv_logit_x;
810  }
811 
842  template <typename T, typename TL, typename TU>
843  inline
844  typename boost::math::tools::promote_args<T,TL,TU>::type
845  lub_free(const T y, TL lb, TU ub) {
846  using stan::math::logit;
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())
851  return ub_free(y,ub);
852  if (ub == std::numeric_limits<double>::infinity())
853  return lb_free(y,lb);
854  return logit((y - lb) / (ub - lb));
855  }
856 
857 
858  // PROBABILITY
859 
873  template <typename T>
874  inline
875  T prob_constrain(const T x) {
876  using stan::math::inv_logit;
877  return inv_logit(x);
878  }
879 
901  template <typename T>
902  inline
903  T prob_constrain(const T x, T& lp) {
904  using stan::math::inv_logit;
905  using stan::math::log1m;
906  T inv_logit_x = inv_logit(x);
907  lp += log(inv_logit_x) + log1m(inv_logit_x);
908  return inv_logit_x;
909  }
910 
925  template <typename T>
926  inline
927  T prob_free(const T y) {
928  using stan::math::logit;
929  stan::math::check_bounded<T,double,double,T>
930  ("stan::prob::prob_free(%1%)",
931  y, 0, 1, "Probability variable",(double*)0);
932  return logit(y);
933  }
934 
935 
936  // CORRELATION
937 
950  template <typename T>
951  inline
952  T corr_constrain(const T x) {
953  return tanh(x);
954  }
955 
968  template <typename T>
969  inline
970  T corr_constrain(const T x, T& lp) {
971  using stan::math::log1m;
972  T tanh_x = tanh(x);
973  lp += log1m(tanh_x * tanh_x);
974  return tanh_x;
975  }
976 
993  template <typename T>
994  inline
995  T corr_free(const T y) {
996  stan::math::check_bounded<T,double,double,double>
997  ("stan::prob::lub_free(%1%)",
998  y, -1, 1, "Correlation variable", (double*)0);
999  return atanh(y);
1000  }
1001 
1002 
1003  // Unit vector
1004 
1013  template <typename T>
1014  Eigen::Matrix<T,Eigen::Dynamic,1>
1015  unit_vector_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& y) {
1016  using Eigen::Matrix;
1017  using Eigen::Dynamic;
1018  using stan::math::index_type;
1019  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1020  int Km1 = y.size();
1021  Matrix<T,Dynamic,1> x(Km1 + 1);
1022  x(0) = 1.0;
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);
1029  }
1030  return x;
1031  }
1032 
1042  template <typename T>
1043  Eigen::Matrix<T,Eigen::Dynamic,1>
1044  unit_vector_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& y, T &lp) {
1045  using Eigen::Matrix;
1046  using Eigen::Dynamic;
1047  using stan::math::index_type;
1048  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1049 
1050  int Km1 = y.size();
1051  Matrix<T,Dynamic,1> x(Km1 + 1);
1052  x(0) = 1.0;
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);
1059  if (k < Km1)
1060  lp += (Km1 - k) * log(fabs(sin_yk_1));
1061  }
1062  return x;
1063  }
1064 
1065  template <typename T>
1066  Eigen::Matrix<T,Eigen::Dynamic,1>
1067  unit_vector_free(const Eigen::Matrix<T,Eigen::Dynamic,1>& x) {
1068  using Eigen::Matrix;
1069  using Eigen::Dynamic;
1070  using stan::math::index_type;
1071  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1072 
1073  stan::math::check_unit_vector("stan::prob::unit_vector_free(%1%)",
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;
1081  sumSq += x(k)*x(k);
1082  }
1083  return y;
1084  }
1085 
1086 
1087  // SIMPLEX
1088 
1101  template <typename T>
1102  Eigen::Matrix<T,Eigen::Dynamic,1>
1103  simplex_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& y) {
1104 
1105  // cut & paste simplex_constrain(Eigen::Matrix,T) w/o Jacobian
1106  using Eigen::Matrix;
1107  using Eigen::Dynamic;
1108  using stan::math::index_type;
1109  using stan::math::inv_logit;
1110  using stan::math::logit;
1111  using stan::math::log1m;
1112  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1113 
1114 
1115  int Km1 = y.size();
1116  Matrix<T,Dynamic,1> x(Km1 + 1);
1117  T stick_len(1.0);
1118  for (size_type k = 0; k < Km1; ++k) {
1119  T z_k(inv_logit(y(k) - log(Km1 - k)));
1120  x(k) = stick_len * z_k;
1121  stick_len -= x(k);
1122  }
1123  x(Km1) = stick_len;
1124  return x;
1125  }
1126 
1140  template <typename T>
1141  Eigen::Matrix<T,Eigen::Dynamic,1>
1142  simplex_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& y,
1143  T& lp) {
1144 
1145  using Eigen::Dynamic;
1146  using Eigen::Matrix;
1147  using stan::math::index_type;
1148  using stan::math::inv_logit;
1149  using stan::math::logit;
1150  using stan::math::log1m;
1151  using stan::math::log1p_exp;
1152 
1153  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1154 
1155  int Km1 = y.size(); // K = Km1 + 1
1156  Matrix<T,Dynamic,1> x(Km1 + 1);
1157  T stick_len(1.0);
1158  for (size_type k = 0; k < Km1; ++k) {
1159  double eq_share = -log(Km1 - k); // = logit(1.0/(Km1 + 1 - k));
1160  T adj_y_k(y(k) + eq_share);
1161  T z_k(inv_logit(adj_y_k));
1162  x(k) = stick_len * z_k;
1163  lp += log(stick_len);
1164  lp -= log1p_exp(-adj_y_k);
1165  lp -= log1p_exp(adj_y_k);
1166  stick_len -= x(k); // equivalently *= (1 - z_k);
1167  }
1168  x(Km1) = stick_len; // no Jacobian contrib for last dim
1169  return x;
1170  }
1171 
1186  template <typename T>
1187  Eigen::Matrix<T,Eigen::Dynamic,1>
1188  simplex_free(const Eigen::Matrix<T,Eigen::Dynamic,1>& x) {
1189  using Eigen::Dynamic;
1190  using Eigen::Matrix;
1191  using stan::math::index_type;
1192  using stan::math::logit;
1193 
1194  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1195 
1196  stan::math::check_simplex("stan::prob::simplex_free(%1%)", x, "Simplex variable",
1197  (double*)0);
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; ) {
1202  stick_len += x(k);
1203  T z_k(x(k) / stick_len);
1204  y(k) = logit(z_k) + log(Km1 - k);
1205  // note: log(Km1 - k) = logit(1.0 / (Km1 + 1 - k));
1206  }
1207  return y;
1208  }
1209 
1210 
1211  // ORDERED
1212 
1222  template <typename T>
1223  Eigen::Matrix<T,Eigen::Dynamic,1>
1224  ordered_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& x) {
1225  using Eigen::Matrix;
1226  using Eigen::Dynamic;
1227  using stan::math::index_type;
1228 
1229  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1230 
1231  size_type k = x.size();
1232  Matrix<T,Dynamic,1> y(k);
1233  if (k == 0)
1234  return y;
1235  y[0] = x[0];
1236  for (size_type i = 1; i < k; ++i)
1237  y[i] = y[i-1] + exp(x[i]);
1238  return y;
1239  }
1240 
1253  template <typename T>
1254  inline
1255  Eigen::Matrix<T,Eigen::Dynamic,1>
1256  ordered_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& x, T& lp) {
1257  using Eigen::Matrix;
1258  using Eigen::Dynamic;
1259  using stan::math::index_type;
1260 
1261  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1262 
1263  for (size_type i = 1; i < x.size(); ++i)
1264  lp += x(i);
1265  return ordered_constrain(x);
1266  }
1267 
1268 
1269 
1283  template <typename T>
1284  Eigen::Matrix<T,Eigen::Dynamic,1>
1285  ordered_free(const Eigen::Matrix<T,Eigen::Dynamic,1>& y) {
1286  stan::math::check_ordered("stan::prob::ordered_free(%1%)",
1287  y, "Ordered variable");
1288  using Eigen::Matrix;
1289  using Eigen::Dynamic;
1290  using stan::math::index_type;
1291  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1292 
1293  size_type k = y.size();
1294  Matrix<T,Dynamic,1> x(k);
1295  if (k == 0)
1296  return x;
1297  x[0] = y[0];
1298  for (size_type i = 1; i < k; ++i)
1299  x[i] = log(y[i] - y[i-1]);
1300  return x;
1301  }
1302 
1303 
1304  // POSITIVE ORDERED
1305 
1315  template <typename T>
1316  Eigen::Matrix<T,Eigen::Dynamic,1>
1317  positive_ordered_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& x) {
1318  using Eigen::Matrix;
1319  using Eigen::Dynamic;
1320  using stan::math::index_type;
1321  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1322 
1323  size_type k = x.size();
1324  Matrix<T,Dynamic,1> y(k);
1325  if (k == 0)
1326  return y;
1327  y[0] = exp(x[0]);
1328  for (size_type i = 1; i < k; ++i)
1329  y[i] = y[i-1] + exp(x[i]);
1330  return y;
1331  }
1332 
1345  template <typename T>
1346  inline
1347  Eigen::Matrix<T,Eigen::Dynamic,1>
1348  positive_ordered_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& x, T& lp) {
1349  using Eigen::Matrix;
1350  using Eigen::Dynamic;
1351  using stan::math::index_type;
1352  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1353 
1354  for (size_type i = 0; i < x.size(); ++i)
1355  lp += x(i);
1356  return positive_ordered_constrain(x);
1357  }
1358 
1359 
1360 
1374  template <typename T>
1375  Eigen::Matrix<T,Eigen::Dynamic,1>
1376  positive_ordered_free(const Eigen::Matrix<T,Eigen::Dynamic,1>& y) {
1377  using Eigen::Matrix;
1378  using Eigen::Dynamic;
1379  using stan::math::index_type;
1380 
1381  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1382 
1383  stan::math::check_positive_ordered("stan::prob::positive_ordered_free(%1%)",
1384  y, "Positive ordered variable", (double*)0);
1385 
1386  size_type k = y.size();
1387  Matrix<T,Dynamic,1> x(k);
1388  if (k == 0)
1389  return x;
1390  x[0] = log(y[0]);
1391  for (size_type i = 1; i < k; ++i)
1392  x[i] = log(y[i] - y[i-1]);
1393  return x;
1394  }
1395 
1396 
1397  // CHOLESKY FACTOR
1398 
1410  template <typename T>
1411  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1412  cholesky_factor_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
1413  int M,
1414  int N) {
1415  using std::exp;
1416  if (M < N)
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);
1422  T zero(0);
1423  int pos = 0;
1424  // upper square
1425  for (int m = 0; m < N; ++m) {
1426  for (int n = 0; n < m; ++n)
1427  y(m,n) = x(pos++);
1428  y(m,m) = exp(x(pos++));
1429  for (int n = m + 1; n < N; ++n)
1430  y(m,n) = zero;
1431  }
1432  // lower rectangle
1433  for (int m = N; m < M; ++m)
1434  for (int n = 0; n < N; ++n)
1435  y(m,n) = x(pos++);
1436  return y;
1437  }
1438 
1453  template <typename T>
1454  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1455  cholesky_factor_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
1456  int M,
1457  int N,
1458  T& lp) {
1459  // cut-and-paste from above, so checks twice
1460 
1461  using stan::math::sum;
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");
1464  int pos = 0;
1465  std::vector<T> log_jacobians(N);
1466  for (int n = 0; n < N; ++n) {
1467  pos += n;
1468  log_jacobians[n] = x(pos++);
1469  }
1470  lp += sum(log_jacobians); // optimized for autodiff vs. direct lp +=
1471  return cholesky_factor_constrain(x,M,N);
1472  }
1473 
1483  template <typename T>
1484  Eigen::Matrix<T,Eigen::Dynamic,1>
1485  cholesky_factor_free(const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& y) {
1486  using std::log;
1487  double result;
1488  if (!stan::math::check_cholesky_factor("cholesky_factor_free(%1%)",y,"y",&result))
1489  throw std::domain_error("cholesky_factor_free: y is not a Cholesky factor");
1490  int M = y.rows();
1491  int N = y.cols();
1492  Eigen::Matrix<T,Eigen::Dynamic,1> x((N * (N + 1)) / 2 + (M - N) * N);
1493  int pos = 0;
1494  // lower triangle of upper square
1495  for (int m = 0; m < N; ++m) {
1496  for (int n = 0; n < m; ++n)
1497  x(pos++) = y(m,n);
1498  // diagonal of upper square
1499  x(pos++) = log(y(m,m));
1500  }
1501  // lower rectangle
1502  for (int m = N; m < M; ++m)
1503  for (int n = 0; n < N; ++n)
1504  x(pos++) = y(m,n);
1505  return x;
1506  }
1507 
1508  // CHOLESKY CORRELATION MATRIX
1509 
1510  template <typename T>
1511  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1512  cholesky_corr_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& y,
1513  int K) {
1514  using std::sqrt;
1515  using Eigen::Matrix;
1516  using Eigen::Dynamic;
1517  using stan::math::square;
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.");
1523  }
1524  Matrix<T,Dynamic,1> z(k_choose_2);
1525  for (int i = 0; i < k_choose_2; ++i)
1526  z(i) = corr_constrain(y(i));
1527  Matrix<T,Dynamic,Dynamic> x(K,K);
1528  if (K == 0) return x;
1529  T zero(0);
1530  for (int j = 1; j < K; ++j)
1531  for (int i = 0; i < j; ++i)
1532  x(i,j) = zero;
1533  x(0,0) = 1;
1534  int k = 0;
1535  for (int i = 1; i < K; ++i) {
1536  x(i,0) = z(k++);
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));
1541  }
1542  x(i,i) = sqrt(1.0 - sum_sqs);
1543  }
1544  return x;
1545  }
1546 
1547  // FIXME to match above after debugged
1548  template <typename T>
1549  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1550  cholesky_corr_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& y,
1551  int K,
1552  T& lp) {
1553  using std::sqrt;
1554  using Eigen::Matrix;
1555  using Eigen::Dynamic;
1556  using stan::math::log1m;
1557  using stan::math::square;
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.");
1563  }
1564  Matrix<T,Dynamic,1> z(k_choose_2);
1565  for (int i = 0; i < k_choose_2; ++i)
1566  z(i) = corr_constrain(y(i),lp);
1567  Matrix<T,Dynamic,Dynamic> x(K,K);
1568  if (K == 0) return x;
1569  T zero(0);
1570  for (int j = 1; j < K; ++j)
1571  for (int i = 0; i < j; ++i)
1572  x(i,j) = zero;
1573  x(0,0) = 1;
1574  int k = 0;
1575  for (int i = 1; i < K; ++i) {
1576  x(i,0) = z(k++);
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));
1582  }
1583  x(i,i) = sqrt(1.0 - sum_sqs);
1584  }
1585  return x;
1586  }
1587 
1588 
1589 
1590  template <typename T>
1591  Eigen::Matrix<T,Eigen::Dynamic,1>
1592  cholesky_corr_free(const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& x) {
1593  using std::sqrt;
1594  using Eigen::Matrix;
1595  using Eigen::Dynamic;
1596  using stan::math::square;
1597 
1598  double bad_result; // won't be used
1599  stan::math::check_square("cholesky_corr_free",x,"x",&bad_result);
1600  // should validate lower-triangular, unit lengths
1601 
1602  int K = (x.rows() * (x.rows() - 1)) / 2;
1603  Matrix<T,Dynamic,1> z(K);
1604  int k = 0;
1605  for (int i = 1; i < x.rows(); ++i) {
1606  z(k++) = corr_free(x(i,0));
1607  double sum_sqs = square(x(i,0));
1608  for (int j = 1; j < i; ++j) {
1609  z(k++) = corr_free(x(i,j) / sqrt(1.0 - sum_sqs));
1610  sum_sqs += square(x(i,j));
1611  }
1612  }
1613  return z;
1614  }
1615 
1616  // CORRELATION MATRIX
1617 
1642  template <typename T>
1643  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1644  corr_matrix_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
1645  typename math::index_type<Eigen::Matrix<T,Eigen::Dynamic,1> >::type k) {
1646 
1647  using Eigen::Dynamic;
1648  using Eigen::Matrix;
1649  using stan::math::index_type;
1650  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1651 
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)
1657  cpcs[i] = corr_constrain(x[i]);
1658  return read_corr_matrix(cpcs,k);
1659  }
1660 
1680  template <typename T>
1681  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1682  corr_matrix_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
1683  typename math::index_type<Eigen::Matrix<T,Eigen::Dynamic,1> >::type k,
1684  T& lp) {
1685  using Eigen::Array;
1686  using Eigen::Dynamic;
1687  using Eigen::Matrix;
1688  using stan::math::index_type;
1689  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1690 
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)
1696  cpcs[i] = corr_constrain(x[i],lp);
1697  return read_corr_matrix(cpcs,k,lp);
1698  }
1699 
1720  template <typename T>
1721  Eigen::Matrix<T,Eigen::Dynamic,1>
1722  corr_matrix_free(const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& y) {
1723  using Eigen::Array;
1724  using Eigen::Dynamic;
1725  using Eigen::Matrix;
1726  using stan::math::index_type;
1727  typedef typename index_type<Matrix<T,Dynamic,1> >::type size_type;
1728 
1729  size_type k = y.rows();
1730  if (y.cols() != k)
1731  throw std::domain_error("y is not a square matrix");
1732  if (k == 0)
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);
1737  bool successful = factor_cov_matrix(y,x,sds);
1738  if (!successful)
1739  throw std::runtime_error("factor_cov_matrix failed on y");
1740  for (size_type i = 0; i < k; ++i) {
1741  // sds on log scale unconstrained
1742  if (fabs(sds[i] - 0.0) >= CONSTRAINT_TOLERANCE) {
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()));
1747  }
1748  }
1749  return x.matrix();
1750  }
1751 
1752 
1753  // COVARIANCE MATRIX
1754 
1767  template <typename T>
1768  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1769  cov_matrix_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
1770  typename math::index_type<Eigen::Matrix<T,Eigen::Dynamic,1> >::type K) {
1771  using std::exp;
1772 
1773  using Eigen::Dynamic;
1774  using Eigen::Matrix;
1775  using stan::math::index_type;
1777  typedef typename index_type<Matrix<T,Dynamic,Dynamic> >::type size_type;
1778 
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)");
1782  int i = 0;
1783  for (size_type m = 0; m < K; ++m) {
1784  for (int n = 0; n < m; ++n)
1785  L(m,n) = x(i++);
1786  L(m,m) = exp(x(i++));
1787  for (size_type n = m + 1; n < K; ++n)
1788  L(m,n) = 0.0;
1789  }
1791  }
1792 
1793 
1806  template <typename T>
1807  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1808  cov_matrix_constrain(const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
1809  typename math::index_type<Eigen::Matrix<T,
1810  Eigen::Dynamic,
1811  Eigen::Dynamic> >::type K,
1812  T& lp) {
1813  using std::exp;
1814 
1815  using Eigen::Dynamic;
1816  using Eigen::Matrix;
1817  using stan::math::index_type;
1818  typedef typename index_type<Matrix<T,Dynamic,Dynamic> >::type size_type;
1819 
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);
1823  int i = 0;
1824  for (size_type m = 0; m < K; ++m) {
1825  for (size_type n = 0; n < m; ++n)
1826  L(m,n) = x(i++);
1827  L(m,m) = exp(x(i++));
1828  for (size_type n = m + 1; n < K; ++n)
1829  L(m,n) = 0.0;
1830  }
1831  // Jacobian for complete transform, including exp() above
1832  lp += (K * stan::math::LOG_2); // needless constant; want propto
1833  for (int k = 0; k < K; ++k)
1834  lp += (K - k + 1) * log(L(k,k)); // only +1 because index from 0
1835  return L * L.transpose();
1836  // return tri_multiply_transpose(L);
1837  }
1838 
1861  template <typename T>
1862  Eigen::Matrix<T,Eigen::Dynamic,1>
1863  cov_matrix_free(const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& y) {
1864  using std::log;
1865  int K = y.rows();
1866  if (y.cols() != K)
1867  throw std::domain_error("y is not a square matrix");
1868  if (K == 0)
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);
1874  // FIXME: see Eigen LDLT for rank-revealing version -- use that
1875  // even if less efficient?
1876  Eigen::LLT<Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> >
1877  llt(y.rows());
1878  llt.compute(y);
1879  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> L = llt.matrixL();
1880  int i = 0;
1881  for (int m = 0; m < K; ++m) {
1882  for (int n = 0; n < m; ++n)
1883  x(i++) = L(m,n);
1884  x(i++) = log(L(m,m));
1885  }
1886  return x;
1887  }
1888 
1908  template <typename T>
1909  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1910  cov_matrix_constrain_lkj(const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
1911  size_t k) {
1912  size_t k_choose_2 = (k * (k - 1)) / 2;
1913  Eigen::Array<T,Eigen::Dynamic,1> cpcs(k_choose_2);
1914  int pos = 0;
1915  for (size_t i = 0; i < k_choose_2; ++i)
1916  cpcs[i] = corr_constrain(x[pos++]);
1917  Eigen::Array<T,Eigen::Dynamic,1> sds(k);
1918  for (size_t i = 0; i < k; ++i)
1919  sds[i] = positive_constrain(x[pos++]);
1920  return read_cov_matrix(cpcs, sds);
1921  }
1922 
1947  template <typename T>
1948  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>
1949  cov_matrix_constrain_lkj(const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
1950  size_t k,
1951  T& lp) {
1952  size_t k_choose_2 = (k * (k - 1)) / 2;
1953  Eigen::Array<T,Eigen::Dynamic,1> cpcs(k_choose_2);
1954  int pos = 0;
1955  for (size_t i = 0; i < k_choose_2; ++i)
1956  cpcs[i] = corr_constrain(x[pos++], lp);
1957  Eigen::Array<T,Eigen::Dynamic,1> sds(k);
1958  for (size_t i = 0; i < k; ++i)
1959  sds[i] = positive_constrain(x[pos++],lp);
1960  return read_cov_matrix(cpcs, sds, lp);
1961  }
1962 
1981  template <typename T>
1982  Eigen::Matrix<T,Eigen::Dynamic,1>
1984  const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& y) {
1985 
1986  using Eigen::Array;
1987  using Eigen::Dynamic;
1988  using Eigen::Matrix;
1989  using stan::math::index_type;
1990  typedef typename index_type<Matrix<T,Dynamic,Dynamic> >::type size_type;
1991 
1992  size_type k = y.rows();
1993  if (y.cols() != k)
1994  throw std::domain_error("y is not a square matrix");
1995  if (k == 0)
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);
2000  bool successful = factor_cov_matrix(y,cpcs,sds);
2001  if (!successful)
2002  throw std::runtime_error ("factor_cov_matrix failed on y");
2003  Matrix<T,Dynamic,1> x(k_choose_2 + k);
2004  size_type pos = 0;
2005  for (size_type i = 0; i < k_choose_2; ++i)
2006  x[pos++] = cpcs[i];
2007  for (size_type i = 0; i < k; ++i)
2008  x[pos++] = sds[i];
2009  return x;
2010  }
2011 
2012  }
2013 
2014 }
2015 
2016 #endif
fvar< T > log1m(const fvar< T > &x)
Definition: log1m.hpp:16
T square(const T x)
Return the square of the specified argument.
Definition: square.hpp:22
const double LOG_2
The natural logarithm of 2, .
Definition: constants.hpp:32
fvar< T > fabs(const fvar< T > &x)
Definition: fabs.hpp:18
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...
Definition: transform.hpp:1376
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)
Definition: log1p_exp.hpp:15
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.
Definition: log1p_exp.hpp:44
T corr_constrain(const T x)
Return the result of transforming the specified scalar to have a valid correlation value between -1 a...
Definition: transform.hpp:952
fvar< T > sin(const fvar< T > &x)
Definition: sin.hpp:14
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...
Definition: transform.hpp:43
T identity_free(const T y)
Returns the result of applying the inverse of the identity constraint transform to the input...
Definition: transform.hpp:427
boost::math::tools::promote_args< T >::type logit(const T a)
Returns the logit function applied to the argument.
Definition: logit.hpp:44
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.
Definition: transform.hpp:1188
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.
Definition: transform.hpp:1412
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 ...
Definition: transform.hpp:1910
T prob_free(const T y)
Return the free scalar that when transformed to a probability produces the specified scalar...
Definition: transform.hpp:927
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...
Definition: transform.hpp:668
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.
Definition: transform.hpp:1224
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 ...
Definition: transform.hpp:1722
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.
Definition: transform.hpp:1103
fvar< T > square(const fvar< T > &x)
Definition: square.hpp:15
boost::math::tools::promote_args< T >::type inv_logit(const T a)
Returns the inverse logit function applied to the argument.
Definition: inv_logit.hpp:52
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...
Definition: transform.hpp:300
T prob_constrain(const T x)
Return a probability value constrained to fall between 0 and 1 (inclusive) for the specified free sca...
Definition: transform.hpp:875
fvar< T > sum(const Eigen::Matrix< fvar< T >, R, C > &m)
Definition: sum.hpp:14
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)
Definition: inv_logit.hpp:15
fvar< T > tanh(const fvar< T > &x)
Definition: tanh.hpp:14
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...
Definition: transform.hpp:1983
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 ...
Definition: transform.hpp:343
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.
Definition: transform.hpp:604
Primary template class for the metaprogram to compute the index type of a container.
Definition: index_type.hpp:21
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic >::Index size_type
Definition: typedefs.hpp:14
T corr_free(const T y)
Return the unconstrained scalar that when transformed to a valid correlation produces the specified v...
Definition: transform.hpp:995
T lb_constrain(const T x, const TL lb)
Return the lower-bounded value for the specified unconstrained input and specified lower bound...
Definition: transform.hpp:520
fvar< T > sqrt(const fvar< T > &x)
Definition: sqrt.hpp:15
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)
Definition: transform.hpp:1067
Eigen::Matrix< T, Eigen::Dynamic, 1 > cholesky_corr_free(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &x)
Definition: transform.hpp:1592
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.
Definition: transform.hpp:1317
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.
Definition: transform.hpp:1015
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...
Definition: transform.hpp:845
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...
Definition: transform.hpp:1863
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...
Definition: transform.hpp:139
fvar< T > atan2(const fvar< T > &x1, const fvar< T > &x2)
Definition: atan2.hpp:16
fvar< T > atanh(const fvar< T > &x)
Definition: atanh.hpp:16
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.
Definition: transform.hpp:391
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...
Definition: transform.hpp:1769
fvar< T > cos(const fvar< T > &x)
Definition: cos.hpp:15
fvar< T > logit(const fvar< T > &x)
Definition: logit.hpp:17
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...
Definition: transform.hpp:1485
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...
Definition: transform.hpp:1644
double sum(std::vector< double > &x)
Definition: sum.hpp:10
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.
Definition: constants.hpp:86
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...
Definition: transform.hpp:181
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...
Definition: transform.hpp:710
const double E
The base of the natural logarithm, .
Definition: constants.hpp:14
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.
Definition: transform.hpp:88
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...
Definition: transform.hpp:1285
T positive_free(const T y)
Return the unconstrained value corresponding to the specified positive-constrained value...
Definition: transform.hpp:491
T positive_constrain(const T x)
Return the positive value for the specified unconstrained input.
Definition: transform.hpp:446
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...
Definition: transform.hpp:280
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)
Definition: check_less.hpp:52
const double CONSTRAINT_TOLERANCE
Definition: transform.hpp:32
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > cholesky_corr_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &y, int K)
Definition: transform.hpp:1512
fvar< T > log(const fvar< T > &x)
Definition: log.hpp:15
fvar< T > exp(const fvar< T > &x)
Definition: exp.hpp:16
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)
Definition: log1p.hpp:16
boost::math::tools::promote_args< T >::type log1m(T x)
Return the natural logarithm of one minus the specified value.
Definition: log1m.hpp:40
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...
Definition: transform.hpp:570

     [ Stan Home Page ] © 2011–2014, Stan Development Team.