MUQ  0.4.3
muq::SamplingAlgorithms::SampleCollection Class Reference

A class to hold and analyze a collection of SamplingState objects. More...

#include <SampleCollection.h>

Inheritance diagram for muq::SamplingAlgorithms::SampleCollection:

Detailed Description

A class to hold and analyze a collection of SamplingState objects.

Definition at line 71 of file SampleCollection.h.

Public Member Functions

 SampleCollection ()=default
 
virtual ~SampleCollection ()=default
 
virtual void Add (std::shared_ptr< SamplingState > newSamp)
 
virtual std::shared_ptr< SamplingStateat (unsigned i)
 
virtual const std::shared_ptr< SamplingStateat (unsigned i) const
 
virtual unsigned int size () const
 
virtual std::shared_ptr< SampleCollectionhead (unsigned int N) const
 
virtual std::shared_ptr< SampleCollectiontail (unsigned int N) const
 
virtual std::shared_ptr< SampleCollectionsegment (unsigned int startInd, unsigned int length, unsigned int skipBy=1) const
 
virtual const std::shared_ptr< SamplingStateback () const
 
virtual Eigen::VectorXd CentralMoment (unsigned order, Eigen::VectorXd const &mean, int blockNum=-1) const override
 Computes the componentwise central moments (e.g., variance, skewness, kurtosis, etc..) of a specific order given that we already know the mean. More...
 
virtual Eigen::VectorXd CentralMoment (unsigned int order, int blockNum=-1) const override
 
virtual Eigen::VectorXd Mean (int blockInd=-1) const override
 
virtual Eigen::MatrixXd Covariance (int blockInd=-1) const override
 
virtual Eigen::MatrixXd Covariance (Eigen::VectorXd const &mean, int blockInd=-1) const override
 
virtual std::vector< Eigen::MatrixXd > RunningCovariance (int blockInd=-1) const
 
virtual std::vector< Eigen::MatrixXd > RunningCovariance (Eigen::VectorXd const &mean, int blockInd=-1) const
 
virtual Eigen::VectorXd ESS (std::string const &method="Batch") const override
 Returns the effective sample size of the samples. More...
 
virtual Eigen::VectorXd ESS (int blockDim) const override
 
virtual Eigen::VectorXd ESS (int blockDim, std::string const &method) const override
 
Eigen::VectorXd BatchESS (int blockInd=-1, int batchSize=-1, int overlap=-1) const
 
double MultiBatchESS (int blockInd=-1, int batchSize=-1, int overlap=-1) const
 
virtual Eigen::VectorXd StandardError (int blockInd, std::string const &method) const override
 
virtual Eigen::VectorXd StandardError (int blockInd) const override
 
virtual Eigen::VectorXd StandardError (std::string const &method="Batch") const override
 
virtual Eigen::VectorXd BatchError (int blockInd=-1, int batchSize=-1, int overlap=-1) const
 
virtual Eigen::VectorXd MultiBatchError (int blockInd=-1, int batchSize=-1, int overlap=-1) const
 
virtual Eigen::MatrixXd AsMatrix (int blockDim=-1) const
 
virtual Eigen::VectorXd Weights () const
 
virtual Eigen::VectorXd Rhat (int blockDim, unsigned int numSegments=4, boost::property_tree::ptree options=boost::property_tree::ptree()) const
 
virtual Eigen::VectorXd Rhat (unsigned int numSegments=4, boost::property_tree::ptree options=boost::property_tree::ptree()) const
 
- Public Member Functions inherited from muq::SamplingAlgorithms::SampleEstimator
virtual ~SampleEstimator ()=default
 
virtual unsigned int BlockSize (int blockInd) const =0
 
virtual unsigned int NumBlocks () const =0
 
virtual Eigen::VectorXd CentralMoment (unsigned int order, Eigen::VectorXd const &mean, int blockNum=-1) const
 
virtual Eigen::VectorXd StandardizedMoment (unsigned int order, int blockInd=-1) const
 
virtual Eigen::VectorXd StandardizedMoment (unsigned int order, Eigen::VectorXd const &mean, int blockInd=-1) const
 
