Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
chains.hpp
Go to the documentation of this file.
1 #ifndef STAN__MCMC__CHAINS_HPP
2 #define STAN__MCMC__CHAINS_HPP
3 
4 #include <algorithm>
5 #include <cmath>
6 #include <iostream>
7 #include <map>
8 #include <stdexcept>
9 #include <string>
10 #include <sstream>
11 #include <utility>
12 #include <vector>
13 #include <fstream>
14 #include <cstdlib>
15 
16 #include <boost/accumulators/accumulators.hpp>
17 #include <boost/accumulators/statistics/stats.hpp>
18 #include <boost/accumulators/statistics/mean.hpp>
19 #include <boost/accumulators/statistics/tail_quantile.hpp>
20 #include <boost/accumulators/statistics/p_square_quantile.hpp>
21 #include <boost/accumulators/statistics/variance.hpp>
22 #include <boost/accumulators/statistics/covariance.hpp>
23 #include <boost/accumulators/statistics/variates/covariate.hpp>
24 
25 #include <boost/random/uniform_int_distribution.hpp>
26 #include <boost/random/additive_combine.hpp>
27 
29 #include <stan/math/matrix.hpp>
34 
35 
36 namespace stan {
37 
38  namespace mcmc {
39 
54  template <typename RNG = boost::random::ecuyer1988>
55  class chains {
56  private:
57  Eigen::Matrix<std::string, Eigen::Dynamic, 1> param_names_;
58  Eigen::Matrix<Eigen::MatrixXd, Eigen::Dynamic, 1> samples_;
59  Eigen::VectorXi warmup_;
60 
61  static double mean(const Eigen::VectorXd& x) {
62  return (x.array() / x.size()).sum();
63  }
64 
65  static double variance(const Eigen::VectorXd& x) {
66  double m = mean(x);
67  return ((x.array() - m) / std::sqrt((x.size() - 1.0))).square().sum();
68  }
69 
70  static double sd(const Eigen::VectorXd& x) {
71  return std::sqrt(variance(x));
72  }
73 
74 
75  static double covariance(const Eigen::VectorXd& x,
76  const Eigen::VectorXd& y) {
77  if (x.rows() != y.rows())
78  std::cerr << "warning: covariance of different length chains";
79  using boost::accumulators::accumulator_set;
80  using boost::accumulators::stats;
82  using boost::accumulators::tag::covariance;
83  using boost::accumulators::tag::covariate1;
84 
85  accumulator_set<double, stats<covariance<double, covariate1> > > acc;
86 
87  int M = std::min(x.size(), y.size());
88  for (int i = 0; i < M; i++)
89  acc(x(i), boost::accumulators::covariate1=y(i));
90 
91  return boost::accumulators::covariance(acc) * M / (M-1);
92  }
93 
94  static double correlation(const Eigen::VectorXd& x,
95  const Eigen::VectorXd& y) {
96  if (x.rows() != y.rows())
97  std::cerr << "warning: covariance of different length chains";
98  using boost::accumulators::accumulator_set;
99  using boost::accumulators::stats;
101  using boost::accumulators::tag::covariance;
102  using boost::accumulators::tag::covariate1;
103 
104  accumulator_set<double, stats<variance,
105  covariance<double, covariate1> > > acc_xy;
106  accumulator_set<double, stats<variance> > acc_y;
107 
108  int M = std::min(x.size(), y.size());
109  for (int i = 0; i < M; i++) {
110  acc_xy(x(i), boost::accumulators::covariate1=y(i));
111  acc_y(y(i));
112  }
113 
114  double cov = boost::accumulators::covariance(acc_xy);
115  if (cov > -1e-8 && cov < 1e-8)
116  return cov;
117  return cov / std::sqrt(boost::accumulators::variance(acc_xy)
119  }
120 
121  static double quantile(const Eigen::VectorXd& x, const double prob) {
122  using boost::accumulators::accumulator_set;
123  using boost::accumulators::left;
124  using boost::accumulators::quantile_probability;
125  using boost::accumulators::right;
126  using boost::accumulators::stats;
128  using boost::accumulators::tag::tail_quantile;
129  double M = x.rows();
130  //size_t cache_size = std::min(prob, 1-prob)*M + 2;
131  size_t cache_size = M;
132 
133  if (prob < 0.5) {
134  accumulator_set<double, stats<tail_quantile<left> > >
135  acc(tail<left>::cache_size = cache_size);
136  for (int i = 0; i < M; i++)
137  acc(x(i));
138  return boost::accumulators::quantile(acc, quantile_probability=prob);
139  }
140  accumulator_set<double, stats<tail_quantile<right> > >
141  acc(tail<right>::cache_size = cache_size);
142  for (int i = 0; i < M; i++)
143  acc(x(i));
144  return boost::accumulators::quantile(acc, quantile_probability=prob);
145  }
146 
147  static Eigen::VectorXd
148  quantiles(const Eigen::VectorXd& x, const Eigen::VectorXd& probs)
149  {
150  using boost::accumulators::accumulator_set;
151  using boost::accumulators::left;
152  using boost::accumulators::quantile_probability;
153  using boost::accumulators::right;
154  using boost::accumulators::stats;
156  using boost::accumulators::tag::tail_quantile;
157  double M = x.rows();
158 
159  //size_t cache_size = M/2 + 2;
160  size_t cache_size = M;
161 
162  accumulator_set<double, stats<tail_quantile<left> > >
163  acc_left(tail<left>::cache_size = cache_size);
164  accumulator_set<double, stats<tail_quantile<right> > >
165  acc_right(tail<right>::cache_size = cache_size);
166 
167  for (int i = 0; i < M; i++) {
168  acc_left(x(i));
169  acc_right(x(i));
170  }
171 
172  Eigen::VectorXd q(probs.size());
173  for (int i = 0; i < probs.size(); i++) {
174  if (probs(i) < 0.5)
175  q(i) = boost::accumulators::quantile(acc_left,
176  quantile_probability=probs(i));
177  else
178  q(i) = boost::accumulators::quantile(acc_right,
179  quantile_probability=probs(i));
180  }
181  return q;
182  }
183 
184  static Eigen::VectorXd autocorrelation(const Eigen::VectorXd& x) {
185  using std::vector;
187  typedef typename index_type<vector<double> >::type idx_t;
188 
189  std::vector<double> ac;
190  std::vector<double> sample(x.size());
191  for (int i = 0; i < x.size(); i++)
192  sample[i] = x(i);
194 
195  Eigen::VectorXd ac2(ac.size());
196  for (idx_t i = 0; i < ac.size(); i++)
197  ac2(i) = ac[i];
198  return ac2;
199  }
200 
201  static Eigen::VectorXd autocovariance(const Eigen::VectorXd& x) {
202  using std::vector;
204  typedef typename index_type<vector<double> >::type idx_t;
205 
206  std::vector<double> ac;
207  std::vector<double> sample(x.size());
208  for (int i = 0; i < x.size(); i++)
209  sample[i] = x(i);
211 
212  Eigen::VectorXd ac2(ac.size());
213  for (idx_t i = 0; i < ac.size(); i++)
214  ac2(i) = ac[i];
215  return ac2;
216  }
217 
233  double effective_sample_size(const Eigen::Matrix<Eigen::VectorXd,
234  Eigen::Dynamic, 1> &samples) const {
235  int chains = samples.size();
236 
237  // need to generalize to each jagged samples per chain
238  int n_samples = samples(0).size();
239  for (int chain = 1; chain < chains; chain++) {
240  n_samples = std::min(n_samples, int(samples(chain).size()));
241  }
242 
243  Eigen::Matrix<Eigen::VectorXd, Eigen::Dynamic, 1> acov(chains);
244  for (int chain = 0; chain < chains; chain++) {
245  acov(chain) = autocovariance(samples(chain));
246  }
247 
248  Eigen::VectorXd chain_mean(chains);
249  Eigen::VectorXd chain_var(chains);
250  for (int chain = 0; chain < chains; chain++) {
251  double n_kept_samples = num_kept_samples(chain);
252  chain_mean(chain) = mean(samples(chain));
253  chain_var(chain) = acov(chain)(0)*n_kept_samples/(n_kept_samples-1);
254  }
255 
256  double mean_var = mean(chain_var);
257  double var_plus = mean_var*(n_samples-1)/n_samples;
258  if (chains > 1)
259  var_plus += variance(chain_mean);
260  Eigen::VectorXd rho_hat_t(n_samples);
261  rho_hat_t.setZero();
262  double rho_hat = 0;
263  int max_t = 0;
264  for (int t = 1; (t < n_samples && rho_hat >= 0); t++) {
265  Eigen::VectorXd acov_t(chains);
266  for (int chain = 0; chain < chains; chain++) {
267  acov_t(chain) = acov(chain)(t);
268  }
269  rho_hat = 1 - (mean_var - mean(acov_t)) / var_plus;
270  if (rho_hat >= 0)
271  rho_hat_t(t) = rho_hat;
272  max_t = t;
273  }
274  double ess = chains * n_samples;
275  if (max_t > 1) {
276  ess /= 1 + 2 * rho_hat_t.sum();
277  }
278  return ess;
279  }
280 
294  double
295  split_potential_scale_reduction(
296  const Eigen::Matrix<Eigen::VectorXd,
297  Eigen::Dynamic, 1> &samples) const {
298  int chains = samples.size();
299  int n_samples = samples(0).size();
300  for (int chain = 1; chain < chains; chain++) {
301  n_samples = std::min(n_samples, int(samples(chain).size()));
302  }
303  if (n_samples % 2 == 1)
304  n_samples--;
305  int n = n_samples / 2;
306 
307  Eigen::VectorXd split_chain_mean(2*chains);
308  Eigen::VectorXd split_chain_var(2*chains);
309 
310  for (int chain = 0; chain < chains; chain++) {
311  split_chain_mean(2*chain) = mean(samples(chain).topRows(n));
312  split_chain_mean(2*chain+1) = mean(samples(chain).bottomRows(n));
313 
314  split_chain_var(2*chain) = variance(samples(chain).topRows(n));
315  split_chain_var(2*chain+1) = variance(samples(chain).bottomRows(n));
316  }
317 
318  double var_between = n * variance(split_chain_mean);
319  double var_within = mean(split_chain_var);
320 
321  // rewrote [(n-1)*W/n + B/n]/W as (n-1+ B/W)/n
322  return sqrt((var_between/var_within + n-1)/n);
323  }
324 
325  public:
326  chains(const Eigen::Matrix<std::string, Eigen::Dynamic, 1>& param_names)
327  : param_names_(param_names) { }
328 
329  chains(const std::vector<std::string>& param_names) : param_names_(param_names.size()) {
330  for (size_t i = 0; i < param_names.size(); i++)
331  param_names_(i) = param_names[i];
332  }
333 
334  chains(const stan::io::stan_csv& stan_csv)
335  : param_names_(stan_csv.header) {
336  if (stan_csv.samples.rows() > 0)
337  add(stan_csv);
338  }
339 
340  inline int num_chains() const {
341  return samples_.size();
342  }
343 
344  inline int num_params() const {
345  return param_names_.size();
346  }
347 
348  const Eigen::Matrix<std::string, Eigen::Dynamic, 1>& param_names() const {
349  return param_names_;
350  }
351 
352  const std::string& param_name(int j) const {
353  return param_names_(j);
354  }
355 
356  int index(const std::string& name) const {
357  int index = -1;
358  for (int i = 0; i < param_names_.size(); i++)
359  if (param_names_(i) == name)
360  return i;
361  return index;
362  }
363 
364  void set_warmup(const int chain, const int warmup) {
365  warmup_(chain) = warmup;
366  }
367 
368  void set_warmup(const int warmup) {
369  warmup_.setConstant(warmup);
370  }
371 
372  const Eigen::VectorXi& warmup() const {
373  return warmup_;
374  }
375 
376  int warmup(const int chain) const {
377  return warmup_(chain);
378  }
379 
380  int num_samples(const int chain) const {
381  return samples_(chain).rows();
382  }
383 
384  int num_samples() const {
385  int n = 0;
386  for (int chain = 0; chain < num_chains(); chain++)
387  n += num_samples(chain);
388  return n;
389  }
390 
391  int num_kept_samples(const int chain) const {
392  return num_samples(chain) - warmup(chain);
393  }
394 
395  int num_kept_samples() const {
396  int n = 0;
397  for (int chain = 0; chain < num_chains(); chain++)
398  n += num_kept_samples(chain);
399  return n;
400  }
401 
402  void add(const int chain,
403  const Eigen::MatrixXd& sample) {
404  if (sample.cols() != num_params())
405  throw std::invalid_argument("add(chain,sample): number of columns"
406  " in sample does not match chains");
407  if (num_chains() == 0 || chain >= num_chains()) {
408  int n = num_chains();
409 
410  // Need this block for Windows. conservativeResize
411  // does not keep the references.
412  Eigen::Matrix<Eigen::MatrixXd, Eigen::Dynamic, 1>
413  samples_copy(num_chains());
414  Eigen::VectorXi warmup_copy(num_chains());
415  for (int i = 0; i < n; i++) {
416  samples_copy(i) = samples_(i);
417  warmup_copy(i) = warmup_(i);
418  }
419 
420  samples_.resize(chain+1);
421  warmup_.resize(chain+1);
422  for (int i = 0; i < n; i++) {
423  samples_(i) = samples_copy(i);
424  warmup_(i) = warmup_copy(i);
425  }
426  for (int i = n; i < chain+1; i++) {
427  samples_(i) = Eigen::MatrixXd(0, num_params());
428  warmup_(i) = 0;
429  }
430  }
431  int row = samples_(chain).rows();
432  Eigen::MatrixXd new_samples(row+sample.rows(), num_params());
433  new_samples << samples_(chain), sample;
434  samples_(chain) = new_samples;
435  }
436 
437  void add(const Eigen::MatrixXd& sample) {
438  if (sample.rows() == 0)
439  return;
440  if (sample.cols() != num_params())
441  throw std::invalid_argument("add(sample): number of columns in"
442  " sample does not match chains");
443  add(num_chains(), sample);
444  }
445 
453  void add(const std::vector<std::vector<double> >& sample) {
454  int n_row = sample.size();
455  if (n_row == 0)
456  return;
457  int n_col = sample[0].size();
458  Eigen::MatrixXd sample_copy(n_row, n_col);
459  for (int i = 0; i < n_row; i++) {
460  sample_copy.row(i)
461  = Eigen::VectorXd::Map(&sample[i][0], sample[0].size());
462  }
463  add(sample_copy);
464  }
465 
466  void add(const stan::io::stan_csv& stan_csv) {
467  if (stan_csv.header.size() != num_params())
468  throw std::invalid_argument("add(stan_csv): number of columns in"
469  " sample does not match chains");
470  if (!param_names_.cwiseEqual(stan_csv.header).all()) {
471  throw std::invalid_argument("add(stan_csv): header does not match"
472  " chain's header");
473  }
474  add(stan_csv.samples);
475  if (stan_csv.metadata.save_warmup)
476  set_warmup(num_chains()-1, stan_csv.metadata.num_warmup);
477  }
478 
479  Eigen::VectorXd samples(const int chain, const int index) const {
480  return samples_(chain).col(index).bottomRows(num_kept_samples(chain));
481  }
482 
483  Eigen::VectorXd samples(const int index) const {
484  Eigen::VectorXd s(num_kept_samples());
485  int start = 0;
486  for (int chain = 0; chain < num_chains(); chain++) {
487  int n = num_kept_samples(chain);
488  s.middleRows(start,n) = samples_(chain).col(index).bottomRows(n);
489  start += n;
490  }
491  return s;
492  }
493 
494  Eigen::VectorXd samples(const int chain, const std::string& name) const {
495  return samples(chain,index(name));
496  }
497 
498  Eigen::VectorXd samples(const std::string& name) const {
499  return samples(index(name));
500  }
501 
502  double mean(const int chain, const int index) const {
503  return mean(samples(chain,index));
504  }
505 
506  double mean(const int index) const {
507  return mean(samples(index));
508  }
509 
510  double mean(const int chain, const std::string& name) const {
511  return mean(chain, index(name));
512  }
513 
514  double mean(const std::string& name) const {
515  return mean(index(name));
516  }
517 
518  double sd(const int chain, const int index) const {
519  return sd(samples(chain,index));
520  }
521 
522  double sd(const int index) const {
523  return sd(samples(index));
524  }
525 
526  double sd(const int chain, const std::string& name) const {
527  return sd(chain, index(name));
528  }
529 
530  double sd(const std::string& name) const {
531  return sd(index(name));
532  }
533 
534  double variance(const int chain, const int index) const {
535  return variance(samples(chain,index));
536  }
537 
538  double variance(const int index) const {
539  return variance(samples(index));
540  }
541 
542  double variance(const int chain, const std::string& name) const {
543  return variance(chain, index(name));
544  }
545 
546  double variance(const std::string& name) const {
547  return variance(index(name));
548  }
549 
550  double
551  covariance(const int chain, const int index1, const int index2) const {
552  return covariance(samples(chain,index1), samples(chain,index2));
553  }
554 
555  double covariance(const int index1, const int index2) const {
556  return covariance(samples(index1), samples(index2));
557  }
558 
559  double covariance(const int chain, const std::string& name1,
560  const std::string& name2) const {
561  return covariance(chain, index(name1), index(name2));
562  }
563 
564  double
565  covariance(const std::string& name1, const std::string& name2) const {
566  return covariance(index(name1), index(name2));
567  }
568 
569  double
570  correlation(const int chain, const int index1, const int index2) const {
571  return correlation(samples(chain,index1),samples(chain,index2));
572  }
573 
574  double correlation(const int index1, const int index2) const {
575  return correlation(samples(index1),samples(index2));
576  }
577 
578  double correlation(const int chain, const std::string& name1,
579  const std::string& name2) const {
580  return correlation(chain, index(name1), index(name2));
581  }
582 
583  double
584  correlation(const std::string& name1, const std::string& name2) const {
585  return correlation(index(name1), index(name2));
586  }
587 
588  double
589  quantile(const int chain, const int index, const double prob) const {
590  return quantile(samples(chain,index), prob);
591  }
592 
593  double quantile(const int index, const double prob) const {
594  return quantile(samples(index), prob);
595  }
596 
597  double quantile(int chain, const std::string& name, double prob) const {
598  return quantile(chain, index(name), prob);
599  }
600 
601  double quantile(const std::string& name, const double prob) const {
602  return quantile(index(name), prob);
603  }
604 
605  Eigen::VectorXd
606  quantiles(int chain, int index, const Eigen::VectorXd& probs) const {
607  return quantiles(samples(chain,index), probs);
608  }
609 
610  Eigen::VectorXd quantiles(int index, const Eigen::VectorXd& probs) const {
611  return quantiles(samples(index), probs);
612  }
613 
614  Eigen::VectorXd
615  quantiles(int chain, const std::string& name,
616  const Eigen::VectorXd& probs) const {
617  return quantiles(chain, index(name), probs);
618  }
619 
620  Eigen::VectorXd
621  quantiles(const std::string& name, const Eigen::VectorXd& probs) const {
622  return quantiles(index(name), probs);
623  }
624 
625  Eigen::Vector2d central_interval(int chain, int index,
626  double prob) const {
627  double low_prob = (1-prob)/2;
628  double high_prob = 1-low_prob;
629 
630  Eigen::Vector2d interval;
631  interval
632  << quantile(chain,index,low_prob), quantile(chain,index,high_prob);
633  return interval;
634  }
635 
636  Eigen::Vector2d central_interval(int index, double prob) const {
637  double low_prob = (1-prob)/2;
638  double high_prob = 1-low_prob;
639 
640  Eigen::Vector2d interval;
641  interval << quantile(index,low_prob), quantile(index,high_prob);
642  return interval;
643  }
644 
645  Eigen::Vector2d
646  central_interval(int chain, const std::string& name,
647  double prob) const {
648  return central_interval(chain, index(name), prob);
649  }
650 
651  Eigen::Vector2d central_interval(const std::string& name,
652  double prob) const {
653  return central_interval(index(name), prob);
654  }
655 
656  Eigen::VectorXd autocorrelation(const int chain, const int index) const {
657  return autocorrelation(samples(chain, index));
658  }
659 
660  Eigen::VectorXd autocorrelation(int chain,
661  const std::string& name) const {
662  return autocorrelation(chain, index(name));
663  }
664 
665  Eigen::VectorXd autocovariance(const int chain, const int index) const {
666  return autocovariance(samples(chain,index));
667  }
668 
669  Eigen::VectorXd autocovariance(int chain, const std::string& name) const {
670  return autocovariance(chain, index(name));
671  }
672 
673  // FIXME: reimplement using autocorrelation.
674  double effective_sample_size(const int index) const {
675  Eigen::Matrix<Eigen::VectorXd, Eigen::Dynamic, 1>
676  samples(num_chains());
677  for (int chain = 0; chain < num_chains(); chain++) {
678  samples(chain) = this->samples(chain, index);
679  }
680  return effective_sample_size(samples);
681  }
682 
683  double effective_sample_size(const std::string& name) const {
684  return effective_sample_size(index(name));
685  }
686 
687  double split_potential_scale_reduction(const int index) const
688  {
689  Eigen::Matrix<Eigen::VectorXd, Eigen::Dynamic, 1>
690  samples(num_chains());
691  for (int chain = 0; chain < num_chains(); chain++) {
692  samples(chain) = this->samples(chain, index);
693  }
694  return split_potential_scale_reduction(samples);
695  }
696 
697  double split_potential_scale_reduction(const std::string& name) const {
698  return split_potential_scale_reduction(index(name));
699  }
700  };
701 
702  }
703 }
704 
705 #endif
int warmup(const int chain) const
Definition: chains.hpp:376
int num_kept_samples() const
Definition: chains.hpp:395
double sd(const std::string &name) const
Definition: chains.hpp:530
void set_warmup(const int warmup)
Definition: chains.hpp:368
Eigen::VectorXd samples(const int chain, const std::string &name) const
Definition: chains.hpp:494
double covariance(const std::string &name1, const std::string &name2) const
Definition: chains.hpp:565
const std::string & param_name(int j) const
Definition: chains.hpp:352
double correlation(const int index1, const int index2) const
Definition: chains.hpp:574
Eigen::MatrixXd samples
double quantile(const std::string &name, const double prob) const
Definition: chains.hpp:601
double mean(const int chain, const int index) const
Definition: chains.hpp:502
Eigen::VectorXd quantiles(int chain, const std::string &name, const Eigen::VectorXd &probs) const
Definition: chains.hpp:615
double quantile(const int index, const double prob) const
Definition: chains.hpp:593
Eigen::Matrix< T, 1, Eigen::Dynamic > row(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &m, size_t i)
Return the specified row of the specified matrix, using start-at-1 indexing.
Definition: row.hpp:26
chains(const stan::io::stan_csv &stan_csv)
Definition: chains.hpp:334
void add(const stan::io::stan_csv &stan_csv)
Definition: chains.hpp:466
Eigen::VectorXd autocovariance(int chain, const std::string &name) const
Definition: chains.hpp:669
Eigen::VectorXd samples(const int chain, const int index) const
Definition: chains.hpp:479
double sd(const int chain, const int index) const
Definition: chains.hpp:518
int index(const std::string &name) const
Definition: chains.hpp:356
Eigen::Vector2d central_interval(int index, double prob) const
Definition: chains.hpp:636
double correlation(const std::string &name1, const std::string &name2) const
Definition: chains.hpp:584
double mean(const std::string &name) const
Definition: chains.hpp:514
double variance(const std::string &name) const
Definition: chains.hpp:546
double sd(const int index) const
Definition: chains.hpp:522
fvar< T > sum(const Eigen::Matrix< fvar< T >, R, C > &m)
Definition: sum.hpp:14
double split_potential_scale_reduction(const int index) const
Definition: chains.hpp:687
double variance(const int chain, const int index) const
Definition: chains.hpp:534
double quantile(int chain, const std::string &name, double prob) const
Definition: chains.hpp:597
double covariance(const int chain, const int index1, const int index2) const
Definition: chains.hpp:551
double effective_sample_size(const int index) const
Definition: chains.hpp:674
Primary template class for the metaprogram to compute the index type of a container.
Definition: index_type.hpp:21
double correlation(const int chain, const int index1, const int index2) const
Definition: chains.hpp:570
double mean(const int chain, const std::string &name) const
Definition: chains.hpp:510
double sd(const int chain, const std::string &name) const
Definition: chains.hpp:526
void autocovariance(const std::vector< T > &y, std::vector< T > &acov, Eigen::FFT< T > &fft)
Write autocovariance estimates for every lag for the specified input sequence into the specified resu...
double covariance(const int index1, const int index2) const
Definition: chains.hpp:555
void add(const std::vector< std::vector< double > > &sample)
Convert a vector of vector<double> to Eigen::MatrixXd.
Definition: chains.hpp:453
const Eigen::Matrix< std::string, Eigen::Dynamic, 1 > & param_names() const
Definition: chains.hpp:348
Eigen::VectorXd autocorrelation(const int chain, const int index) const
Definition: chains.hpp:656
fvar< T > sqrt(const fvar< T > &x)
Definition: sqrt.hpp:15
double variance(const int chain, const std::string &name) const
Definition: chains.hpp:542
chains(const std::vector< std::string > &param_names)
Definition: chains.hpp:329
Eigen::VectorXd samples(const int index) const
Definition: chains.hpp:483
Eigen::Matrix< T, Eigen::Dynamic, 1 > tail(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &v, size_t n)
Return the specified number of elements as a vector from the back of the specified vector...
Definition: tail.hpp:22
int num_samples(const int chain) const
Definition: chains.hpp:380
double correlation(const int chain, const std::string &name1, const std::string &name2) const
Definition: chains.hpp:578
void add(const Eigen::MatrixXd &sample)
Definition: chains.hpp:437
double effective_sample_size(const std::string &name) const
Definition: chains.hpp:683
double variance(const int index) const
Definition: chains.hpp:538
double mean(const int index) const
Definition: chains.hpp:506
An mcmc::chains object stores parameter names and dimensionalities along with samples from multiple c...
Definition: chains.hpp:55
int num_params() const
Definition: chains.hpp:344
Eigen::Vector2d central_interval(int chain, const std::string &name, double prob) const
Definition: chains.hpp:646
Eigen::VectorXd autocovariance(const int chain, const int index) const
Definition: chains.hpp:665
double split_potential_scale_reduction(const std::string &name) const
Definition: chains.hpp:697
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:86
void add(const int chain, const Eigen::MatrixXd &sample)
Definition: chains.hpp:402
Eigen::Vector2d central_interval(const std::string &name, double prob) const
Definition: chains.hpp:651
Eigen::VectorXd quantiles(int index, const Eigen::VectorXd &probs) const
Definition: chains.hpp:610
int size(const std::vector< T > &x)
Definition: size.hpp:11
Eigen::VectorXd autocorrelation(int chain, const std::string &name) const
Definition: chains.hpp:660
int num_chains() const
Definition: chains.hpp:340
int num_samples() const
Definition: chains.hpp:384
Eigen::VectorXd quantiles(int chain, int index, const Eigen::VectorXd &probs) const
Definition: chains.hpp:606
double min(const double a, const double b)
Definition: min.hpp:7
chains(const Eigen::Matrix< std::string, Eigen::Dynamic, 1 > &param_names)
Definition: chains.hpp:326
var variance(const std::vector< var > &v)
Return the sample variance of the specified standard vector.
Definition: variance.hpp:51
Eigen::VectorXd quantiles(const std::string &name, const Eigen::VectorXd &probs) const
Definition: chains.hpp:621
stan_csv_metadata metadata
double covariance(const int chain, const std::string &name1, const std::string &name2) const
Definition: chains.hpp:559
double quantile(const int chain, const int index, const double prob) const
Definition: chains.hpp:589
Eigen::Matrix< std::string, Eigen::Dynamic, 1 > header
void sample(stan::mcmc::base_mcmc *sampler, int num_warmup, int num_samples, int num_thin, int refresh, bool save, stan::io::mcmc_writer< Model, SampleRecorder, DiagnosticRecorder, MessageRecorder > &writer, stan::mcmc::sample &init_s, Model &model, RNG &base_rng, const std::string &prefix, const std::string &suffix, std::ostream &o, StartTransitionCallback &callback)
Definition: sample.hpp:13
Eigen::VectorXd samples(const std::string &name) const
Definition: chains.hpp:498
Eigen::Vector2d central_interval(int chain, int index, double prob) const
Definition: chains.hpp:625
const Eigen::VectorXi & warmup() const
Definition: chains.hpp:372
void autocorrelation(const std::vector< T > &y, std::vector< T > &ac, Eigen::FFT< T > &fft)
Write autocorrelation estimates for every lag for the specified input sequence into the specified res...
int num_kept_samples(const int chain) const
Definition: chains.hpp:391
void set_warmup(const int chain, const int warmup)
Definition: chains.hpp:364

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