MUQ  0.4.3
MultiIndexEstimator.h
Go to the documentation of this file.
1 #ifndef MULTIINDEXESTIMATOR_H
2 #define MULTIINDEXESTIMATOR_H
3 
7 
8 namespace muq{
9 namespace SamplingAlgorithms{
10 
16  {
17  public:
18 
27  MultiIndexEstimator(std::vector<std::shared_ptr<MIMCMCBox>> const& boxesIn,
28  bool useQoisIn = false);
29 
30  virtual ~MultiIndexEstimator() = default;
31 
32  virtual unsigned int BlockSize(int blockInd) const override;
33 
34  virtual unsigned int NumBlocks() const override;
35 
36  virtual Eigen::VectorXd ExpectedValue(std::shared_ptr<muq::Modeling::ModPiece> const& f,
37  std::vector<std::string> const& metains = std::vector<std::string>()) const override;
38 
39  virtual Eigen::MatrixXd Covariance(int blockInd=-1) const override { return SampleEstimator::Covariance(blockInd);};
40  virtual Eigen::MatrixXd Covariance(Eigen::VectorXd const& mean,
41  int blockInd=-1) const override;
42 
43 
44 
54  virtual Eigen::VectorXd StandardError(int blockDim, std::string const& method) const override;
55  virtual Eigen::VectorXd StandardError(std::string const& method="Batch") const override{return StandardError(-1,method);};
56  virtual Eigen::VectorXd StandardError(int blockDim) const override{return StandardError(blockDim,"Batch");};
57 
68  virtual Eigen::VectorXd ESS(std::string const& method="Batch") const override{return ESS(-1,method);};
69  virtual Eigen::VectorXd ESS(int blockDim) const override{return ESS(blockDim,"Batch");};
70  virtual Eigen::VectorXd ESS(int blockDim, std::string const& method) const override;
71 
72  private:
73 
79  std::vector<std::shared_ptr<MarkovChain>> GetDiffChains(int blockInd=-1,
80  std::shared_ptr<muq::Modeling::ModPiece> const& f=nullptr) const;
81 
82 
83 
84  const Eigen::VectorXi blockSizes;
85  const bool useQois;
86 
87  std::vector<std::shared_ptr<MIMCMCBox>> boxes;
88  std::vector<std::shared_ptr<MarkovChain>> diffChains;
89 
90  };
91 }
92 }
93 
94 
95 
96 #endif
Class for estimating expectations using multi-index samples from MC or MCMC.
virtual unsigned int BlockSize(int blockInd) const override
std::vector< std::shared_ptr< MarkovChain > > GetDiffChains(int blockInd=-1, std::shared_ptr< muq::Modeling::ModPiece > const &f=nullptr) const
MultiIndexEstimator(std::vector< std::shared_ptr< MIMCMCBox >> const &boxesIn, bool useQoisIn=false)
virtual unsigned int NumBlocks() const override
virtual Eigen::MatrixXd Covariance(int blockInd=-1) const override
virtual Eigen::VectorXd StandardError(std::string const &method="Batch") 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::VectorXd ESS(std::string const &method="Batch") const override
virtual Eigen::VectorXd StandardError(int blockDim) const override
virtual Eigen::VectorXd StandardError(int blockDim, std::string const &method) const override
virtual Eigen::VectorXd ESS(int blockDim) const override
std::vector< std::shared_ptr< MarkovChain > > diffChains
std::vector< std::shared_ptr< MIMCMCBox > > boxes
Abstract base class for computing sample-based approximations of expectations.
virtual Eigen::MatrixXd Covariance(int blockInd=-1) const