virtual Eigen::VectorXd StandardizedMoment (unsigned int order, Eigen::VectorXd const &mean, Eigen::VectorXd const &stdDev, int blockInd=-1) const
 
virtual Eigen::VectorXd Variance (int blockInd=-1) const
 
virtual Eigen::VectorXd Variance (Eigen::VectorXd const &mean, int blockInd=-1) const
 
virtual Eigen::VectorXd Skewness (int blockInd=-1) const
 
virtual Eigen::VectorXd Skewness (Eigen::VectorXd const &mean, int blockInd=-1) const
 
virtual Eigen::VectorXd Skewness (Eigen::VectorXd const &mean, Eigen::VectorXd const &stdDev, int blockInd=-1) const
 
virtual Eigen::VectorXd Kurtosis (int blockInd=-1) const
 
virtual Eigen::VectorXd Kurtosis (Eigen::VectorXd const &mean, int blockInd=-1) const
 
virtual Eigen::VectorXd Kurtosis (Eigen::VectorXd const &mean, Eigen::VectorXd const &stdDev, int blockInd=-1) const
 
virtual Eigen::VectorXd ExpectedValue (std::shared_ptr< muq::Modeling::ModPiece > const &f, std::vector< std::string > const &metains=std::vector< std::string >()) const =0
 

Static Public Member Functions

static std::shared_ptr< SampleCollectionFromMatrix (Eigen::MatrixXd const &samps)
 

Constructor & Destructor Documentation

◆ SampleCollection()

muq::SamplingAlgorithms::SampleCollection::SampleCollection ( )
default

◆ ~SampleCollection()

virtual muq::SamplingAlgorithms::SampleCollection::~SampleCollection ( )
virtualdefault

Member Function Documentation

◆ Add()

void SampleCollection::Add ( std::shared_ptr< SamplingState newSamp)
virtual

Reimplemented in muq::SamplingAlgorithms::DistributedCollection.

Definition at line 107 of file SampleCollection.cpp.

References samples.

◆ AsMatrix()

Eigen::MatrixXd SampleCollection::AsMatrix ( int  blockDim = -1) const
virtual

Returns the samples in this collection as a matrix. Each column of the matrix will correspond to a single state in the chain.

Reimplemented in muq::SamplingAlgorithms::DistributedCollection.

Definition at line 438 of file SampleCollection.cpp.

References samples.

Referenced by muq::SamplingAlgorithms::MarkovChain::WolffESS().

◆ at() [1/2]

std::shared_ptr< SamplingState > SampleCollection::at ( unsigned  i)
virtual

◆ at() [2/2]

const std::shared_ptr< SamplingState > SampleCollection::at ( unsigned  i) const
virtual

Reimplemented in muq::SamplingAlgorithms::DistributedCollection.

Definition at line 117 of file SampleCollection.cpp.

References samples.

◆ back()

const std::shared_ptr< SamplingState > SampleCollection::back ( ) const
virtual

Definition at line 134 of file SampleCollection.cpp.

References samples.

◆ BatchError()

Eigen::VectorXd SampleCollection::BatchError ( int  blockInd = -1,
int  batchSize = -1,
int  overlap = -1 
) const
virtual

Estimates the Monte Carlo standard error (square root of estimator variance) using overlapping batch means (OBM). For details, consult [Flegal2010] .

Parameters
[in]blockIndSpecifies the block of the sampling state we're interested in. Defaults to -1, which will result in all blocks of the sampling state being concatenated in the MCSE estimate.
[in]batchSizeThe size of the batches \(b_n\) to use. Defaults to \(\sqrt{N}\). Note that to ensure the convergence of the MCSE estimate as \(N\rightarrow\infty\), the batch size must increase with \(N\).
[in]overlapHow much the batches overlap. Defaults to \(0.75*b_n\), where \(b_n\) is the batch size. A nonoverlapping batch estimate will be use is overlap=0. This value bust be smaller than \(b_n\).
Returns
A vector containing either the MCSE for each component of the chain.

