MUQ  0.4.3
muq::SamplingAlgorithms::SampleEstimator Class Referenceabstract

Abstract base class for computing sample-based approximations of expectations. More...

#include <SampleEstimator.h>

Inheritance diagram for muq::SamplingAlgorithms::SampleEstimator:

Detailed Description

Abstract base class for computing sample-based approximations of expectations.

Consider a vector-valued random variable $x$ taking values in \(\mathbb{R}^N\). Now consider a decomposition of this random variable into \(M\) "blocks" \(x_1, x_2, \ldots, x_M\) with lengths \(N_1,\ldots, N_M\), where \(\sum N_i = N\).
This class provides an abstract interface for methods that approximate expectations of the form

\[ \mathbb{E}[f(x)] = \int f(x) p(x) dx \]

or

\[ \mathbb{E}[f(x_i)] = \int f(x_i) p(x_i) dx, \]

where \(p(x)\) and \(p(x_i)\) are the distributions of \(x\) and \(x_i\), respectively. This class of expectations include the mean \(\mu=\mathbb{E}[x]\) and variance \(\sigma^2 = \mathbb{E}[(x-\mu)^2]\).

In Monte Carlo, samples \(x^{(k)}\) of the random variable are used to approximate the expectation via the law of large numbers:

\[ \mathbb{E}[f(x)] \approx \frac{1}{K} \sum_{k=1}^K f\left(x^{(k)}\right). \]

The muq::SamplingAlgorithms::SampleCollection class is a child of this class that provides such Monte Carlo estimators. Of course, it is possible to construct other estimators of the expection, using multilevel Monte Carlo methods for example. This class aims to provide a common interface for all such approaches.

Definition at line 37 of file SampleEstimator.h.

Public Member Functions

virtual ~SampleEstimator ()=default
 
virtual unsigned int BlockSize (int blockInd) const =0
 
virtual unsigned int NumBlocks () const =0
 
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
 
virtual Eigen::MatrixXd Covariance (int blockInd=-1) const
 
virtual Eigen::MatrixXd Covariance (Eigen::VectorXd const &mean, int blockInd=-1) const =0
 
virtual Eigen::VectorXd ExpectedValue (std::shared_ptr< muq::Modeling::ModPiece > const &f, std::vector< std::string > const &metains=std::vector< std::string >()) const =0
 
virtual Eigen::VectorXd StandardError (int blockInd, std::string const &method) const =0
 
virtual Eigen::VectorXd StandardError (std::string const &method="Batch") const
 
virtual Eigen::VectorXd StandardError (int blockDim) const
 
virtual Eigen::VectorXd ESS (int blockInd, std::string const &method) const =0
 
virtual Eigen::VectorXd ESS (int blockDim) const
 
virtual Eigen::VectorXd ESS (std::string const &method="Batch") const
 

Constructor & Destructor Documentation

◆ ~SampleEstimator()

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

Member Function Documentation

◆ BlockSize()

virtual unsigned int muq::SamplingAlgorithms::SampleEstimator::BlockSize ( int  blockInd) const
pure virtual

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

Implemented in muq::SamplingAlgorithms::MultiIndexEstimator.

Referenced by muq::SamplingAlgorithms::SampleCollection::BatchError(), and muq::SamplingAlgorithms::SampleCollection::MultiBatchESS().

◆ CentralMoment() [1/2]

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

Compute the central moment, as in the other SampleEstimator::CentralMoment function, but use a precomputed (or known) mean. Note that using a vector of zeros for the mean allows non-central moments to be computed.

Parameters
[in]orderThe order \(p\) of the central moment. \(p=2\) yields the variance.
[in]meanA vector containing the mean of \(x\) (if blockNum==-1) or \(x_i\) (if blockNum==i).
[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.

Definition at line 16 of file SampleEstimator.cpp.

◆ CentralMoment() [2/2]

Eigen::VectorXd SampleEstimator::CentralMoment ( unsigned int  order,
int  blockNum = -1 
) const
virtual

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 in muq::SamplingAlgorithms::SampleCollection.

Definition at line 9 of file SampleEstimator.cpp.

◆ Covariance() [1/2]

virtual Eigen::MatrixXd muq::SamplingAlgorithms::SampleEstimator::Covariance ( Eigen::VectorXd const &  mean,
int  blockInd = -1 
) const
pure virtual

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.

Implemented in muq::SamplingAlgorithms::SampleCollection, and muq::SamplingAlgorithms::MultiIndexEstimator.

◆ Covariance() [2/2]

Eigen::MatrixXd SampleEstimator::Covariance ( int  blockInd = -1) const
virtual

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 in muq::SamplingAlgorithms::SampleCollection, muq::SamplingAlgorithms::MultiIndexEstimator, and muq::SamplingAlgorithms::DistributedCollection.

Definition at line 152 of file SampleEstimator.cpp.

Referenced by muq::SamplingAlgorithms::MultiIndexEstimator::Covariance(), and muq::SamplingAlgorithms::SampleCollection::Covariance().

◆ ESS() [1/3]

virtual Eigen::VectorXd muq::SamplingAlgorithms::SampleEstimator::ESS ( int  blockDim) const
inlinevirtual

◆ ESS() [2/3]

virtual Eigen::VectorXd muq::SamplingAlgorithms::SampleEstimator::ESS ( int  blockInd,
std::string const &  method 
) const
pure virtual

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").

Implemented in muq::SamplingAlgorithms::SampleCollection, muq::SamplingAlgorithms::MultiIndexEstimator, and muq::SamplingAlgorithms::MarkovChain.

◆ ESS() [3/3]

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

◆ ExpectedValue()

virtual Eigen::VectorXd muq::SamplingAlgorithms::SampleEstimator::ExpectedValue ( std::shared_ptr< muq::Modeling::ModPiece > const &  f,
std::vector< std::string > const &  metains = std::vector< std::string >() 
) const
pure virtual

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\).

