kde_evaluate.cpp
// (C) Copyright Renaud Detry 2007-2015.
// Distributed under the GNU General Public License and under the
// BSD 3-Clause License (See accompanying file LICENSE.txt).
// This example demonstrates how Nulkei can be used to estimate a density
// function at a given point.
//
// This example can be compiled with
//
// g++ `pkg-config --cflags --libs nuklei` -O3 kde_evaluate.cpp -o kde_evaluate
//
// (It may be necessary to add -I/path/to/boost/include -L/path/to/boost/lib if
// Boost doesn't reside at a standard place.)
int main(int argc, char ** argv)
{
try {
// ----------- //
// Parameters: //
// ----------- //
// Set of datapoints that represent a density
std::string densityFilename = "data/points1.txt";
// Set of points at which the density will be evaluated
std::string pointsFilename = "data/points2.txt";
// Kernel widths, for position and orientation:
double locH = 40; // in the same unit as the datapoints forming the density
double oriH = .4; // in radians
// ------------- //
// Read-in data: //
// ------------- //
nuklei::KernelCollection density, points;
nuklei::readObservations(densityFilename, density);
nuklei::readObservations(pointsFilename, points);
// ------------------------------- //
// Prepare density for evaluation: //
// ------------------------------- //
density.setKernelLocH(locH);
density.setKernelOriH(oriH);
density.normalizeWeights();
density.buildKdTree();
// At this point, "density" should not be modified anymore. (Modifying it will
// destroy the kd-tree, and the kernel statistics.)
// --------------------- //
// Evaluate the density: //
// --------------------- //
i != points.end(); ++i)
{
std::cout << density.evaluationAt(*i) << std::endl;
}
return 0;
}
catch (std::exception &e) {
std::cerr << "Exception caught: ";
std::cerr << e.what() << std::endl;
return EXIT_FAILURE;
}
catch (...) {
std::cerr << "Caught unknown exception." << std::endl;
return EXIT_FAILURE;
}
}
Container::const_iterator const_iterator
KernelCollection iterator.
This class acts as a vector-like container for kernels. It also provides methods related to kernel de...
void setKernelLocH(coord_t h)
Sets the location bandwidth of all kernels.
void setKernelOriH(coord_t h)
Sets the orientation bandwidth of all kernels.
Container::iterator end()
Returns an iterator pointing to the last kernel.
Container::iterator begin()
Returns an iterator pointing to the first kernel.
void buildKdTree()
Builds a kd-tree of the kernel positions and stores the tree internally. See Intermediary Results.
weight_t evaluationAt(const kernel::base &f, const EvaluationStrategy strategy=WEIGHTED_SUM_EVAL) const
Evaluates the density represented by *this at f.
void normalizeWeights()
Divides all weights by the total weight.
void readObservations(ObservationReader &r, KernelCollection &kc)
Reads the data available from the reader r and stores it into kc.
© Copyright 2007-2013 Renaud Detry.
Distributed under the terms of the GNU General Public License (GPL).
(See accompanying file LICENSE.txt or copy at http://www.gnu.org/copyleft/gpl.html.)
Revised Sun Sep 13 2020 19:10:06.