MUQ  0.4.3
LocalRegression.h
Go to the documentation of this file.
1 #ifndef LOCALREGRESSION_H_
2 #define LOCALREGRESSION_H_
3 
4 #include "MUQ/config.h"
5 
6 #if MUQ_HAS_PARCER
7 #include <parcer/Communicator.h>
8 #endif
9 
10 #include "MUQ/Modeling/ModPiece.h"
12 
14 
15 namespace muq {
16  namespace Approximation {
18  public:
19 
24  LocalRegression(std::shared_ptr<muq::Modeling::ModPiece> function, boost::property_tree::ptree& pt);
25 
26 #if MUQ_HAS_PARCER
32  LocalRegression(std::shared_ptr<muq::Modeling::ModPiece> function, boost::property_tree::ptree& pt, std::shared_ptr<parcer::Communicator> comm);
33 #endif
34 
35 #if MUQ_HAS_PARCER
37 #else
38  ~LocalRegression() = default;
39 #endif
40 
42 
45  void Add(std::vector<Eigen::VectorXd> const& inputs) const;
46 
48 
52  Eigen::VectorXd Add(Eigen::VectorXd const& input) const;
53 
55 
58  unsigned int CacheSize() const;
59 
61 
65  bool InCache(Eigen::VectorXd const& point) const;
66 
68 
72  Eigen::VectorXd CachePoint(unsigned int const index) const;
73 
75 
80  std::tuple<Eigen::VectorXd, double, unsigned int> PoisednessConstant(Eigen::VectorXd const& input) const;
81 
83 
89  std::tuple<Eigen::VectorXd, double, unsigned int> PoisednessConstant(Eigen::VectorXd const& input, std::vector<Eigen::VectorXd> const& neighbors) const;
90 
92 
97  std::pair<double, double> ErrorIndicator(Eigen::VectorXd const& input) const;
98 
100 
106  std::pair<double, double> ErrorIndicator(Eigen::VectorXd const& input, std::vector<Eigen::VectorXd> const& neighbors) const;
107 
109 
113  void NearestNeighbors(Eigen::VectorXd const& input, std::vector<Eigen::VectorXd>& neighbors) const;
114 
116 
121  void NearestNeighbors(Eigen::VectorXd const& input, std::vector<Eigen::VectorXd>& neighbors, std::vector<Eigen::VectorXd>& result) const;
122 
123  Eigen::VectorXd EvaluateRegressor(Eigen::VectorXd const& input, std::vector<Eigen::VectorXd> const& neighbors, std::vector<Eigen::VectorXd> const& result) const;
124 
125 #if MUQ_HAS_PARCER
126  void Probe() const;
127 #endif
128 
130  const unsigned int kn;
131 
133 
136  Eigen::VectorXd CacheCentroid() const;
137 
139  void ClearCache();
140 
142  unsigned int Order() const;
143 
144  private:
145 
146  virtual void EvaluateImpl(muq::Modeling::ref_vector<Eigen::VectorXd> const& inputs) override;
147 
149  void FitRegression(Eigen::VectorXd const& input) const;
150 
152  void SetUp(std::shared_ptr<muq::Modeling::ModPiece> function, boost::property_tree::ptree& pt);
153 
155  std::shared_ptr<muq::Modeling::FlannCache> cache;
156 
158  std::shared_ptr<Regression> reg;
159 
160 #if MUQ_HAS_PARCER
161  std::shared_ptr<parcer::Communicator> comm = nullptr;
162  const int tagSingle = 0;
163 #endif
164  };
165  } // namespace Approximation
166 } // namespace muq
167 
168 #endif
std::pair< double, double > ErrorIndicator(Eigen::VectorXd const &input) const
Get the error indicator.
Eigen::VectorXd CacheCentroid() const
Get the centroid of the cache.
virtual void EvaluateImpl(muq::Modeling::ref_vector< Eigen::VectorXd > const &inputs) override
void SetUp(std::shared_ptr< muq::Modeling::ModPiece > function, boost::property_tree::ptree &pt)
Set up the regressor.
unsigned int Order() const
Polynomial order.
Eigen::VectorXd EvaluateRegressor(Eigen::VectorXd const &input, std::vector< Eigen::VectorXd > const &neighbors, std::vector< Eigen::VectorXd > const &result) const
Eigen::VectorXd CachePoint(unsigned int const index) const
A point in the cache.
unsigned int CacheSize() const
Get the total size of the cache.
std::shared_ptr< muq::Modeling::FlannCache > cache
A cache containing previous model evaluations.
std::tuple< Eigen::VectorXd, double, unsigned int > PoisednessConstant(Eigen::VectorXd const &input) const
Get the poisedness constant.
std::shared_ptr< parcer::Communicator > comm
void Add(std::vector< Eigen::VectorXd > const &inputs) const
Add some points to the cache.
std::shared_ptr< Regression > reg
A regressor.
LocalRegression(std::shared_ptr< muq::Modeling::ModPiece > function, boost::property_tree::ptree &pt)
const unsigned int kn
The number of nearest neighbors to use by the regressor.
bool InCache(Eigen::VectorXd const &point) const
Is a point in the cache?
void FitRegression(Eigen::VectorXd const &input) const
Fit the regression to the nearest neighbors.
void NearestNeighbors(Eigen::VectorXd const &input, std::vector< Eigen::VectorXd > &neighbors) const
Get the number of nearest neighbors.
Provides an abstract interface for defining vector-valued model components.
Definition: ModPiece.h:148
std::vector< std::reference_wrapper< const T > > ref_vector
A vector of references to something ...
Definition: WorkPiece.h:37