Definition at line 384 of file SampleCollection.cpp.

References muq::SamplingAlgorithms::SampleEstimator::BlockSize(), Mean(), segment(), and size().

Referenced by BatchESS(), muq::SamplingAlgorithms::MarkovChain::StandardError(), and StandardError().

◆ BatchESS()

Eigen::VectorXd SampleCollection::BatchESS ( int  blockInd = -1,
int  batchSize = -1,
int  overlap = -1 
) const

Computes the effective sample size using the overlapping or nonoverlapping batch method. (See e.g., [Flegal2010]).

Parameters
[in]blockIndThe block of the SampleState that we want to use in the ESS computation.
[in]batchSizeThe size of the batches \(b_n\) to use. Defaults to \(N^{1/3}\). Note that to ensure the convergence of the ESS estimate as \(N\rightarrow\infty\), the batch size must increase with \(N\).
[in]overlapHow much the batches overlap. Defaults to \(b_n/2\), where \(b_n\) is the batch size. A nonoverlapping batch estimate will be use is overlap=0. This value bust be smaller than \(b_n\).
Returns
A vector of length \(D\) containing an ESS estimate for each component

Definition at line 356 of file SampleCollection.cpp.

References BatchError(), size(), and muq::SamplingAlgorithms::SampleEstimator::Variance().

Referenced by muq::SamplingAlgorithms::MarkovChain::ESS(), and ESS().

◆ CentralMoment() [1/2]

virtual Eigen::VectorXd muq::SamplingAlgorithms::SampleCollection::CentralMoment ( unsigned int  order,
int  blockNum = -1 
) const
inlineoverridevirtual

The central moment of order $p$ is given by

\[ \mathbb{E}\left[ \left(X-\mathbb{E}[X]\right)^p \right] \]

This function returns a sample-based estimate of the central moment. In the default implementation of this function, the SampleEstimator::ExpectedValue function is called twice. Once to compute \(\mathbb{E}[X]\) and then again to compute \(\mathbb{E}\left[ \left(X-\mathbb{E}[X]\right)^p \right]\). Child classes might provide alternative implementations.

Parameters
[in]orderThe order \(p\) of the central moment. \(p=2\) yields the variance.
[in]blockNum(Optional) The block of the random variable \(x\) to use in the expectation. By default, blockNum=-1 and the expectation is computed with respect to the entire random variable \(x\).
Returns
A vector with the same size as \(x\) or \(x_i\) containing an estimate of the central moment.

Reimplemented from muq::SamplingAlgorithms::SampleEstimator.

Definition at line 108 of file SampleCollection.h.

References CentralMoment(), and Mean().

◆ CentralMoment() [2/2]

Eigen::VectorXd SampleCollection::CentralMoment ( unsigned  order,
Eigen::VectorXd const &  mean,
int  blockNum = -1 
) const
overridevirtual

Computes the componentwise central moments (e.g., variance, skewness, kurtosis, etc..) of a specific order given that we already know the mean.

Definition at line 139 of file SampleCollection.cpp.

References RecursiveSum(), and samples.

Referenced by CentralMoment().

◆ Covariance() [1/2]

Eigen::MatrixXd SampleCollection::Covariance ( Eigen::VectorXd const &  mean,
int  blockInd = -1 
) const
overridevirtual

Computes the sample covariance using a precomputed (or known) mean.

If blockInd is non-negative, only the mean of one block of the samples is computed.

Implements muq::SamplingAlgorithms::SampleEstimator.

Definition at line 161 of file SampleCollection.cpp.

References samples.

◆ Covariance() [2/2]

virtual Eigen::MatrixXd muq::SamplingAlgorithms::SampleCollection::Covariance ( int  blockInd = -1) const
inlineoverridevirtual

Computes the sample covariance of \(x\) with itself (if blockInd==-1) or \(x_i\) with itself (if blockInd==i), i.e.,

\[ \text{Cov}[x] = \mathbb{E}\left[ \left(x-\mathbb{E}[x]\right)\left(x-\mathbb{E}[x]\right)^T\right] \]

or

