1 #ifndef STAN__MCMC__CHAINS_HPP
2 #define STAN__MCMC__CHAINS_HPP
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>
25 #include <boost/random/uniform_int_distribution.hpp>
26 #include <boost/random/additive_combine.hpp>
54 template <
typename RNG = boost::random::ecuyer1988>
57 Eigen::Matrix<std::string, Eigen::Dynamic, 1> param_names_;
58 Eigen::Matrix<Eigen::MatrixXd, Eigen::Dynamic, 1> samples_;
59 Eigen::VectorXi warmup_;
61 static double mean(
const Eigen::VectorXd& x) {
62 return (x.array() / x.size()).
sum();
65 static double variance(
const Eigen::VectorXd& x) {
67 return ((x.array() - m) /
std::sqrt((x.size() - 1.0))).square().sum();
70 static double sd(
const Eigen::VectorXd& x) {
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;
85 accumulator_set<double, stats<covariance<double, covariate1> > > acc;
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));
91 return boost::accumulators::covariance(acc) * M / (M-1);
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;
104 accumulator_set<double, stats<variance,
105 covariance<double, covariate1> > > acc_xy;
106 accumulator_set<double, stats<variance> > acc_y;
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));
114 double cov = boost::accumulators::covariance(acc_xy);
115 if (cov > -1
e-8 && cov < 1
e-8)
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;
131 size_t cache_size = M;
134 accumulator_set<double, stats<tail_quantile<left> > >
135 acc(tail<left>::cache_size = cache_size);
136 for (
int i = 0; i < M; i++)
138 return boost::accumulators::quantile(acc, quantile_probability=prob);
140 accumulator_set<double, stats<tail_quantile<right> > >
141 acc(tail<right>::cache_size = cache_size);
142 for (
int i = 0; i < M; i++)
144 return boost::accumulators::quantile(acc, quantile_probability=prob);
147 static Eigen::VectorXd
148 quantiles(
const Eigen::VectorXd& x,
const Eigen::VectorXd& probs)
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;
160 size_t cache_size = M;
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);
167 for (
int i = 0; i < M; i++) {
172 Eigen::VectorXd q(probs.size());
173 for (
int i = 0; i < probs.size(); i++) {
175 q(i) = boost::accumulators::quantile(acc_left,
176 quantile_probability=probs(i));
178 q(i) = boost::accumulators::quantile(acc_right,
179 quantile_probability=probs(i));
184 static Eigen::VectorXd autocorrelation(
const Eigen::VectorXd& x) {
187 typedef typename index_type<vector<double> >::type idx_t;
189 std::vector<double> ac;
190 std::vector<double>
sample(x.size());
191 for (
int i = 0; i < x.size(); i++)
195 Eigen::VectorXd ac2(ac.size());
196 for (idx_t i = 0; i < ac.size(); i++)
201 static Eigen::VectorXd autocovariance(
const Eigen::VectorXd& x) {
204 typedef typename index_type<vector<double> >::type idx_t;
206 std::vector<double> ac;
207 std::vector<double>
sample(x.size());
208 for (
int i = 0; i < x.size(); i++)
212 Eigen::VectorXd ac2(ac.size());
213 for (idx_t i = 0; i < ac.size(); i++)
233 double effective_sample_size(
const Eigen::Matrix<Eigen::VectorXd,
234 Eigen::Dynamic, 1> &
samples)
const {
238 int n_samples =
samples(0).size();
239 for (
int chain = 1; chain <
chains; chain++) {
243 Eigen::Matrix<Eigen::VectorXd, Eigen::Dynamic, 1> acov(chains);
244 for (
int chain = 0; chain <
chains; chain++) {
245 acov(chain) = autocovariance(
samples(chain));
248 Eigen::VectorXd chain_mean(chains);
249 Eigen::VectorXd chain_var(chains);
250 for (
int chain = 0; chain <
chains; chain++) {
252 chain_mean(chain) = mean(
samples(chain));
253 chain_var(chain) = acov(chain)(0)*n_kept_samples/(n_kept_samples-1);
256 double mean_var = mean(chain_var);
257 double var_plus = mean_var*(n_samples-1)/n_samples;
259 var_plus += variance(chain_mean);
260 Eigen::VectorXd rho_hat_t(n_samples);
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);
269 rho_hat = 1 - (mean_var - mean(acov_t)) / var_plus;
271 rho_hat_t(t) = rho_hat;
274 double ess = chains * n_samples;
276 ess /= 1 + 2 * rho_hat_t.sum();
295 split_potential_scale_reduction(
296 const Eigen::Matrix<Eigen::VectorXd,
297 Eigen::Dynamic, 1> &
samples)
const {
299 int n_samples =
samples(0).size();
300 for (
int chain = 1; chain <
chains; chain++) {
303 if (n_samples % 2 == 1)
305 int n = n_samples / 2;
307 Eigen::VectorXd split_chain_mean(2*chains);
308 Eigen::VectorXd split_chain_var(2*chains);
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));
314 split_chain_var(2*chain) = variance(
samples(chain).topRows(n));
315 split_chain_var(2*chain+1) = variance(
samples(chain).bottomRows(n));
318 double var_between = n * variance(split_chain_mean);
319 double var_within = mean(split_chain_var);
322 return sqrt((var_between/var_within + n-1)/n);
327 : param_names_(param_names) { }
330 for (
size_t i = 0; i < param_names.size(); i++)
331 param_names_(i) = param_names[i];
335 : param_names_(stan_csv.header) {
336 if (stan_csv.
samples.rows() > 0)
341 return samples_.size();
345 return param_names_.size();
348 const Eigen::Matrix<std::string, Eigen::Dynamic, 1>&
param_names()
const {
353 return param_names_(j);
356 int index(
const std::string& name)
const {
358 for (
int i = 0; i < param_names_.size(); i++)
359 if (param_names_(i) == name)
369 warmup_.setConstant(warmup);
377 return warmup_(chain);
381 return samples_(chain).rows();
386 for (
int chain = 0; chain <
num_chains(); chain++)
397 for (
int chain = 0; chain <
num_chains(); chain++)
403 const Eigen::MatrixXd&
sample) {
405 throw std::invalid_argument(
"add(chain,sample): number of columns"
406 " in sample does not match chains");
412 Eigen::Matrix<Eigen::MatrixXd, Eigen::Dynamic, 1>
415 for (
int i = 0; i < n; i++) {
416 samples_copy(i) = samples_(i);
417 warmup_copy(i) = warmup_(i);
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);
426 for (
int i = n; i < chain+1; i++) {
427 samples_(i) = Eigen::MatrixXd(0,
num_params());
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;
438 if (sample.rows() == 0)
441 throw std::invalid_argument(
"add(sample): number of columns in"
442 " sample does not match chains");
453 void add(
const std::vector<std::vector<double> >&
sample) {
454 int n_row =
sample.size();
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++) {
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"
486 for (
int chain = 0; chain <
num_chains(); chain++) {
488 s.middleRows(start,n) = samples_(chain).col(index).bottomRows(n);
494 Eigen::VectorXd
samples(
const int chain,
const std::string& name)
const {
498 Eigen::VectorXd
samples(
const std::string& name)
const {
503 return mean(
samples(chain,index));
510 double mean(
const int chain,
const std::string& name)
const {
511 return mean(chain,
index(name));
514 double mean(
const std::string& name)
const {
515 return mean(
index(name));
518 double sd(
const int chain,
const int index)
const {
519 return sd(
samples(chain,index));
526 double sd(
const int chain,
const std::string& name)
const {
527 return sd(chain,
index(name));
530 double sd(
const std::string& name)
const {
531 return sd(
index(name));
535 return variance(
samples(chain,index));
539 return variance(
samples(index));
542 double variance(
const int chain,
const std::string& name)
const {
543 return variance(chain,
index(name));
547 return variance(
index(name));
551 covariance(
const int chain,
const int index1,
const int index2)
const {
555 double covariance(
const int index1,
const int index2)
const {
560 const std::string& name2)
const {
561 return covariance(chain,
index(name1),
index(name2));
565 covariance(
const std::string& name1,
const std::string& name2)
const {
566 return covariance(
index(name1),
index(name2));
570 correlation(
const int chain,
const int index1,
const int index2)
const {
579 const std::string& name2)
const {
580 return correlation(chain,
index(name1),
index(name2));
584 correlation(
const std::string& name1,
const std::string& name2)
const {
585 return correlation(
index(name1),
index(name2));
590 return quantile(
samples(chain,index), prob);
594 return quantile(
samples(index), prob);
597 double quantile(
int chain,
const std::string& name,
double prob)
const {
598 return quantile(chain,
index(name), prob);
601 double quantile(
const std::string& name,
const double prob)
const {
602 return quantile(
index(name), prob);
607 return quantiles(
samples(chain,index), probs);
611 return quantiles(
samples(index), probs);
616 const Eigen::VectorXd& probs)
const {
617 return quantiles(chain,
index(name), probs);
621 quantiles(
const std::string& name,
const Eigen::VectorXd& probs)
const {
622 return quantiles(
index(name), probs);
627 double low_prob = (1-prob)/2;
628 double high_prob = 1-low_prob;
630 Eigen::Vector2d interval;
632 << quantile(chain,index,low_prob), quantile(chain,index,high_prob);
637 double low_prob = (1-prob)/2;
638 double high_prob = 1-low_prob;
640 Eigen::Vector2d interval;
641 interval << quantile(index,low_prob), quantile(index,high_prob);
657 return autocorrelation(
samples(chain, index));
661 const std::string& name)
const {
662 return autocorrelation(chain,
index(name));
666 return autocovariance(
samples(chain,index));
670 return autocovariance(chain,
index(name));
675 Eigen::Matrix<Eigen::VectorXd, Eigen::Dynamic, 1>
677 for (
int chain = 0; chain <
num_chains(); chain++) {
680 return effective_sample_size(samples);
684 return effective_sample_size(
index(name));
689 Eigen::Matrix<Eigen::VectorXd, Eigen::Dynamic, 1>
691 for (
int chain = 0; chain <
num_chains(); chain++) {
694 return split_potential_scale_reduction(samples);
698 return split_potential_scale_reduction(
index(name));
int warmup(const int chain) const
int num_kept_samples() const
double sd(const std::string &name) const
void set_warmup(const int warmup)
Eigen::VectorXd samples(const int chain, const std::string &name) const
double covariance(const std::string &name1, const std::string &name2) const
const std::string & param_name(int j) const
double correlation(const int index1, const int index2) const
double quantile(const std::string &name, const double prob) const
double mean(const int chain, const int index) const
Eigen::VectorXd quantiles(int chain, const std::string &name, const Eigen::VectorXd &probs) const
double quantile(const int index, const double prob) const
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.
chains(const stan::io::stan_csv &stan_csv)
void add(const stan::io::stan_csv &stan_csv)
Eigen::VectorXd autocovariance(int chain, const std::string &name) const
Eigen::VectorXd samples(const int chain, const int index) const
double sd(const int chain, const int index) const
int index(const std::string &name) const
Eigen::Vector2d central_interval(int index, double prob) const
double correlation(const std::string &name1, const std::string &name2) const
double mean(const std::string &name) const
double variance(const std::string &name) const
double sd(const int index) const
fvar< T > sum(const Eigen::Matrix< fvar< T >, R, C > &m)
double split_potential_scale_reduction(const int index) const
double variance(const int chain, const int index) const
double quantile(int chain, const std::string &name, double prob) const
double covariance(const int chain, const int index1, const int index2) const
double effective_sample_size(const int index) const
Primary template class for the metaprogram to compute the index type of a container.
double correlation(const int chain, const int index1, const int index2) const
double mean(const int chain, const std::string &name) const
double sd(const int chain, const std::string &name) const
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
void add(const std::vector< std::vector< double > > &sample)
Convert a vector of vector<double> to Eigen::MatrixXd.
const Eigen::Matrix< std::string, Eigen::Dynamic, 1 > & param_names() const
Eigen::VectorXd autocorrelation(const int chain, const int index) const
fvar< T > sqrt(const fvar< T > &x)
double variance(const int chain, const std::string &name) const
chains(const std::vector< std::string > ¶m_names)
Eigen::VectorXd samples(const int index) const
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...
int num_samples(const int chain) const
double correlation(const int chain, const std::string &name1, const std::string &name2) const
void add(const Eigen::MatrixXd &sample)
double effective_sample_size(const std::string &name) const
double variance(const int index) const
double mean(const int index) const
An mcmc::chains object stores parameter names and dimensionalities along with samples from multiple c...
Eigen::Vector2d central_interval(int chain, const std::string &name, double prob) const
Eigen::VectorXd autocovariance(const int chain, const int index) const
double split_potential_scale_reduction(const std::string &name) const
double e()
Return the base of the natural logarithm.
void add(const int chain, const Eigen::MatrixXd &sample)
Eigen::Vector2d central_interval(const std::string &name, double prob) const
Eigen::VectorXd quantiles(int index, const Eigen::VectorXd &probs) const
int size(const std::vector< T > &x)
Eigen::VectorXd autocorrelation(int chain, const std::string &name) const
Eigen::VectorXd quantiles(int chain, int index, const Eigen::VectorXd &probs) const
double min(const double a, const double b)
chains(const Eigen::Matrix< std::string, Eigen::Dynamic, 1 > ¶m_names)
var variance(const std::vector< var > &v)
Return the sample variance of the specified standard vector.
Eigen::VectorXd quantiles(const std::string &name, const Eigen::VectorXd &probs) const
stan_csv_metadata metadata
double covariance(const int chain, const std::string &name1, const std::string &name2) const
double quantile(const int chain, const int index, const double prob) const
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)
Eigen::VectorXd samples(const std::string &name) const
Eigen::Vector2d central_interval(int chain, int index, double prob) const
const Eigen::VectorXi & warmup() const
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
void set_warmup(const int chain, const int warmup)