1 #ifndef LOCALREGRESSION_H_
2 #define LOCALREGRESSION_H_
7 #include <parcer/Communicator.h>
16 namespace Approximation {
24 LocalRegression(std::shared_ptr<muq::Modeling::ModPiece>
function, boost::property_tree::ptree& pt);
32 LocalRegression(std::shared_ptr<muq::Modeling::ModPiece>
function, boost::property_tree::ptree& pt, std::shared_ptr<parcer::Communicator>
comm);
45 void Add(std::vector<Eigen::VectorXd>
const& inputs)
const;
52 Eigen::VectorXd
Add(Eigen::VectorXd
const& input)
const;
65 bool InCache(Eigen::VectorXd
const& point)
const;
72 Eigen::VectorXd
CachePoint(
unsigned int const index)
const;
80 std::tuple<Eigen::VectorXd, double, unsigned int>
PoisednessConstant(Eigen::VectorXd
const& input)
const;
89 std::tuple<Eigen::VectorXd, double, unsigned int>
PoisednessConstant(Eigen::VectorXd
const& input, std::vector<Eigen::VectorXd>
const& neighbors)
const;
97 std::pair<double, double>
ErrorIndicator(Eigen::VectorXd
const& input)
const;
106 std::pair<double, double>
ErrorIndicator(Eigen::VectorXd
const& input, std::vector<Eigen::VectorXd>
const& neighbors)
const;
113 void NearestNeighbors(Eigen::VectorXd
const& input, std::vector<Eigen::VectorXd>& neighbors)
const;
121 void NearestNeighbors(Eigen::VectorXd
const& input, std::vector<Eigen::VectorXd>& neighbors, std::vector<Eigen::VectorXd>& result)
const;
123 Eigen::VectorXd
EvaluateRegressor(Eigen::VectorXd
const& input, std::vector<Eigen::VectorXd>
const& neighbors, std::vector<Eigen::VectorXd>
const& result)
const;
130 const unsigned int kn;
142 unsigned int Order()
const;
152 void SetUp(std::shared_ptr<muq::Modeling::ModPiece>
function, boost::property_tree::ptree& pt);
155 std::shared_ptr<muq::Modeling::FlannCache>
cache;
158 std::shared_ptr<Regression>
reg;
161 std::shared_ptr<parcer::Communicator>
comm =
nullptr;
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
void ClearCache()
Clear the cache.
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)
~LocalRegression()=default
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.
std::vector< std::reference_wrapper< const T > > ref_vector
A vector of references to something ...