MUQ  0.4.3
ExpensiveSamplingProblem.h
Go to the documentation of this file.
1 #ifndef EXPENSIVESAMPLINGPROBLEM_H_
2 #define EXPENSIVESAMPLINGPROBLEM_H_
3 
4 #include "MUQ/config.h"
5 
6 #if MUQ_HAS_PARCER
7 #include <parcer/Communicator.h>
8 #endif
9 
11 
13 
14 namespace muq {
15  namespace SamplingAlgorithms {
17  public:
18 
25  ExpensiveSamplingProblem(std::shared_ptr<muq::Modeling::ModPiece> const& target, std::shared_ptr<muq::Modeling::ModPiece> const& qoi, Eigen::VectorXd const& centroid, boost::property_tree::ptree pt);
26 
32  ExpensiveSamplingProblem(std::shared_ptr<muq::Modeling::ModPiece> const& target, std::shared_ptr<muq::Modeling::ModPiece> const& qoi, boost::property_tree::ptree pt);
33 
38  ExpensiveSamplingProblem(std::shared_ptr<muq::Modeling::ModPiece> const& target, boost::property_tree::ptree pt);
39 
45  ExpensiveSamplingProblem(std::shared_ptr<muq::Modeling::ModPiece> const& target, Eigen::VectorXd const& centroid, boost::property_tree::ptree pt);
46 
47 #if MUQ_HAS_PARCER
53  ExpensiveSamplingProblem(std::shared_ptr<muq::Modeling::ModPiece> const& target, boost::property_tree::ptree pt, std::shared_ptr<parcer::Communicator> comm);
54 
61  ExpensiveSamplingProblem(std::shared_ptr<muq::Modeling::ModPiece> const& target, Eigen::VectorXd const& centroid, boost::property_tree::ptree pt, std::shared_ptr<parcer::Communicator> comm);
62 #endif
63 
64  virtual ~ExpensiveSamplingProblem() = default;
65 
66  virtual double LogDensity(std::shared_ptr<SamplingState> const& state) override;
67 
68  unsigned int CacheSize() const;
69 
70 
71  virtual void AddOptions(boost::property_tree::ptree & pt) const override;
72 
73  virtual std::shared_ptr<SamplingState> QOI() override;
74 
75  private:
76 
78 
81  void SetUp(boost::property_tree::ptree& pt);
82 
90  double RefineSurrogate(
91  unsigned int const step,
92  std::shared_ptr<SamplingState> const& state,
93  std::vector<Eigen::VectorXd>& neighbors,
94  std::vector<Eigen::VectorXd>& results);
95 
103  double RefineSurrogate(
104  std::shared_ptr<SamplingState> const& state,
105  double const radius,
106  std::vector<Eigen::VectorXd>& neighbors,
107  std::vector<Eigen::VectorXd>& results);
108 
116  void RefineSurrogate(
117  Eigen::VectorXd const& point,
118  unsigned int const index,
119  std::vector<Eigen::VectorXd>& neighbors,
120  std::vector<Eigen::VectorXd>& results);
121 
123 
129  void CheckNeighbors(
130  std::shared_ptr<SamplingState> const& state,
131  std::vector<Eigen::VectorXd>& neighbors,
132  std::vector<Eigen::VectorXd>& results) const;
133 
135 
138  void CheckNumNeighbors(std::shared_ptr<SamplingState> const& state);
139 
141 
145  double ErrorThreshold(Eigen::VectorXd const& x) const;
146 
148 
153  double LogLyapunovFunction(Eigen::VectorXd const& x) const;
154 
156  std::shared_ptr<muq::Approximation::LocalRegression> reg;
157 
159  std::shared_ptr<muq::Approximation::LocalRegression> regQoI;
160 
162 
165  std::pair<double, double> beta;
166 
168 
171  double tau0;
172 
174 
177  std::pair<double, double> gamma;
178 
180  std::pair<double, double> lyapunovPara;
181 
183  const Eigen::VectorXd centroid;
184 
186  double eta;
187 
189  unsigned int level = 1;
190 
192  unsigned int step = 0;
193 public:
195  unsigned int cumbeta = 0;
196 
198  unsigned int cumgamma = 0;
199 
201  unsigned int cumkappa = 0;
202 
204  double lastLyapunov = 0.0;
205 
207  double lastThreshold = 0.0;
208  };
209  } // namespace SamplingAlgorithms
210 } // namespace muq
211 
212 #endif
std::pair< double, double > beta
Parameters for random refinement.
ExpensiveSamplingProblem(std::shared_ptr< muq::Modeling::ModPiece > const &target, std::shared_ptr< muq::Modeling::ModPiece > const &qoi, boost::property_tree::ptree pt)
unsigned int step
The current step in the MCMC chain.
std::pair< double, double > gamma
Parameters for structural refinement.
unsigned int cumkappa
Cumulative kappa refinements.
virtual void AddOptions(boost::property_tree::ptree &pt) const override
virtual std::shared_ptr< SamplingState > QOI() override
double lastLyapunov
The Lyapunov function at the most recently accepted point.
std::shared_ptr< muq::Approximation::LocalRegression > regQoI
The regression for the quantity of interest.
std::shared_ptr< muq::Approximation::LocalRegression > reg
The regression for the log-likelihood.
unsigned int cumbeta
Cumulative beta refinements.
double lastThreshold
The error threshold at the most recently accepted point.
const Eigen::VectorXd centroid
The centroid of the distribution (input parameter)
virtual double LogDensity(std::shared_ptr< SamplingState > const &state) override
ExpensiveSamplingProblem(std::shared_ptr< muq::Modeling::ModPiece > const &target, Eigen::VectorXd const &centroid, boost::property_tree::ptree pt)
double RefineSurrogate(unsigned int const step, std::shared_ptr< SamplingState > const &state, std::vector< Eigen::VectorXd > &neighbors, std::vector< Eigen::VectorXd > &results)
void CheckNumNeighbors(std::shared_ptr< SamplingState > const &state)
Check to make sure we have enough model evaluations.
double LogLyapunovFunction(Eigen::VectorXd const &x) const
Estimate the Lyapunov function.
ExpensiveSamplingProblem(std::shared_ptr< muq::Modeling::ModPiece > const &target, std::shared_ptr< muq::Modeling::ModPiece > const &qoi, Eigen::VectorXd const &centroid, boost::property_tree::ptree pt)
unsigned int level
The current error threshold level.
double ErrorThreshold(Eigen::VectorXd const &x) const
Compute the error threshold.
void SetUp(boost::property_tree::ptree &pt)
Set up the sampling problem.
void CheckNeighbors(std::shared_ptr< SamplingState > const &state, std::vector< Eigen::VectorXd > &neighbors, std::vector< Eigen::VectorXd > &results) const
Check to make sure we have enough model evaluations.
unsigned int cumgamma
Cumulative gamma refinements.
ExpensiveSamplingProblem(std::shared_ptr< muq::Modeling::ModPiece > const &target, boost::property_tree::ptree pt)
std::pair< double, double > lyapunovPara
The scaling and exponent for the Lyapunov function.
Class for sampling problems based purely on a density function.
std::shared_ptr< muq::Modeling::ModPiece > qoi
std::shared_ptr< muq::Modeling::ModPiece > target
The target distribution (the prior in the inference case)