\[ \text{Cov}[x_i] = \mathbb{E}\left[ \left(x_i-\mathbb{E}[x_i]\right)\left(x_i-\mathbb{E}[x_i]\right)^T\right] \]

Note that it is only possible to compute the cross covariance of \(x_i\) with \(x_j\) by setting blockInd=-1 and computing the entire covariance matrix.

Parameters
[in]blockInd(Optional) The block of the random variable $x$ to use in the expectation. By default, blockInd=-1 and the expectation is computed with respect to the entire random variable $x$.
Returns
A matrix containing an estimate of \(\text{Cov}[x]\) or \(\text{Cov}[x_i]\).

Reimplemented from muq::SamplingAlgorithms::SampleEstimator.

Reimplemented in muq::SamplingAlgorithms::DistributedCollection.

Definition at line 112 of file SampleCollection.h.

References muq::SamplingAlgorithms::SampleEstimator::Covariance().

Referenced by MultiBatchESS().

◆ CreateDataset()

bool SampleCollection::CreateDataset ( std::shared_ptr< muq::Utilities::HDF5File hdf5file,
std::string const &  dataname,
int const  dataSize,
int const  totSamps 
) const
private
Parameters
[in]hdf5fileThe hdf5 file where the data will be written
[in]datanameThe name of the data set we (may) want to create
[in]dataSizeThe number of rows (the size of the data in one sample)
[in]totSampsThe total number of samples we need to write to the file (max. number of samples)
Returns
true: the data set exists and is the right size, false: the data set does not exist or is the wrong size

Definition at line 474 of file SampleCollection.cpp.

References size().

◆ ESS() [1/3]

virtual Eigen::VectorXd muq::SamplingAlgorithms::SampleCollection::ESS ( int  blockDim) const
inlineoverridevirtual
Parameters
[in]blockDimThe block of the SampleState that we want to use in the ESS computation.
Returns
Either a vector of length \(D\) containing an ESS estimate for each component or a vector of length 1 containing a single multivariate ESS estimate.

Reimplemented from muq::SamplingAlgorithms::SampleEstimator.

Reimplemented in muq::SamplingAlgorithms::DistributedCollection, and muq::SamplingAlgorithms::MarkovChain.

Definition at line 159 of file SampleCollection.h.

References ESS().

Referenced by ESS().

◆ ESS() [2/3]

Eigen::VectorXd SampleCollection::ESS ( int  blockDim,
std::string const &  method 
) const
overridevirtual
Parameters
[in]blockDimThe block of the SampleState that we want to use in the ESS computation.
[in]methodA string specifying the method used to estimate the ESS. Options in SampleCollection are "Batch" or "MultiBatch". The MarkovChain class also supports a "Wolff" method.
Returns
Either a vector of length \(D\) containing an ESS estimate for each component or a vector of length 1 containing a single multivariate ESS estimate.

Implements muq::SamplingAlgorithms::SampleEstimator.

Reimplemented in muq::SamplingAlgorithms::MarkovChain.

Definition at line 288 of file SampleCollection.cpp.

References BatchESS(), and MultiBatchESS().

◆ ESS() [3/3]

virtual Eigen::VectorXd muq::SamplingAlgorithms::SampleCollection::ESS ( std::string const &  method = "Batch") const
inlineoverridevirtual

Returns the effective sample size of the samples.

For almost all random variables of interest, the central limit theorem states that the variance of a Monte Carlo estimator

\[\hat{\mu}_N = \frac{1}{N}\sum_{i=1}^N \theta^{(i)},\]

satisfies

\[ \sqrt{N}\mathbb{V}\left[\hat{\mu}_N-\mu\right] \rightarrow N(0,\sigma^2),\]

as the number of samples \(N\) increases. Here, \(mu\) is the true mean of the random variable, \(\sigma^2\) is the true variance of the random variable. This assumes that each \(\theta^{(i)}\) is an independent sample with unity weight.

When the samples are not independent or the weights are not unity, the variance of the estimator will generally be larger. The "effective sample size" (ESS) describes how many independent samples would have been needed to obtain the same estimator variance. In particular, let \(\tilde{\mu}_N\) be a Monte Carlo estimator based on \(N\) correlated samples.

