8 #include <nanoflann.hpp>
19 template <
class Distance = nanoflann::metric_L2,
typename IndexType =
size_t>
25 typedef typename Distance::template traits<double,self_t>::distance_t
metric_t;
33 index = std::make_shared<index_t>(dim, *
this , nanoflann::KDTreeSingleIndexAdaptorParams(leaf_max_size ) );
38 index = std::make_shared<index_t>(
m_data.at(0).size(), *
this , nanoflann::KDTreeSingleIndexAdaptorParams(leaf_max_size ) );
44 inline void add(Eigen::VectorXd
const& newPt) {
54 inline std::pair<std::vector<IndexType>, std::vector<double>>
query(Eigen::VectorXd
const& query_point,
const size_t num_closest,
const int nChecks_IGNORED = 10)
const {
55 std::vector<IndexType> out_indices(num_closest);
56 std::vector<double> out_distances_sq(num_closest);
58 nanoflann::KNNResultSet<double,IndexType> resultSet(num_closest);
59 resultSet.init(&out_indices[0], &out_distances_sq[0]);
60 #if MUQ_NANOFLAN_PARAMS_COMPILES
61 index->findNeighbors(resultSet, query_point.data(), nanoflann::SearchParams());
63 index->findNeighbors(resultSet, query_point.data(), nanoflann::SearchParameters());
65 return std::make_pair(out_indices, out_distances_sq);
84 assert(dim<
m_data[idx].size());
110 FlannCache(std::shared_ptr<ModPiece>
function);
119 int InCache(Eigen::VectorXd
const& input)
const;
126 std::vector<Eigen::VectorXd>
Add(std::vector<Eigen::VectorXd>
const& inputs);
133 void Add(std::vector<Eigen::VectorXd>
const& inputs, std::vector<Eigen::VectorXd>
const& results);
140 const Eigen::VectorXd
at(
unsigned int const index)
const;
147 Eigen::VectorXd
at(
unsigned int const index);
150 Eigen::VectorXd
const&
OutputValue(
unsigned int index)
const;
157 Eigen::VectorXd
Add(Eigen::VectorXd
const& input);
165 unsigned int Add(Eigen::VectorXd
const& input, Eigen::VectorXd
const& result);
171 void Remove(Eigen::VectorXd
const& input);
188 unsigned int const k,
189 std::vector<Eigen::VectorXd>& neighbors,
190 std::vector<Eigen::VectorXd>& result)
const;
199 unsigned int const k,
200 std::vector<Eigen::VectorXd>& neighbors)
const;
206 unsigned int Size()
const;
218 std::shared_ptr<ModPiece>
Function()
const;
234 std::shared_ptr<ModPiece>
function;
237 std::shared_ptr<DynamicKDTreeAdaptor<>>
kdTree;
Create a cache of model evaluations (input/output pairs)
std::vector< Eigen::VectorXd > outputCache
virtual void EvaluateImpl(ref_vector< Eigen::VectorXd > const &inputs) override
void UpdateCentroid(Eigen::VectorXd const &point)
Update the centroid when a new point is added to the cache.
std::vector< Eigen::VectorXd > Add(std::vector< Eigen::VectorXd > const &inputs)
Add new points to the cache.
std::shared_ptr< DynamicKDTreeAdaptor<> > kdTree
The nearest neighbor index, used to perform searches.
void NearestNeighbors(Eigen::VectorXd const &point, unsigned int const k, std::vector< Eigen::VectorXd > &neighbors) const
Find the nearest neighbors.
FlannCache(std::shared_ptr< ModPiece > function)
Eigen::VectorXd Centroid() const
Get the centroid of the cache.
unsigned int Size() const
Get the size of the cache.
const Eigen::VectorXd at(unsigned int const index) const
Get an input point from the cache.
Eigen::VectorXd centroid
The centroid of the input addPoints.
size_t NearestNeighborIndex(Eigen::VectorXd const &point) const
The index of the nearest neighbor.
void Remove(Eigen::VectorXd const &input)
Remove point from the cache.
int InCache(Eigen::VectorXd const &input) const
Determine if an entry is in the cache.
std::shared_ptr< ModPiece > Function() const
Get the underlying function.
Eigen::VectorXd const & OutputValue(unsigned int index) const
Returns the model for a specific cache index.
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 ...
bool kdtree_get_bbox(BBOX &) const
size_t kdtree_get_point_count() const
void UpdateIndex(const int leaf_max_size=10)
nanoflann::KDTreeSingleIndexDynamicAdaptor< metric_t, self_t,-1, IndexType > index_t
Distance::template traits< double, self_t >::distance_t metric_t
DynamicKDTreeAdaptor(const int dim, const int leaf_max_size=10)
Constructor: takes a const ref to the vector of vectors object with the data points.
double kdtree_get_pt(const size_t idx, int dim) const
std::shared_ptr< index_t > index
std::deque< Eigen::VectorXd > m_data
The kd-tree index for the user to call its methods as usual with any other FLANN index.
DynamicKDTreeAdaptor< Distance, IndexType > self_t
const self_t & derived() const
std::pair< std::vector< IndexType >, std::vector< double > > query(Eigen::VectorXd const &query_point, const size_t num_closest, const int nChecks_IGNORED=10) const
virtual ~DynamicKDTreeAdaptor()
void add(Eigen::VectorXd const &newPt)