Implemented in muq::SamplingAlgorithms::MultiIndexEstimator, and muq::SamplingAlgorithms::DistributedCollection.

◆ Kurtosis() [1/3]

Eigen::VectorXd SampleEstimator::Kurtosis ( Eigen::VectorXd const &  mean,
Eigen::VectorXd const &  stdDev,
int  blockInd = -1 
) const
virtual

Definition at line 144 of file SampleEstimator.cpp.

◆ Kurtosis() [2/3]

Eigen::VectorXd SampleEstimator::Kurtosis ( Eigen::VectorXd const &  mean,
int  blockInd = -1 
) const
virtual

Evaluate the kurtosis using a precomputed (or known) mean vector.

Definition at line 138 of file SampleEstimator.cpp.

◆ Kurtosis() [3/3]

Eigen::VectorXd SampleEstimator::Kurtosis ( int  blockInd = -1) const
virtual

The marginal kurtosis of a random variable is given by

\[ \tilde{\mu}_4 = \mathbb{E}\left[ \left(\frac{x-\mu}{\sigma}\right)^4 \right], \]

where \(\mu\) is the mean and \(\sigma\) is the standard deviation. This function returns an estimate of this quantity by using sample approximations of \(\mu\), \(\sigma\), and the outer expectation.

In the default implementation, three calls to the ExpectedValue function are used to compute this quantity. Children of this class, like the SampleCollection, may provide more efficient implementations.

Just like the variance is the diagonal of the covariance matrix, this function returns a vector representing the "diagonal" of the full fourth order kurtosis tensor.

Definition at line 133 of file SampleEstimator.cpp.

◆ Mean()

Eigen::VectorXd SampleEstimator::Mean ( int  blockInd = -1) const
virtual

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 in muq::SamplingAlgorithms::SampleCollection, and muq::SamplingAlgorithms::DistributedCollection.

Definition at line 76 of file SampleEstimator.cpp.

◆ NumBlocks()

virtual unsigned int muq::SamplingAlgorithms::SampleEstimator::NumBlocks ( ) const
pure virtual

Returns the nubmer of block \(M\).

Implemented in muq::SamplingAlgorithms::MultiIndexEstimator.

◆ Skewness() [1/3]

Eigen::VectorXd SampleEstimator::Skewness ( Eigen::VectorXd const &  mean,
Eigen::VectorXd const &  stdDev,
int  blockInd = -1 
) const
virtual

Definition at line 126 of file SampleEstimator.cpp.

◆ Skewness() [2/3]

Eigen::VectorXd SampleEstimator::Skewness ( Eigen::VectorXd const &  mean,
int  blockInd = -1 
) const
virtual

Evaluate the skewness using a precomputed (or known) mean vector.

Definition at line 120 of file SampleEstimator.cpp.

◆ Skewness() [3/3]

Eigen::VectorXd SampleEstimator::Skewness ( int  blockInd = -1) const
virtual

The marginal skewness of a random variable is given by

\[ \tilde{\mu}_3 = \mathbb{E}\left[ \left(\frac{x-\mu}{\sigma}\right)^3 \right], \]

where \(\mu\) is the mean and \(\sigma\) is the standard deviation. This function returns an estimate of this quantity by using sample approximations of \(\mu\), \(\sigma\), and the outer expectation.

In the default implementation, three calls to the ExpectedValue function are used to compute this quantity. Children of this class, like the SampleCollection, may provide more efficient implementations.

Just like the variance is the diagonal of the covariance matrix, this function returns a vector representing the "diagonal" of the full third order skewness tensor.

Definition at line 115 of file SampleEstimator.cpp.