If method=="Batch" (default) The overlapping batch method (OBM) described in [Flegal2010] is used. This method is also applied to each component independently, resulting in an ESS estimate for each component.

If method=="MultiBatch", The multivariate method of [Vats2019] is employed. This method takes into account the joint correlation of all components of the chain and returns a single ESS. This approach is preferred in high dimensional settings.

Parameters
[in]methodA string specifying the method used to estimate the ESS. Options in SampleCollection are "Batch" or "MultiBatch". The MarkovChain class also supports a "Wolff" method.
Returns
Either a vector of length \(D\) containing an ESS estimate for each component or a vector of length 1 containing a single multivariate ESS estimate.

Reimplemented from muq::SamplingAlgorithms::SampleEstimator.

Reimplemented in muq::SamplingAlgorithms::MarkovChain.

Definition at line 153 of file SampleCollection.h.

References ESS().

Referenced by ESS().

◆ FromMatrix()

std::shared_ptr< SampleCollection > SampleCollection::FromMatrix ( Eigen::MatrixXd const &  samps)
static

Constructs a SampleCollection from a matrix of samples. Each column of the matrix is a sample.

Definition at line 427 of file SampleCollection.cpp.

◆ GetMeta()

std::unordered_map< std::string, Eigen::MatrixXd > SampleCollection::GetMeta ( ) const
protected
Returns
A map from meta data name to a matrix where each column corresponds to a sample

Definition at line 540 of file SampleCollection.cpp.

References at(), and size().

◆ head()

virtual std::shared_ptr<SampleCollection> muq::SamplingAlgorithms::SampleCollection::head ( unsigned int  N) const
inlinevirtual

Returns a new sample collection with the first \(N\) states of this collection. Equivalent to a python command like list[0:N]

Definition at line 88 of file SampleCollection.h.

References segment().

◆ Mean()

Eigen::VectorXd SampleCollection::Mean ( int  blockInd = -1) const
overridevirtual

Computes the sample mean of \(x\) (if blockInd==-1) or \(x_i\) (if blockInd==i).
If blockInd is non-negative, only the mean of one block of the samples is computed.

Reimplemented from muq::SamplingAlgorithms::SampleEstimator.

Reimplemented in muq::SamplingAlgorithms::DistributedCollection.

Definition at line 149 of file SampleCollection.cpp.

References RecursiveSum(), and samples.

Referenced by BatchError(), CentralMoment(), MultiBatchESS(), and RunningCovariance().

◆ MultiBatchError()

Eigen::VectorXd SampleCollection::MultiBatchError ( int  blockInd = -1,
int  batchSize = -1,
int  overlap = -1 
) const
virtual

Estimates the Monte Carlo standard error using the effective sample size coming from the multivariate overlapping batch method of [Vats2019] .

Parameters
[in]blockIndSpecifies the block of the sampling state we're interested in. Defaults to -1, which will result in all blocks of the sampling state being concatenated in the MCSE estimate.
[in]batchSizeThe size of the batches \(b_n\) to use. Defaults to \(N^{1/2}\). Note that to ensure the convergence of the MCSE estimate as \(N\rightarrow\infty\), the batch size must increase with \(N\).
[in]overlapHow much the batches overlap. Defaults to \(0.75*b_n\), where \(b_n\) is the batch size. A nonoverlapping batch estimate will be use is overlap=0. This value bust be smaller than \(b_n\).
Returns
A vector containing \(\sqrt{\sigma_i / N_{ess}}\), where \(N_{ess}\) is the multivariate ESS estimate return by MultiBatchESS.

Definition at line 420 of file SampleCollection.cpp.

References MultiBatchESS(), and muq::SamplingAlgorithms::SampleEstimator::Variance().

Referenced by muq::SamplingAlgorithms::MarkovChain::StandardError(), and StandardError().

◆ MultiBatchESS()

double SampleCollection::MultiBatchESS ( int  blockInd = -1,
int  batchSize = -1,
int  overlap = -1 
) const

