MUQ  0.4.3
AMProposal.h
Go to the documentation of this file.
1 #ifndef AMPROPOSAL_H_
2 #define AMPROPOSAL_H_
3 
5 
6 namespace muq {
7  namespace SamplingAlgorithms {
8 
21  class AMProposal : public MCMCProposal {
22  public:
23 
24  AMProposal(boost::property_tree::ptree pt,
25  std::shared_ptr<AbstractSamplingProblem> const& prob);
26 
33  AMProposal(boost::property_tree::ptree pt,
34  std::shared_ptr<AbstractSamplingProblem> const& prob,
35  Eigen::MatrixXd const& initialCov);
36 
37 
38  virtual ~AMProposal() = default;
39 
41 
46  virtual void Adapt(unsigned int const t, std::vector<std::shared_ptr<SamplingState>> const& states) override;
47 
49  virtual Eigen::MatrixXd ProposalCovariance() const;
50 
51  private:
52 
53  static Eigen::MatrixXd ConstructCovariance(boost::property_tree::ptree const& pt ,
54  std::shared_ptr<AbstractSamplingProblem> const& prob);
55 
56 
57  virtual std::shared_ptr<SamplingState> Sample(std::shared_ptr<SamplingState> const& currentState) override;
58 
59  virtual double LogDensity(std::shared_ptr<SamplingState> const& currState,
60  std::shared_ptr<SamplingState> const& propState) override;
61 
63  Eigen::MatrixXd newSamps;
64  unsigned int numNewSamps=0;
65 
67  unsigned int numAdaptSamps=0;
68 
70  Eigen::VectorXd mean;
71 
73  Eigen::MatrixXd propCov;
74 
76  Eigen::LLT<Eigen::MatrixXd, Eigen::Lower> propChol;
77 
79  const unsigned int adaptSteps;
80 
82  const unsigned int adaptStart;
83 
85  const unsigned int adaptEnd;
86 
87  // Multiplicative scaling of the sample covariance
88  const double adaptScale;
89 
90  };
91 
92  } // namespace SamplingAlgorithms
93 } // namespace muq
94 
95 #endif
An implemental of the Adaptive Metropolis algorithm.
Definition: AMProposal.h:21
virtual std::shared_ptr< SamplingState > Sample(std::shared_ptr< SamplingState > const &currentState) override
Definition: AMProposal.cpp:84
Eigen::LLT< Eigen::MatrixXd, Eigen::Lower > propChol
The Cholesky factorization of the current proposal covariance.
Definition: AMProposal.h:76
virtual Eigen::MatrixXd ProposalCovariance() const
Definition: AMProposal.cpp:79
AMProposal(boost::property_tree::ptree pt, std::shared_ptr< AbstractSamplingProblem > const &prob, Eigen::MatrixXd const &initialCov)
AMProposal(boost::property_tree::ptree pt, std::shared_ptr< AbstractSamplingProblem > const &prob)
Eigen::VectorXd mean
The current mean.
Definition: AMProposal.h:70
static Eigen::MatrixXd ConstructCovariance(boost::property_tree::ptree const &pt, std::shared_ptr< AbstractSamplingProblem > const &prob)
Definition: AMProposal.cpp:34
Eigen::MatrixXd propCov
The current proposal covariance.
Definition: AMProposal.h:73
Eigen::MatrixXd newSamps
Samples seen between visits.
Definition: AMProposal.h:63
virtual void Adapt(unsigned int const t, std::vector< std::shared_ptr< SamplingState >> const &states) override
Adapt the proposal after each step.
Definition: AMProposal.cpp:43
const unsigned int adaptStart
When should we start adapting?
Definition: AMProposal.h:82
const unsigned int adaptEnd
After how many iterations should we stop adapting?
Definition: AMProposal.h:85
unsigned int numAdaptSamps
Total number of samples that have been used so far to compute mean and covariance.
Definition: AMProposal.h:67
const unsigned int adaptSteps
How frequently should we update the adaption?
Definition: AMProposal.h:79
virtual double LogDensity(std::shared_ptr< SamplingState > const &currState, std::shared_ptr< SamplingState > const &propState) override
Definition: AMProposal.cpp:97
std::shared_ptr< AbstractSamplingProblem > prob
Definition: MCMCProposal.h:81