Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
reader.hpp
Go to the documentation of this file.
1 #ifndef STAN__IO__READER_HPP
2 #define STAN__IO__READER_HPP
3 
4 #include <boost/throw_exception.hpp>
6 
7 namespace stan {
8 
9  namespace io {
10 
11 
33  template <typename T>
34  class reader {
35 
36  private:
37 
38  std::vector<T>& data_r_;
39  std::vector<int>& data_i_;
40  size_t pos_;
41  size_t int_pos_;
42 
43  inline T& scalar_ptr() {
44  return data_r_.at(pos_);
45  }
46 
47  inline T& scalar_ptr_increment(size_t m) {
48  pos_ += m;
49  return data_r_.at(pos_ - m);
50  }
51 
52  inline int& int_ptr() {
53  return data_i_.at(int_pos_);
54  }
55 
56  inline int& int_ptr_increment(size_t m) {
57  int_pos_ += m;
58  return data_i_.at(int_pos_ - m);
59  }
60 
61  public:
62 
63  typedef Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> matrix_t;
64  typedef Eigen::Matrix<T,Eigen::Dynamic,1> vector_t;
65  typedef Eigen::Matrix<T,1,Eigen::Dynamic> row_vector_t;
66 
67  typedef Eigen::Map<matrix_t> map_matrix_t;
68  typedef Eigen::Map<vector_t> map_vector_t;
69  typedef Eigen::Map<row_vector_t> map_row_vector_t;
70 
71 
83  reader(std::vector<T>& data_r,
84  std::vector<int>& data_i)
85  : data_r_(data_r),
86  data_i_(data_i),
87  pos_(0),
88  int_pos_(0) {
89  }
90 
94  ~reader() { }
95 
101  inline size_t available() {
102  return data_r_.size() - pos_;
103  }
104 
110  inline size_t available_i() {
111  return data_i_.size() - int_pos_;
112  }
113 
119  inline int integer() {
120  if (int_pos_ >= data_i_.size())
121  BOOST_THROW_EXCEPTION(
122  std::runtime_error("no more integers to read."));
123  return data_i_[int_pos_++];
124  }
125 
133  inline int integer_constrain() {
134  return integer();
135  }
136 
144  inline int integer_constrain(T& /*log_prob*/) {
145  return integer();
146  }
147 
148 
149 
155  inline T scalar() {
156  if (pos_ >= data_r_.size())
157  BOOST_THROW_EXCEPTION(std::runtime_error("no more scalars to read"));
158  return data_r_[pos_++];
159  }
160 
167  inline T scalar_constrain() {
168  return scalar();
169  }
170 
182  T scalar_constrain(T& /*log_prob*/) {
183  return scalar();
184  }
185 
186 
194  inline std::vector<T> std_vector(size_t m) {
195  std::vector<T> vec;
196  T& start = scalar_ptr_increment(m);
197  vec.insert(vec.begin(), &start, &scalar_ptr());
198  return vec;
199  }
200 
208  inline vector_t vector(size_t m) {
209  return map_vector_t(&scalar_ptr_increment(m),m);
210  }
218  inline vector_t vector_constrain(size_t m) {
219  return map_vector_t(&scalar_ptr_increment(m),m);
220  }
229  inline vector_t vector_constrain(size_t m, T& /*lp*/) {
230  return map_vector_t(&scalar_ptr_increment(m),m);
231  }
232 
233 
234 
242  inline row_vector_t row_vector(size_t m) {
243  return map_row_vector_t(&scalar_ptr_increment(m),m);
244  }
245 
254  return map_row_vector_t(&scalar_ptr_increment(m),m);
255  }
256 
266  inline row_vector_t row_vector_constrain(size_t m, T& /*lp*/) {
267  return map_row_vector_t(&scalar_ptr_increment(m),m);
268  }
269 
287  inline matrix_t matrix(size_t m, size_t n) {
288  return map_matrix_t(&scalar_ptr_increment(m*n),m,n);
289  }
290 
301  inline matrix_t matrix_constrain(size_t m, size_t n) {
302  return map_matrix_t(&scalar_ptr_increment(m*n),m,n);
303  }
304 
317  inline matrix_t matrix_constrain(size_t m, size_t n, T& /*lp*/) {
318  return map_matrix_t(&scalar_ptr_increment(m*n),m,n);
319  }
320 
321 
331  inline int integer_lb(int lb) {
332  int i = integer();
333  if (!(i >= lb))
334  BOOST_THROW_EXCEPTION(
335  std::runtime_error("required value greater than or equal to lb"));
336  return i;
337  }
347  inline int integer_lb_constrain(int lb) {
348  return integer_lb(lb);
349  }
360  inline int integer_lb_constrain(int lb, T& /*lp*/) {
361  return integer_lb(lb);
362  }
363 
364 
374  inline int integer_ub(int ub) {
375  int i = integer();
376  if (!(i <= ub))
377  BOOST_THROW_EXCEPTION(
378  std::runtime_error("required value less than or equal to ub"));
379  return i;
380  }
390  inline int integer_ub_constrain(int ub) {
391  return integer_ub(ub);
392  }
403  int integer_ub_constrain(int ub, T& /*lp*/) {
404  return integer_ub(ub);
405  }
406 
419  inline int integer_lub(int lb, int ub) {
420  // read first to make position deterministic [arbitrary choice]
421  int i = integer();
422  if (lb > ub)
423  BOOST_THROW_EXCEPTION(
424  std::runtime_error("lower bound must be less than or equal to ub"));
425  if (!(i >= lb))
426  BOOST_THROW_EXCEPTION(
427  std::runtime_error("required value greater than or equal to lb"));
428  if (!(i <= ub))
429  BOOST_THROW_EXCEPTION(
430  std::runtime_error("required value less than or equal to ub"));
431  return i;
432  }
443  inline int integer_lub_constrain(int lb, int ub) {
444  return integer_lub(lb,ub);
445  }
457  inline int integer_lub_constrain(int lb, int ub, T& /*lp*/) {
458  return integer_lub(lb,ub);
459  }
460 
461 
462 
472  inline T scalar_pos() {
473  T x(scalar());
474  stan::math::check_positive("stan::io::scalar_pos(%1%)", x,
475  "Constrained scalar", (double*)0);
476  return x;
477  }
478 
486  inline T scalar_pos_constrain() {
488  }
489 
500  inline T scalar_pos_constrain(T& lp) {
502  }
503 
516  template <typename TL>
517  inline T scalar_lb(const TL lb) {
518  T x(scalar());
519  stan::math::check_greater_or_equal("stan::io::scalar_lb(%1%)",
520  x, lb, "Constrained scalar",
521  (double*)0);
522  return x;
523  }
524 
536  template <typename TL>
537  inline T scalar_lb_constrain(const TL lb) {
538  return stan::prob::lb_constrain(scalar(),lb);
539  }
540 
552  template <typename TL>
553  inline T scalar_lb_constrain(const TL lb, T& lp) {
554  return stan::prob::lb_constrain(scalar(),lb,lp);
555  }
556 
557 
558 
571  template <typename TU>
572  inline T scalar_ub(TU ub) {
573  T x(scalar());
574  stan::math::check_less_or_equal("stan::io::scalar_ub(%1%)", x, ub,
575  "Constrained scalar", (double*)0);
576  return x;
577  }
578 
590  template <typename TU>
591  inline T scalar_ub_constrain(const TU ub) {
592  return stan::prob::ub_constrain(scalar(),ub);
593  }
594 
606  template <typename TU>
607  inline T scalar_ub_constrain(const TU ub, T& lp) {
608  return stan::prob::ub_constrain(scalar(),ub,lp);
609  }
610 
625  template <typename TL, typename TU>
626  inline T scalar_lub(const TL lb, const TU ub) {
627  T x(scalar());
628  stan::math::check_bounded<T,TL,TU,T>
629  ("stan::io::scalar_lub(%1%)", x, lb, ub, "Constrained scalar",0);
630  return x;
631  }
632 
646  template <typename TL, typename TU>
647  inline T scalar_lub_constrain(const TL lb, const TU ub) {
648  return stan::prob::lub_constrain(scalar(),lb,ub);
649  }
650 
664  template <typename TL, typename TU>
665  inline T scalar_lub_constrain(TL lb, TU ub, T& lp) {
666  return stan::prob::lub_constrain(scalar(),lb,ub,lp);
667  }
668 
677  inline T prob() {
678  T x(scalar());
679  stan::math::check_bounded<T,double,double,double>
680  ("stan::io::prob(%1%)", x, 0, 1, "Constrained probability", 0);
681  return x;
682  }
683 
692  inline T prob_constrain() {
694  }
695 
706  inline T prob_constrain(T& lp) {
707  return stan::prob::prob_constrain(scalar(),lp);
708  }
709 
710 
711 
712 
724  inline T corr() {
725  T x(scalar());
726  stan::math::check_bounded<T,double,double,double>
727  ("stan::io::corr(%1%)", x, -1, 1, "Correlation value",0);
728  return x;
729  }
730 
739  inline T corr_constrain() {
741  }
742 
754  inline T corr_constrain(T& lp) {
755  return stan::prob::corr_constrain(scalar(),lp);
756  }
757 
768  inline vector_t unit_vector(size_t k) {
769  vector_t theta(vector(k));
770  stan::math::check_unit_vector("stan::io::unit_vector(%1%)", theta, "Constrained vector",
771  (double*)0);
772  return theta;
773  }
774 
785  inline
786  Eigen::Matrix<T,Eigen::Dynamic,1> unit_vector_constrain(size_t k) {
788  }
789 
802  inline vector_t unit_vector_constrain(size_t k, T& lp) {
804  }
805 
816  inline vector_t simplex(size_t k) {
817  vector_t theta(vector(k));
818  stan::math::check_simplex("stan::io::simplex(%1%)", theta, "Constrained vector",
819  (double*)0);
820  return theta;
821  }
822 
833  inline
834  Eigen::Matrix<T,Eigen::Dynamic,1> simplex_constrain(size_t k) {
836  }
837 
850  inline vector_t simplex_constrain(size_t k, T& lp) {
851  return stan::prob::simplex_constrain(vector(k-1),lp);
852  }
853 
864  inline vector_t ordered(size_t k) {
865  vector_t x(vector(k));
866  stan::math::check_ordered("stan::io::ordered(%1%)", x, "Constrained vector");
867  return x;
868  }
869 
879  inline vector_t ordered_constrain(size_t k) {
881  }
882 
894  inline vector_t ordered_constrain(size_t k, T& lp) {
896  }
897 
908  inline vector_t positive_ordered(size_t k) {
909  vector_t x(vector(k));
910  stan::math::check_positive_ordered("stan::io::positive_ordered(%1%)",
911  x, "Constrained vector",
912  (double*)0);
913  return x;
914  }
915 
927  }
928 
940  inline vector_t positive_ordered_constrain(size_t k, T& lp) {
942  }
943 
944 
945 
956  inline matrix_t cholesky_factor(size_t M, size_t N) {
957  matrix_t y(matrix(M,N));
958  stan::math::check_cholesky_factor("stan::io::cholesky_factor(%1%)",
959  y, "Constrained matrix",
960  (double*)0);
961  return y;
962  }
963 
975  inline matrix_t cholesky_factor_constrain(size_t M, size_t N) {
976  return stan::prob::cholesky_factor_constrain(vector((N * (N + 1)) / 2 + (M - N) * N),
977  M,N);
978  }
979 
993  inline matrix_t cholesky_factor_constrain(size_t M, size_t N, T& lp) {
994  return stan::prob::cholesky_factor_constrain(vector((N * (N + 1)) / 2 + (M - N) * N),
995  M,N,lp);
996  }
997 
998 
999 
1010  inline matrix_t cholesky_corr(size_t K) {
1011  matrix_t y(matrix(K,K));
1012  stan::math::check_cholesky_factor_corr("stan::io::cholesky_factor_corr(%1%)", y, "Constrained matrix");
1013  return y;
1014  }
1015 
1027  return stan::prob::cholesky_corr_constrain(vector((K * (K - 1)) / 2),
1028  K);
1029  }
1030 
1043  inline matrix_t cholesky_corr_constrain(size_t K, T& lp) {
1044  return stan::prob::cholesky_corr_constrain(vector((K * (K - 1)) / 2),
1045  K,lp);
1046  }
1047 
1048 
1049 
1061  inline matrix_t cov_matrix(size_t k) {
1062  matrix_t y(matrix(k,k));
1063  stan::math::check_cov_matrix("stan::io::cov_matrix(%1%)", y, "Constrained matrix",
1064  (double*)0);
1065  return y;
1066  }
1067 
1076  inline matrix_t cov_matrix_constrain(size_t k) {
1077  return stan::prob::cov_matrix_constrain(vector(k + (k * (k - 1)) / 2),
1078  k);
1079  }
1080 
1092  inline matrix_t cov_matrix_constrain(size_t k, T& lp) {
1093  return stan::prob::cov_matrix_constrain(vector(k + (k * (k - 1)) / 2),
1094  k,lp);
1095  }
1096 
1097 
1107  inline matrix_t corr_matrix(size_t k) {
1108  matrix_t x(matrix(k,k));
1109  stan::math::check_corr_matrix("stan::math::corr_matrix(%1%)",
1110  x, "Constrained matrix",
1111  (double*)0);
1112  return x;
1113  }
1114 
1123  inline matrix_t corr_matrix_constrain(size_t k) {
1124  return stan::prob::corr_matrix_constrain(vector((k * (k - 1)) / 2),k);
1125  }
1126 
1138  inline matrix_t corr_matrix_constrain(size_t k, T& lp) {
1139  return stan::prob::corr_matrix_constrain(vector((k * (k - 1)) / 2),
1140  k,lp);
1141  }
1142 
1143  template <typename TL>
1144  inline vector_t vector_lb(const TL lb, size_t m) {
1145  vector_t v(m);
1146  for (size_t i = 0; i < m; ++i)
1147  v(i) = scalar_lb(lb);
1148  return v;
1149  }
1150  template <typename TL>
1151  inline vector_t vector_lb_constrain(const TL lb, size_t m) {
1152  vector_t v(m);
1153  for (size_t i = 0; i < m; ++i)
1154  v(i) = scalar_lb_constrain(lb);
1155  return v;
1156  }
1157  template <typename TL>
1158  inline vector_t vector_lb_constrain(const TL lb, size_t m, T& lp) {
1159  vector_t v(m);
1160  for (size_t i = 0; i < m; ++i)
1161  v(i) = scalar_lb_constrain(lb,lp);
1162  return v;
1163  }
1164 
1165  template <typename TL>
1166  inline row_vector_t row_vector_lb(const TL lb, size_t m) {
1167  row_vector_t v(m);
1168  for (size_t i = 0; i < m; ++i)
1169  v(i) = scalar_lb(lb);
1170  return v;
1171  }
1172  template <typename TL>
1173  inline row_vector_t row_vector_lb_constrain(const TL lb, size_t m) {
1174  row_vector_t v(m);
1175  for (size_t i = 0; i < m; ++i)
1176  v(i) = scalar_lb_constrain(lb);
1177  return v;
1178  }
1179  template <typename TL>
1180  inline row_vector_t row_vector_lb_constrain(const TL lb, size_t m, T& lp) {
1181  row_vector_t v(m);
1182  for (size_t i = 0; i < m; ++i)
1183  v(i) = scalar_lb_constrain(lb,lp);
1184  return v;
1185  }
1186 
1187  template <typename TL>
1188  inline matrix_t matrix_lb(const TL lb, size_t m, size_t n) {
1189  matrix_t v(m,n);
1190  for (size_t j = 0; j < n; ++j)
1191  for (size_t i = 0; i < m; ++i)
1192  v(i,j) = scalar_lb(lb);
1193  return v;
1194  }
1195  template <typename TL>
1196  inline matrix_t matrix_lb_constrain(const TL lb, size_t m, size_t n) {
1197  matrix_t v(m,n);
1198  for (size_t j = 0; j < n; ++j)
1199  for (size_t i = 0; i < m; ++i)
1200  v(i,j) = scalar_lb_constrain(lb);
1201  return v;
1202  }
1203  template <typename TL>
1204  inline matrix_t matrix_lb_constrain(const TL lb, size_t m, size_t n, T& lp) {
1205  matrix_t v(m,n);
1206  for (size_t j = 0; j < n; ++j)
1207  for (size_t i = 0; i < m; ++i)
1208  v(i,j) = scalar_lb_constrain(lb,lp);
1209  return v;
1210  }
1211 
1212  template <typename TU>
1213  inline vector_t vector_ub(const TU ub, size_t m) {
1214  vector_t v(m);
1215  for (size_t i = 0; i < m; ++i)
1216  v(i) = scalar_ub(ub);
1217  return v;
1218  }
1219  template <typename TU>
1220  inline vector_t vector_ub_constrain(const TU ub, size_t m) {
1221  vector_t v(m);
1222  for (size_t i = 0; i < m; ++i)
1223  v(i) = scalar_ub_constrain(ub);
1224  return v;
1225  }
1226  template <typename TU>
1227  inline vector_t vector_ub_constrain(const TU ub, size_t m, T& lp) {
1228  vector_t v(m);
1229  for (size_t i = 0; i < m; ++i)
1230  v(i) = scalar_ub_constrain(ub,lp);
1231  return v;
1232  }
1233 
1234  template <typename TU>
1235  inline row_vector_t row_vector_ub(const TU ub, size_t m) {
1236  row_vector_t v(m);
1237  for (size_t i = 0; i < m; ++i)
1238  v(i) = scalar_ub(ub);
1239  return v;
1240  }
1241  template <typename TU>
1242  inline row_vector_t row_vector_ub_constrain(const TU ub, size_t m) {
1243  row_vector_t v(m);
1244  for (size_t i = 0; i < m; ++i)
1245  v(i) = scalar_ub_constrain(ub);
1246  return v;
1247  }
1248  template <typename TU>
1249  inline row_vector_t row_vector_ub_constrain(const TU ub, size_t m, T& lp) {
1250  row_vector_t v(m);
1251  for (size_t i = 0; i < m; ++i)
1252  v(i) = scalar_ub_constrain(ub,lp);
1253  return v;
1254  }
1255 
1256  template <typename TU>
1257  inline matrix_t matrix_ub(const TU ub, size_t m, size_t n) {
1258  matrix_t v(m,n);
1259  for (size_t j = 0; j < n; ++j)
1260  for (size_t i = 0; i < m; ++i)
1261  v(i,j) = scalar_ub(ub);
1262  return v;
1263  }
1264  template <typename TU>
1265  inline matrix_t matrix_ub_constrain(const TU ub, size_t m, size_t n) {
1266  matrix_t v(m,n);
1267  for (size_t j = 0; j < n; ++j)
1268  for (size_t i = 0; i < m; ++i)
1269  v(i,j) = scalar_ub_constrain(ub);
1270  return v;
1271  }
1272  template <typename TU>
1273  inline matrix_t matrix_ub_constrain(const TU ub, size_t m, size_t n, T& lp) {
1274  matrix_t v(m,n);
1275  for (size_t j = 0; j < n; ++j)
1276  for (size_t i = 0; i < m; ++i)
1277  v(i,j) = scalar_ub_constrain(ub,lp);
1278  return v;
1279  }
1280 
1281 
1282  template <typename TL, typename TU>
1283  inline vector_t vector_lub(const TL lb, const TU ub, size_t m) {
1284  vector_t v(m);
1285  for (size_t i = 0; i < m; ++i)
1286  v(i) = scalar_lub(lb,ub);
1287  return v;
1288  }
1289  template <typename TL, typename TU>
1290  inline vector_t vector_lub_constrain(const TL lb, const TU ub, size_t m) {
1291  vector_t v(m);
1292  for (size_t i = 0; i < m; ++i)
1293  v(i) = scalar_lub_constrain(lb,ub);
1294  return v;
1295  }
1296  template <typename TL, typename TU>
1297  inline vector_t vector_lub_constrain(const TL lb, const TU ub, size_t m, T& lp) {
1298  vector_t v(m);
1299  for (size_t i = 0; i < m; ++i)
1300  v(i) = scalar_lub_constrain(lb,ub,lp);
1301  return v;
1302  }
1303 
1304  template <typename TL, typename TU>
1305  inline row_vector_t row_vector_lub(const TL lb, const TU ub, size_t m) {
1306  row_vector_t v(m);
1307  for (size_t i = 0; i < m; ++i)
1308  v(i) = scalar_lub(lb,ub);
1309  return v;
1310  }
1311  template <typename TL, typename TU>
1312  inline row_vector_t row_vector_lub_constrain(const TL lb, const TU ub, size_t m) {
1313  row_vector_t v(m);
1314  for (size_t i = 0; i < m; ++i)
1315  v(i) = scalar_lub_constrain(lb,ub);
1316  return v;
1317  }
1318  template <typename TL, typename TU>
1319  inline row_vector_t row_vector_lub_constrain(const TL lb, const TU ub, size_t m, T& lp) {
1320  row_vector_t v(m);
1321  for (size_t i = 0; i < m; ++i)
1322  v(i) = scalar_lub_constrain(lb,ub,lp);
1323  return v;
1324  }
1325 
1326  template <typename TL, typename TU>
1327  inline matrix_t matrix_lub(const TL lb, const TU ub, size_t m, size_t n) {
1328  matrix_t v(m,n);
1329  for (size_t j = 0; j < n; ++j)
1330  for (size_t i = 0; i < m; ++i)
1331  v(i,j) = scalar_lub(lb,ub);
1332  return v;
1333  }
1334  template <typename TL, typename TU>
1335  inline matrix_t matrix_lub_constrain(const TL lb, const TU ub, size_t m, size_t n) {
1336  matrix_t v(m,n);
1337  for (size_t j = 0; j < n; ++j)
1338  for (size_t i = 0; i < m; ++i)
1339  v(i,j) = scalar_lub_constrain(lb,ub);
1340  return v;
1341  }
1342  template <typename TL, typename TU>
1343  inline matrix_t matrix_lub_constrain(const TL lb, const TU ub, size_t m, size_t n, T& lp) {
1344  matrix_t v(m,n);
1345  for (size_t j = 0; j < n; ++j)
1346  for (size_t i = 0; i < m; ++i)
1347  v(i,j) = scalar_lub_constrain(lb,ub,lp);
1348  return v;
1349  }
1350 
1351 
1352 
1353  };
1354 
1355  }
1356 
1357 }
1358 
1359 #endif
bool check_cov_matrix(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 covariance matrix.
matrix_t cholesky_corr_constrain(size_t K, T &lp)
Return the next Cholesky factor for a correlation matrix with the specified dimensionality, reading from an unconstrained vector of the appropriate size, and increment the log probability reference with the log Jacobian adjustment for the transform.
Definition: reader.hpp:1043
vector_t positive_ordered_constrain(size_t k)
Return the next positive ordered vector of the specified length.
Definition: reader.hpp:925
row_vector_t row_vector_ub_constrain(const TU ub, size_t m)
Definition: reader.hpp:1242
vector_t simplex_constrain(size_t k, T &lp)
Return the next simplex of the specified size (using one fewer unconstrained scalars), incrementing the specified reference with the log absolute Jacobian determinant.
Definition: reader.hpp:850
T scalar_lub_constrain(const TL lb, const TU ub)
Return the next scalar transformed to be between the specified lower and upper bounds.
Definition: reader.hpp:647
matrix_t matrix_lub(const TL lb, const TU ub, size_t m, size_t n)
Definition: reader.hpp:1327
T scalar_ub(TU ub)
Return the next scalar, checking that it is less than or equal to the specified upper bound...
Definition: reader.hpp:572
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.
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
T scalar_pos()
Return the next scalar, checking that it is positive.
Definition: reader.hpp:472
vector_t unit_vector_constrain(size_t k, T &lp)
Return the next unit_vector of the specified size (using one fewer unconstrained scalars), incrementing the specified reference with the log absolute Jacobian determinant.
Definition: reader.hpp:802
int integer_ub_constrain(int ub)
Return the next integer, checking that it is less than or equal to the specified upper bound...
Definition: reader.hpp:390
matrix_t matrix_ub(const TU ub, size_t m, size_t n)
Definition: reader.hpp:1257
matrix_t matrix_lub_constrain(const TL lb, const TU ub, size_t m, size_t n)
Definition: reader.hpp:1335
std::vector< T > std_vector(size_t m)
Return a standard library vector of the specified dimensionality made up of the next scalars...
Definition: reader.hpp:194
T scalar_lb(const TL lb)
Return the next scalar, checking that it is greater than or equal to the specified lower bound...
Definition: reader.hpp:517
row_vector_t row_vector_constrain(size_t m)
Return a row vector of specified dimensionality made up of the next scalars.
Definition: reader.hpp:253
int integer_ub_constrain(int ub, T &)
Return the next integer, checking that it is less than or equal to the specified upper bound...
Definition: reader.hpp:403
T scalar_constrain()
Return the next scalar.
Definition: reader.hpp:167
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
vector_t vector_lub_constrain(const TL lb, const TU ub, size_t m, T &lp)
Definition: reader.hpp:1297
Eigen::Map< row_vector_t > map_row_vector_t
Definition: reader.hpp:69
int integer_lb_constrain(int lb, T &)
Return the next integer, checking that it is greater than or equal to the specified lower bound...
Definition: reader.hpp:360
vector_t vector_constrain(size_t m)
Return a column vector of specified dimensionality made up of the next scalars.
Definition: reader.hpp:218
row_vector_t row_vector_ub(const TU ub, size_t m)
Definition: reader.hpp:1235
matrix_t cov_matrix(size_t k)
Return the next covariance matrix with the specified dimensionality.
Definition: reader.hpp:1061
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
matrix_t matrix_lb_constrain(const TL lb, size_t m, size_t n, T &lp)
Definition: reader.hpp:1204
T corr_constrain()
Return the next scalar transformed to be a correlation between -1 and 1.
Definition: reader.hpp:739
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.
T scalar_constrain(T &)
Return the next scalar in the sequence, incrementing the specified reference with the log absolute Ja...
Definition: reader.hpp:182
matrix_t cov_matrix_constrain(size_t k)
Return the next covariance matrix of the specified dimensionality.
Definition: reader.hpp:1076
T prob_constrain(T &lp)
Return the next scalar transformed to be a probability between 0 and 1, incrementing the specified re...
Definition: reader.hpp:706
matrix_t corr_matrix_constrain(size_t k, T &lp)
Return the next correlation matrix of the specified dimensionality, incrementing the specified refere...
Definition: reader.hpp:1138
int integer_constrain(T &)
Return the next integer in the integer sequence.
Definition: reader.hpp:144
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
T corr()
Return the next scalar, checking that it is a valid value for a correlation, between -1 (inclusive) a...
Definition: reader.hpp:724
int integer_ub(int ub)
Return the next integer, checking that it is less than or equal to the specified upper bound...
Definition: reader.hpp:374
int integer_lb_constrain(int lb)
Return the next integer, checking that it is greater than or equal to the specified lower bound...
Definition: reader.hpp:347
T scalar_lb_constrain(const TL lb)
Return the next scalar transformed to have the specified lower bound.
Definition: reader.hpp:537
matrix_t cholesky_corr(size_t K)
Return the next Cholesky factor for a correlation matrix with the specified dimensionality, reading it directly without transforms.
Definition: reader.hpp:1010
T scalar_pos_constrain()
Return the next scalar, transformed to be positive.
Definition: reader.hpp:486
Eigen::Matrix< T, Eigen::Dynamic, 1 > vector_t
Definition: reader.hpp:64
vector_t vector_ub(const TU ub, size_t m)
Definition: reader.hpp:1213
row_vector_t row_vector(size_t m)
Return a row vector of specified dimensionality made up of the next scalars.
Definition: reader.hpp:242
vector_t vector(size_t m)
Return a column vector of specified dimensionality made up of the next scalars.
Definition: reader.hpp:208
vector_t vector_lub(const TL lb, const TU ub, size_t m)
Definition: reader.hpp:1283
row_vector_t row_vector_constrain(size_t m, T &)
Return a row vector of specified dimensionality made up of the next scalars.
Definition: reader.hpp:266
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
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.
matrix_t matrix_constrain(size_t m, size_t n)
Return a matrix of the specified dimensionality made up of the next scalars arranged in column-major ...
Definition: reader.hpp:301
vector_t unit_vector(size_t k)
Return a unit_vector of the specified size made up of the next scalars.
Definition: reader.hpp:768
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
row_vector_t row_vector_lb_constrain(const TL lb, size_t m)
Definition: reader.hpp:1173
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > matrix_t
Definition: reader.hpp:63
vector_t vector_ub_constrain(const TU ub, size_t m, T &lp)
Definition: reader.hpp:1227
Eigen::Map< vector_t > map_vector_t
Definition: reader.hpp:68
matrix_t cholesky_factor(size_t M, size_t N)
Return the next Cholesky factor with the specified dimensionality, reading it directly without transf...
Definition: reader.hpp:956
vector_t ordered(size_t k)
Return the next vector of specified size containing values in ascending order.
Definition: reader.hpp:864
A stream-based reader for integer, scalar, vector, matrix and array data types, with Jacobian calcula...
Definition: reader.hpp:34
reader(std::vector< T > &data_r, std::vector< int > &data_i)
Construct a variable reader using the specified vectors as the source of scalar and integer values fo...
Definition: reader.hpp:83
T scalar_ub_constrain(const TU ub)
Return the next scalar transformed to have the specified upper bound.
Definition: reader.hpp:591
matrix_t matrix_lub_constrain(const TL lb, const TU ub, size_t m, size_t n, T &lp)
Definition: reader.hpp:1343
int integer_lub(int lb, int ub)
Return the next integer, checking that it is less than or equal to the specified upper bound...
Definition: reader.hpp:419
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
matrix_t matrix(size_t m, size_t n)
Return a matrix of the specified dimensionality made up of the next scalars arranged in column-major ...
Definition: reader.hpp:287
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.
T corr_constrain(T &lp)
Return the next scalar transformed to be a (partial) correlation between -1 and 1, incrementing the specified reference with the log of the absolute Jacobian determinant.
Definition: reader.hpp:754
row_vector_t row_vector_lub_constrain(const TL lb, const TU ub, size_t m, T &lp)
Definition: reader.hpp:1319
T scalar()
Return the next scalar in the sequence.
Definition: reader.hpp:155
T scalar_pos_constrain(T &lp)
Return the next scalar transformed to be positive, incrementing the specified reference with the log ...
Definition: reader.hpp:500
vector_t ordered_constrain(size_t k)
Return the next ordered vector of the specified length.
Definition: reader.hpp:879
matrix_t matrix_lb(const TL lb, size_t m, size_t n)
Definition: reader.hpp:1188
int integer_lub_constrain(int lb, int ub)
Return the next integer, checking that it is less than or equal to the specified upper bound...
Definition: reader.hpp:443
matrix_t cholesky_factor_constrain(size_t M, size_t N, T &lp)
Return the next Cholesky factor with the specified dimensionality, reading from an unconstrained vect...
Definition: reader.hpp:993
T prob_constrain()
Return the next scalar transformed to be a probability between 0 and 1.
Definition: reader.hpp:692
vector_t simplex(size_t k)
Return a simplex of the specified size made up of the next scalars.
Definition: reader.hpp:816
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
matrix_t corr_matrix(size_t k)
Returns the next correlation matrix of the specified dimensionality.
Definition: reader.hpp:1107
vector_t vector_lb(const TL lb, size_t m)
Definition: reader.hpp:1144
T scalar_lb_constrain(const TL lb, T &lp)
Return the next scalar transformed to have the specified lower bound, incrementing the specified refe...
Definition: reader.hpp:553
bool check_positive(const char *function, const T_y &y, const char *name, T_result *result)
bool check_corr_matrix(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 correlation matrix.
vector_t vector_lb_constrain(const TL lb, size_t m)
Definition: reader.hpp:1151
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
vector_t positive_ordered(size_t k)
Return the next vector of specified size containing positive values in ascending order.
Definition: reader.hpp:908
size_t available()
Return the number of scalars remaining to be read.
Definition: reader.hpp:101
matrix_t cov_matrix_constrain(size_t k, T &lp)
Return the next covariance matrix of the specified dimensionality, incrementing the specified referen...
Definition: reader.hpp:1092
row_vector_t row_vector_lub_constrain(const TL lb, const TU ub, size_t m)
Definition: reader.hpp:1312
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
~reader()
Destroy this variable reader.
Definition: reader.hpp:94
row_vector_t row_vector_lub(const TL lb, const TU ub, size_t m)
Definition: reader.hpp:1305
int integer()
Return the next integer in the integer sequence.
Definition: reader.hpp:119
Eigen::Matrix< T, 1, Eigen::Dynamic > row_vector_t
Definition: reader.hpp:65
bool check_less_or_equal(const char *function, const T_y &y, const T_high &high, const char *name, T_result *result)
matrix_t matrix_ub_constrain(const TU ub, size_t m, size_t n, T &lp)
Definition: reader.hpp:1273
matrix_t matrix_ub_constrain(const TU ub, size_t m, size_t n)
Definition: reader.hpp:1265
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
matrix_t cholesky_factor_constrain(size_t M, size_t N)
Return the next Cholesky factor with the specified dimensionality, reading from an unconstrained vect...
Definition: reader.hpp:975
int integer_constrain()
Return the next integer in the integer sequence.
Definition: reader.hpp:133
size_t available_i()
Return the number of integers remaining to be read.
Definition: reader.hpp:110
row_vector_t row_vector_lb(const TL lb, size_t m)
Definition: reader.hpp:1166
vector_t vector_lub_constrain(const TL lb, const TU ub, size_t m)
Definition: reader.hpp:1290
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...
Eigen::Map< matrix_t > map_matrix_t
Definition: reader.hpp:67
int integer_lb(int lb)
Return the next integer, checking that it is greater than or equal to the specified lower bound...
Definition: reader.hpp:331
int integer_lub_constrain(int lb, int ub, T &)
Return the next integer, checking that it is less than or equal to the specified upper bound...
Definition: reader.hpp:457
T scalar_ub_constrain(const TU ub, T &lp)
Return the next scalar transformed to have the specified upper bound, incrementing the specified refe...
Definition: reader.hpp:607
T scalar_lub_constrain(TL lb, TU ub, T &lp)
Return the next scalar transformed to be between the the specified lower and upper bounds...
Definition: reader.hpp:665
matrix_t corr_matrix_constrain(size_t k)
Return the next correlation matrix of the specified dimensionality.
Definition: reader.hpp:1123
Eigen::Matrix< T, Eigen::Dynamic, 1 > unit_vector_constrain(size_t k)
Return the next unit_vector transformed vector of the specified length.
Definition: reader.hpp:786
T positive_constrain(const T x)
Return the positive value for the specified unconstrained input.
Definition: transform.hpp:446
matrix_t matrix_constrain(size_t m, size_t n, T &)
Return a matrix of the specified dimensionality made up of the next scalars arranged in column-major ...
Definition: reader.hpp:317
vector_t vector_lb_constrain(const TL lb, size_t m, T &lp)
Definition: reader.hpp:1158
vector_t ordered_constrain(size_t k, T &lp)
Return the next ordered vector of the specified size, incrementing the specified reference with the l...
Definition: reader.hpp:894
matrix_t cholesky_corr_constrain(size_t K)
Return the next Cholesky factor for a correlation matrix with the specified dimensionality, reading from an unconstrained vector of the appropriate size.
Definition: reader.hpp:1026
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > cholesky_corr_constrain(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &y, int K)
Definition: transform.hpp:1512
T prob()
Return the next scalar, checking that it is a valid value for a probability, between 0 (inclusive) an...
Definition: reader.hpp:677
vector_t positive_ordered_constrain(size_t k, T &lp)
Return the next positive_ordered vector of the specified size, incrementing the specified reference w...
Definition: reader.hpp:940
row_vector_t row_vector_ub_constrain(const TU ub, size_t m, T &lp)
Definition: reader.hpp:1249
T scalar_lub(const TL lb, const TU ub)
Return the next scalar, checking that it is between the specified lower and upper bound...
Definition: reader.hpp:626
row_vector_t row_vector_lb_constrain(const TL lb, size_t m, T &lp)
Definition: reader.hpp:1180
matrix_t matrix_lb_constrain(const TL lb, size_t m, size_t n)
Definition: reader.hpp:1196
bool check_greater_or_equal(const char *function, const T_y &y, const T_low &low, const char *name, T_result *result)
vector_t vector_ub_constrain(const TU ub, size_t m)
Definition: reader.hpp:1220
Eigen::Matrix< T, Eigen::Dynamic, 1 > simplex_constrain(size_t k)
Return the next simplex transformed vector of the specified length.
Definition: reader.hpp:834
bool check_cholesky_factor_corr(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.
vector_t vector_constrain(size_t m, T &)
Return a column vector of specified dimensionality made up of the next scalars.
Definition: reader.hpp:229

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