Computes the multivariate effective sample size of [Vats2019]. This implementation is adapted from the python implementation at https://github.com/Gabriel-p/multiESS .

The standard definition of a univariate effective sample size is given in terms of the the ratio between the target distribution variance and the estimator variance. More precisely, \(\text{ESS}_u\) is given by

\[ \text{ESS}_u = N \frac{\sigma^2}{\hat{\sigma}^2}, \]

where \(N\) is the total number of samples in the Monte Carlo estimate, \(\sigma^2\) is the variance of the target random variable, and \(\hat{\sigma}^2\) is the variance of a Monte Carlo estimator for the mean of the target distribution.

The multivariate ESS is defined similarly, but using generalized variances. Let \(\Sigma\) denote the target distribution covariance and \(\hat{\Sigma}\) the covariance of the Monte Carlo mean estimator. The multivariate ESS of [Vats2019] is then given by

\[ \text{ESS}_m = N \left( \frac{|\Sigma|}{|\hat{\Sigma}|} \right)^{1/D}, \]

where \(|\cdot|\) denotes the determinant of the matrix and \(D\) is the dimension of target random variable.

As proposed in [Vats2019], it is possible to estimate \(\text{ESS}_m\) using overlapping batches of samples to calculate the estimator covariance \(\hat{\Sigma}\).

Parameters
[in]blockIndThe block of the SampleState that we want to use in the ESS computation.
[in]batchSizeThe size of the batches \(b_n\) to use. Defaults to \(N^{1/3}\). Note that to ensure the convergence of the ESS estimate as \(N\rightarrow\infty\), the batch size must increase with \(N\).
[in]overlapHow much the batches overlap. Defaults to \(b_n/2\), where \(b_n\) is the batch size. A nonoverlapping batch estimate will be use is overlap=0. This value bust be smaller than \(b_n\).
Returns
A double containing the multivariate ESS estimate.

Definition at line 303 of file SampleCollection.cpp.

References muq::SamplingAlgorithms::SampleEstimator::BlockSize(), Covariance(), Mean(), segment(), and size().

Referenced by muq::SamplingAlgorithms::MarkovChain::ESS(), ESS(), and MultiBatchError().

◆ RecursiveSum()

template<typename FuncType >
static std::pair<double,Eigen::VectorXd> muq::SamplingAlgorithms::SampleCollection::RecursiveSum ( std::vector< std::shared_ptr< SamplingState >>::const_iterator  start,
std::vector< std::shared_ptr< SamplingState >>::const_iterator  end,
FuncType &  f 
)
inlinestaticprotected

Definition at line 351 of file SampleCollection.h.

Referenced by CentralMoment(), and Mean().

◆ RecursiveWeightSum()

std::pair< double, double > SampleCollection::RecursiveWeightSum ( std::vector< std::shared_ptr< SamplingState >>::const_iterator  start,
std::vector< std::shared_ptr< SamplingState >>::const_iterator  end 
)
staticprotected

Returns the sum of the weights and the sum of the squared weights.

Definition at line 259 of file SampleCollection.cpp.

◆ Rhat() [1/2]

Eigen::VectorXd SampleCollection::Rhat ( int  blockDim,
unsigned int  numSegments = 4,
boost::property_tree::ptree  options = boost::property_tree::ptree() 
) const
virtual