◆ StandardError() [1/3]

virtual Eigen::VectorXd muq::SamplingAlgorithms::SampleEstimator::StandardError ( int  blockDim) const
inlinevirtual

◆ StandardError() [2/3]

virtual Eigen::VectorXd muq::SamplingAlgorithms::SampleEstimator::StandardError ( int  blockInd,
std::string const &  method 
) const
pure virtual

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"
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").

Implemented in muq::SamplingAlgorithms::SampleCollection, muq::SamplingAlgorithms::MultiIndexEstimator, and muq::SamplingAlgorithms::MarkovChain.

◆ StandardError() [3/3]

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

◆ StandardizedMoment() [1/3]

Eigen::VectorXd SampleEstimator::StandardizedMoment ( unsigned int  order,
Eigen::VectorXd const &  mean,
Eigen::VectorXd const &  stdDev,
int  blockInd = -1 
) const
virtual

Computes the standardized moment with precomputed (or known) mean and standard deviation.

Definition at line 105 of file SampleEstimator.cpp.

◆ StandardizedMoment() [2/3]

Eigen::VectorXd SampleEstimator::StandardizedMoment ( unsigned int  order,
Eigen::VectorXd const &  mean,
int  blockInd = -1 
) const
virtual

The standardize moment of order $p$ is similar to the central moment, but also includes a scaling of the random variable \(x\) by the standard deviation. Mathematially, the standardized moment is given by

\[ \mathbb{E}\left[ \left(\frac{x-\mu}{\sigma}\right)^p\right], \]

where \(\mu=\mathbb{E}[x]\) is the mean and $\sigma$ is the standard deviation of $x$.

Parameters
[in]orderThe order \(p\) of the central moment. \(p=2\) yields the variance.
[in]meanA vector containing a precomputed (or known) estimate of \(\mu=\mathbb{E}[x]\).
[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 standardized moment.

Definition at line 98 of file SampleEstimator.cpp.

◆ StandardizedMoment() [3/3]

Eigen::VectorXd SampleEstimator::StandardizedMoment ( unsigned int  order,
int  blockInd = -1 
) const
virtual

The standardize moment of order $p$ is similar to the central moment, but also includes a scaling of the random variable \(x\) by the standard deviation. Mathematially, the standardized moment is given by

\[ \mathbb{E}\left[ \left(\frac{x-\mu}{\sigma}\right)^p\right], \]

where \(\mu=\mathbb{E}[x]\) is the mean and $\sigma$ is the standard deviation of $x$.

In the default implementation of this function, three calls to the ExpectedValue function will be made to compute the mean, compute the variance, and then to compute the outer expectation.

Note that the standardized moment of order \(p=3\) is commonly called the skewness and the standardized moment of order $p=4$ is called the kurtosis. Shortcuts for these common moments can be found in the SampleEstimator::Skewness and SampleEstimator::Kurtosis functions.

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 standardized moment.

Definition at line 93 of file SampleEstimator.cpp.

◆ Variance() [1/2]

Eigen::VectorXd SampleEstimator::Variance ( Eigen::VectorXd const &  mean,
int  blockInd = -1 
) const
virtual

Computes the sample variance using a precomputed (or known) mean vector. If blockInd is non-negative, only the mean of one block of the samples is computed.

Definition at line 87 of file SampleEstimator.cpp.

◆ Variance() [2/2]

Eigen::VectorXd SampleEstimator::Variance ( int  blockInd = -1) const
virtual

Computes an estiamte of the variance of \(x\) (if blockInd==-1) or \(x_i\) (if blockInd==i). The variance is defined as

\[ \text{Var}[x] = \mathbb{E}\left[\left(x-\mathbb{E}[x]\right)^2\right]. \]

If blockInd is non-negative, only the variance of one block, i.e. \(\text{Var}[x_i]\) of the samples is computed.

Note that the default implementation of this function makes two calls to SampleEstimator::ExpectedValue. One call is used to compute the mean \(\mathbb{E}[x]\) and another call is used to evaluate the outer expectation \(\mathbb{E}\left[\left(x-\mathbb{E}[x]\right)^2\right]\). In some cases this can lead to biased estimates of the variance.

The SamplingAlgorithms::SampleCollection child of this class provides an unbiased implementation for use with standard Monte Carlo and Markov chain Monte Carlo samplers.

Reimplemented in muq::SamplingAlgorithms::DistributedCollection.

Definition at line 82 of file SampleEstimator.cpp.

Referenced by muq::SamplingAlgorithms::SampleCollection::BatchESS(), muq::SamplingAlgorithms::MultiIndexEstimator::ESS(), muq::SamplingAlgorithms::SampleCollection::MultiBatchError(), and muq::SamplingAlgorithms::MarkovChain::WolffError().


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