MUQ  0.4.3
muq::SamplingAlgorithms::MultiIndexEstimator Class Reference

Class for estimating expectations using multi-index samples from MC or MCMC. More...

#include <MultiIndexEstimator.h>

Inheritance diagram for muq::SamplingAlgorithms::MultiIndexEstimator:

Detailed Description

Class for estimating expectations using multi-index samples from MC or MCMC.

Definition at line 15 of file MultiIndexEstimator.h.

Public Member Functions

 MultiIndexEstimator (std::vector< std::shared_ptr< MIMCMCBox >> const &boxesIn, bool useQoisIn=false)
 
virtual ~MultiIndexEstimator ()=default
 
virtual unsigned int BlockSize (int blockInd) const override
 
virtual unsigned int NumBlocks () const override
 
virtual Eigen::VectorXd ExpectedValue (std::shared_ptr< muq::Modeling::ModPiece > const &f, std::vector< std::string > const &metains=std::vector< std::string >()) 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 Eigen::VectorXd StandardError (int blockDim, std::string const &method) const override
 
virtual Eigen::VectorXd StandardError (std::string const &method="Batch") const override
 
virtual Eigen::VectorXd StandardError (int blockDim) const override
 
virtual Eigen::VectorXd ESS (std::string const &method="Batch") const override
 
virtual Eigen::VectorXd ESS (int blockDim) const override
 
virtual Eigen::VectorXd ESS (int blockDim, std::string const &method) const override
 
- Public Member Functions inherited from muq::SamplingAlgorithms::SampleEstimator
virtual ~SampleEstimator ()=default
 
virtual Eigen::VectorXd CentralMoment (unsigned int order, int blockNum=-1) const
 
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 Mean (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
 

Constructor & Destructor Documentation

◆ MultiIndexEstimator()

MultiIndexEstimator::MultiIndexEstimator ( std::vector< std::shared_ptr< MIMCMCBox >> const &  boxesIn,
bool  useQoisIn = false 
)

Construct the multiindex estimator using MIMCMC boxes. These boxes are typically constructed by a MIMCMC methods such as the GreedyMLMCMC or MIMCMC classes.

Parameters
[in]boxesIn"Boxes" holding the differences between chains at different indices
[in]useQoisIn(optional) Whether this estimator should use the QOIs in the chains or the parameters themselves. Defaults to false, which implies the parameters will be used in the estimates.

Definition at line 7 of file MultiIndexEstimator.cpp.

◆ ~MultiIndexEstimator()

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

Member Function Documentation

◆ BlockSize()

unsigned int MultiIndexEstimator::BlockSize ( int  blockInd) const
overridevirtual

Returns the size \(N_i\) of each block. If blockInd==-1, the size \(N\) of the joint random variable is returned.

Implements muq::SamplingAlgorithms::SampleEstimator.

Definition at line 15 of file MultiIndexEstimator.cpp.

References blockSizes.

Referenced by Covariance().

◆ Covariance() [1/2]

Eigen::MatrixXd MultiIndexEstimator::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 157 of file MultiIndexEstimator.cpp.

References BlockSize(), boxes, and useQois.

◆ Covariance() [2/2]

virtual Eigen::MatrixXd muq::SamplingAlgorithms::MultiIndexEstimator::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.

Definition at line 39 of file MultiIndexEstimator.h.

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

◆ ESS() [1/3]

virtual Eigen::VectorXd muq::SamplingAlgorithms::MultiIndexEstimator::ESS ( int  blockDim) const
inlineoverridevirtual

Reimplemented from muq::SamplingAlgorithms::SampleEstimator.

Definition at line 69 of file MultiIndexEstimator.h.

References ESS().

Referenced by ESS().

◆ ESS() [2/3]

Eigen::VectorXd MultiIndexEstimator::ESS ( int  blockInd,
std::string const &  method 
) const
overridevirtual

Returns an estimate of the effective sample size (ESS), which is the number of independent samples of the target distribution that would be needed to obtain the same statistical accuracy as this estimator. For independent samples, the estimator variance (squared MCSE) \(\hat{\sigma}^2\f=\sigma^2/ N\). Given the estimator variance \(\hat{\sigma}^2\), the effective sample size is then given by the ratio of the target distribution variance and the estimator variance:

\[ \text{ESS} = \frac{\sigma^2}{\hat{\sigma}^2}. \]

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". Other options include "MultiBatch" or "Wolff". For details, see the SampleCollection class.
Returns
A vector containing either the MCSE for each component (if method!="MultiBatch") or a single component vector containing the square root of the generalized estimator variance (if method=="MultiBatch").

Implements muq::SamplingAlgorithms::SampleEstimator.

Definition at line 43 of file MultiIndexEstimator.cpp.

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

◆ ESS() [3/3]

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

This function returns an estimate of the Effective Sample Size (ESS) of this estimator. The ESS here refers to the number of indepedent samples that would be required in a classic single-level Monte Carlo estimate to achieve the same statistical accuracy as this multi-index estimator.

This function computes the ESS by computing the ratio of the sample variance with the squared MCSE. The MCSE is computed by the MultiIndex::StandardError function.

Parameters
[in]methodSpecifies the type of MCSE estimator used by the single level chains to estimate the variance of each term in the multiindex telescoping series. See MarkovChain::ESS for more details.
Returns
If method=="MultiBatch", this function returns length 1 vector containing a single multivariate effective sample size estimate. If method!="MultiBatch", this function returns an ESS for each component of the chain.

Reimplemented from muq::SamplingAlgorithms::SampleEstimator.

Definition at line 68 of file MultiIndexEstimator.h.

References ESS().

Referenced by ESS().

◆ ExpectedValue()

Eigen::VectorXd MultiIndexEstimator::ExpectedValue ( std::shared_ptr< muq::Modeling::ModPiece > const &  f,
std::vector< std::string > const &  metains = std::vector< std::string >() 
) const
overridevirtual

Using samples of \(x\) stored in this sample collection, this function approximates the expected value of \(f(x)\) for some function \(f\) defined as a muq::Modeling::ModPiece. The output is a vector containing the expected value of \(f\).

Implements muq::SamplingAlgorithms::SampleEstimator.

Definition at line 115 of file MultiIndexEstimator.cpp.

References boxes, and useQois.

◆ GetDiffChains()

std::vector< std::shared_ptr< MarkovChain > > MultiIndexEstimator::GetDiffChains ( int  blockInd = -1,
std::shared_ptr< muq::Modeling::ModPiece > const &  f = nullptr 
) const
private

Creates a Markov chain for each term in the telescoping sum.

Parameters
[in]blockInd(optional) The index of the block we're interested in. If -1, then all blocks in the state are concatenated. This is not used if the f argument is specified.
[in]f(optional), A ModPiece that evaluates a quantity of interest if we're interested in chains over the QOI difference.
Returns
A vector of Markov chains for each term in the series.

Definition at line 50 of file MultiIndexEstimator.cpp.

References boxes, diffChains, and useQois.

Referenced by StandardError().

◆ NumBlocks()

unsigned int MultiIndexEstimator::NumBlocks ( ) const
overridevirtual

Returns the nubmer of block \(M\).

Implements muq::SamplingAlgorithms::SampleEstimator.

Definition at line 25 of file MultiIndexEstimator.cpp.

References blockSizes.

◆ StandardError() [1/3]

virtual Eigen::VectorXd muq::SamplingAlgorithms::MultiIndexEstimator::StandardError ( int  blockDim) const
inlineoverridevirtual

Reimplemented from muq::SamplingAlgorithms::SampleEstimator.

Definition at line 56 of file MultiIndexEstimator.h.

References StandardError().

Referenced by StandardError().

◆ StandardError() [2/3]

Eigen::VectorXd MultiIndexEstimator::StandardError ( int  blockDim,
std::string const &  method 
) const
overridevirtual

Computes the standard deviation mean value returned by theis MultiIndexEstimator. The estimator variances for each term in the telescoping series are summed and the square root of this quantity is returned. This process assumes that the terms in the series are independent. The value of Method is passed on to the underlying SampleCollection classes to compute the variance of each term. Valid options are Batch, MultiBatch, and Wolff.

Parameters
[in]blockDimSpecifies the block that we wish to use in the MCSE estimator. Defaults to -1, which results in all blocks of the chain being concatenated in the MCSE estimate.
[in]methodSpecifies the type of MCSE estimator used by the single level chains to estimate the variance of each term in the multiindex telescoping series. See MarkovChain::StandardError for more details.
Returns
If method=="MultiBatch", this function returns length 1 vector containing a single multivariate MCSE estimate. If method!="MultiBatch", this function returns the MCSE for each component of the chain.

Implements muq::SamplingAlgorithms::SampleEstimator.

Definition at line 30 of file MultiIndexEstimator.cpp.

References GetDiffChains().

Referenced by ESS().

◆ StandardError() [3/3]

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

Reimplemented from muq::SamplingAlgorithms::SampleEstimator.

Definition at line 55 of file MultiIndexEstimator.h.

References StandardError().

Referenced by StandardError().

Member Data Documentation

◆ blockSizes

const Eigen::VectorXi muq::SamplingAlgorithms::MultiIndexEstimator::blockSizes
private

Definition at line 84 of file MultiIndexEstimator.h.

Referenced by BlockSize(), and NumBlocks().

◆ boxes

std::vector<std::shared_ptr<MIMCMCBox> > muq::SamplingAlgorithms::MultiIndexEstimator::boxes
private

Definition at line 87 of file MultiIndexEstimator.h.

Referenced by Covariance(), ExpectedValue(), and GetDiffChains().

◆ diffChains

std::vector<std::shared_ptr<MarkovChain> > muq::SamplingAlgorithms::MultiIndexEstimator::diffChains
private

Definition at line 88 of file MultiIndexEstimator.h.

Referenced by GetDiffChains().

◆ useQois

const bool muq::SamplingAlgorithms::MultiIndexEstimator::useQois
private

Definition at line 85 of file MultiIndexEstimator.h.

Referenced by Covariance(), ExpectedValue(), and GetDiffChains().


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