Writes the samples, weights, and sample metadata to a group in an HDFfile.

  • The samples themselves will be written as a \(DxN\) matrix to the "/samples" dataset.
  • The sample weights will be written as a length \(1xN\) row vector to the "/weights" dataset.
  • Each additional metadata field (e.g., logdensity) will be written to a \(KxN\) matrix, where \(K\) is the dimension of the vector-valued metadata and \(N|f$ is the number of samples. @param[in] filename The name of the HDF file to write to @param[in] group The name of the group within the HDF file. Defaults to the root directory "/". */ virtual void WriteToFile(std::string const& filename, std::string const& group = "/") const; /** @param[in] name Need the this piece of meta data for each sample \return A matrix of meta data associated with in the input string */ Eigen::MatrixXd GetMeta(std::string const& name) const; /** Returns a vector containing the names of metadata that is available. If requireAll is true, then the returned list contains metadata that is available for all samples. If requireAll is false, then this list contains metadata that is stored with any sample even if not all samples have that metadata. */ std::set<std::string> ListMeta(bool requireAll=true) const; /** Using samples of \)@_fakenlx \( stored in this sample collection, this function computes the expected value of \)@_fakenlf(x) \( for some function \)f \( defined as a muq::Modeling::ModPiece. The output is a vector containing the expected value of \)@_fakenlf \(. */ virtual Eigen::VectorXd ExpectedValue(std::shared_ptr<muq::Modeling::ModPiece> const& f, std::vector<std::string> const& metains = std::vector<std::string>()) const override; /** Computes running estimates of the expected value. Returns a std::vector of Eigen::VectorXd instances. Index \)@_fakenln \( of the std::vector contains a Monte Carlo estimate of the expected value of \)@_fakenlf \( using the first \)@_fakenln \( samples of the sample collection. */ std::vector<Eigen::VectorXd> RunningExpectedValue(std::shared_ptr<muq::Modeling::ModPiece> const& f, std::vector<std::string> const& metains = std::vector<std::string>()) const; /** Returns the size \)@_fakenlN_i \( of each block. If `blockInd==-1`, the size \)N \( of the joint random variable is returned. */ virtual unsigned int BlockSize(int blockInd) const override; /** Returns the nubmer of block \)@_fakenlM \(. */ virtual unsigned int NumBlocks() const override; /** To help assess the convergence of a MCMC chain, this function returns either an estimate of the standard scale reduction \)\hat{R} \( diagnostic from \cite Gelman2013, one of the modifications presented in \cite Vehtari2021, the multivariate potential scale reduction factor (MPSRF) of \cite Brooks1998, or multivariate adapations of MPSRF that are similar to the split and ranked methods of \cite Vehtari2021. These scale reduction factors are typically computed using multiple chains. This function instead splits a single chain into multiple segments and then treats each segment as chain in the usual \)\hat{R} \( estimators. This approach is similar to the split techniques of \cite Vehtari2021. For more details on the definition and computation of \)\hat{R} \(, see Diagnostics::Rhat. Parameter Key | Type | Default Value | Description | ------------- | ------------- | ------------- | ------------- | "Transform" | boolean | False | If the parameters should be rank-transformed before computing Rhat, as in \cite Vehtari2021. | "Multivariate" | boolean | False | If the MPSRF value should be returned instead of the componentwise \$\hat{R}\) statistic. If true, the output vector will have a single component. The MPSRF serves as a worse case estimate of \(\hat{R}\) over all linear combinations of the parameters. |
Parameters
[in]blockDimSpecifies the block of the sampling state we're interested in. If blockDim=-1, then all blocks will be concatenated together.
[in]options(optional) A property tree possibly containing settings for the "Multivariate" or "Transform" parameters listed above.
Returns
If Multivariate==False, a vector of \(\hat{R}\) values for each component of the parameters. If Multivariate==true, then a length 1 vector containing the MPSRF.

Definition at line 684 of file SampleCollection.cpp.

Referenced by Rhat().

◆ Rhat() [2/2]

virtual Eigen::VectorXd muq::SamplingAlgorithms::SampleCollection::Rhat ( unsigned int  numSegments = 4,
boost::property_tree::ptree  options = boost::property_tree::ptree() 
) const
inlinevirtual

Same as the other Rhat function, except all blocks are concatenated together.

Definition at line 338 of file SampleCollection.h.

References Rhat().

◆ RunningCovariance() [1/2]

std::vector< Eigen::MatrixXd > SampleCollection::RunningCovariance ( Eigen::VectorXd const &  mean,
int  blockInd = -1 
) const
virtual

Definition at line 203 of file SampleCollection.cpp.

References at(), and size().

◆ RunningCovariance() [2/2]

std::vector< Eigen::MatrixXd > SampleCollection::RunningCovariance ( int  blockInd = -1) const
virtual

Computes running estimates of the covariance.

This function returns a vector of matrices where index \(n\) of the vector contains a Monte Carlo estimate of the covariance that is constructed using the first \(n\) samples in this SampleCollection.

Definition at line 199 of file SampleCollection.cpp.

References Mean().

◆ segment()

std::shared_ptr< SampleCollection > SampleCollection::segment ( unsigned int  startInd,
unsigned int  length,
unsigned int  skipBy = 1 
) const
virtual

Returns a new sample collection containing a segment of this sample collection. The Equivalent to a python command like list[startInd:startInd+length:skipBy]

Reimplemented in muq::SamplingAlgorithms::MarkovChain.

Definition at line 122 of file SampleCollection.cpp.

References at(), and size().

Referenced by BatchError(), head(), MultiBatchESS(), and tail().

◆ size()

◆ StandardError() [1/3]

virtual Eigen::VectorXd muq::SamplingAlgorithms::SampleCollection::StandardError ( int  blockInd) const
inlineoverridevirtual

Reimplemented from muq::SamplingAlgorithms::SampleEstimator.

Reimplemented in muq::SamplingAlgorithms::MarkovChain.

Definition at line 222 of file SampleCollection.h.

References StandardError().

Referenced by StandardError().

◆ StandardError() [2/3]

Eigen::VectorXd SampleCollection::StandardError ( int  blockInd,
std::string const &  method 
) const
overridevirtual

Returns an estimate of the the Monte Carlo standard error (MCSE) \(\hat{\sigma}\) for a Monte Carlo estimate of the mean derived using this SampleEstimator. Recall at the MCSE is the standard deviation of the estimator variance employed in the Central Limit Theorem.

Parameters
[in]blockIndSpecifies the block of the sampling state we're interested in. Defaults to -1, which will result in all blocks of the sampling state being concatenated in the MCSE estimate.
[in]methodA string describing what method should be used to estimate the MCSE. Defaults to "Batch". Options are the same as the ESS function: "Batch", "MultiBatch", and "Wolff". The "Wolff" option is only valid for the "MarkovChain" and "MultiIndexEstimator" classes.
Returns
A vector containing either the MCSE for each component.

Implements muq::SamplingAlgorithms::SampleEstimator.

Reimplemented in muq::SamplingAlgorithms::MarkovChain.

Definition at line 369 of file SampleCollection.cpp.

References BatchError(), and MultiBatchError().

◆ StandardError() [3/3]

virtual Eigen::VectorXd muq::SamplingAlgorithms::SampleCollection::StandardError ( std::string const &  method = "Batch") const
inlineoverridevirtual

Reimplemented from muq::SamplingAlgorithms::SampleEstimator.

Reimplemented in muq::SamplingAlgorithms::MarkovChain.

Definition at line 223 of file SampleCollection.h.

References StandardError().

Referenced by StandardError().

◆ tail()

virtual std::shared_ptr<SampleCollection> muq::SamplingAlgorithms::SampleCollection::tail ( unsigned int  N) const
inlinevirtual

Returns a new sample collection with the last \(N\) states of this collection. Equivalent to a python command like list[-N:]

Definition at line 93 of file SampleCollection.h.

References segment(), and size().

◆ Weights()

Eigen::VectorXd SampleCollection::Weights ( ) const
virtual

Returns a vector of unnormalized sample weights for computing empirical expectations. If the samples were generated with MC or MCMC, the weights will all be

  1. Importance sampling will generate non-unit weights. Note that these weights are unnormalized and will not generally sum to one.

Reimplemented in muq::SamplingAlgorithms::DistributedCollection.

Definition at line 466 of file SampleCollection.cpp.

References samples.

Member Data Documentation

◆ samples

std::vector<std::shared_ptr<SamplingState> > muq::SamplingAlgorithms::SampleCollection::samples
protected

Definition at line 343 of file SampleCollection.h.

Referenced by Add(), AsMatrix(), at(), back(), CentralMoment(), Covariance(), Mean(), size(), and Weights().


The documentation for this class was generated from the following files: