25 #ifndef LSST_AFW_MATH_GAUSSIAN_PROCESS_H 26 #define LSST_AFW_MATH_GAUSSIAN_PROCESS_H 28 #include <Eigen/Dense> 31 #include "ndarray/eigen.h" 114 double _before, _beginning;
115 double _eigenTime, _iterationTime, _searchTime, _varianceTime, _totalTime;
116 int _interpolationCount;
125 template <
typename T>
150 virtual T operator()(ndarray::Array<const T, 1, 1>
const &p1,
151 ndarray::Array<const T, 1, 1>
const &p2)
const;
163 template <
typename T>
174 void setEllSquared(
double ellSquared);
176 T operator()(ndarray::Array<const T, 1, 1>
const &, ndarray::Array<const T, 1, 1>
const &)
const override;
192 template <
typename T>
202 void setSigma0(
double sigma0);
207 void setSigma1(
double sigma1);
209 T operator()(ndarray::Array<const T, 1, 1>
const &, ndarray::Array<const T, 1, 1>
const &)
const override;
212 double _sigma0, _sigma1;
223 template <
typename T>
245 void Initialize(ndarray::Array<T, 2, 2>
const &dt);
262 void findNeighbors(ndarray::Array<int, 1, 1> neighdex, ndarray::Array<double, 1, 1> dd,
263 ndarray::Array<const T, 1, 1>
const &v,
int n_nn)
const;
272 T getData(
int ipt,
int idim)
const;
295 ndarray::Array<T, 1, 1> getData(
int ipt)
const;
304 void addPoint(ndarray::Array<const T, 1, 1>
const &v);
314 void removePoint(
int dex);
319 int getNPoints()
const;
328 void getTreeNode(ndarray::Array<int, 1, 1>
const &v,
int dex)
const;
331 ndarray::Array<int, 2, 2> _tree;
332 ndarray::Array<int, 1, 1> _inn;
333 ndarray::Array<T, 2, 2> _data;
335 enum { DIMENSION, LT, GEQ,
PARENT };
350 int _npts, _dimensions, _room, _roomStep, _masterParent;
351 mutable int _neighborsFound, _neighborsWanted;
357 mutable ndarray::Array<double, 1, 1> _neighborDistances;
358 mutable ndarray::Array<int, 1, 1> _neighborCandidates;
372 void _organize(ndarray::Array<int, 1, 1>
const &use,
int ct,
int parent,
int dir);
379 int _findNode(ndarray::Array<const T, 1, 1>
const &v)
const;
397 void _lookForNeighbors(ndarray::Array<const T, 1, 1>
const &v,
int consider,
int from)
const;
402 int _testTree()
const;
417 int _walkUpTree(
int target,
int dir,
int root)
const;
427 void _count(
int where,
int *ct)
const;
434 void _descend(
int root);
441 void _reassign(
int target);
446 double _distance(ndarray::Array<const T, 1, 1>
const &p1, ndarray::Array<const T, 1, 1>
const &p2)
const;
470 template <
typename T>
493 GaussianProcess(ndarray::Array<T, 2, 2>
const &dataIn, ndarray::Array<T, 1, 1>
const &ff,
515 GaussianProcess(ndarray::Array<T, 2, 2>
const &dataIn, ndarray::Array<T, 1, 1>
const &mn,
516 ndarray::Array<T, 1, 1>
const &mx, ndarray::Array<T, 1, 1>
const &ff,
530 GaussianProcess(ndarray::Array<T, 2, 2>
const &dataIn, ndarray::Array<T, 2, 2>
const &ff,
548 GaussianProcess(ndarray::Array<T, 2, 2>
const &dataIn, ndarray::Array<T, 1, 1>
const &mn,
549 ndarray::Array<T, 1, 1>
const &mx, ndarray::Array<T, 2, 2>
const &ff,
555 int getNPoints()
const;
572 void getData(ndarray::Array<T, 2, 2> pts, ndarray::Array<T, 1, 1> fn,
573 ndarray::Array<int, 1, 1> indices)
const;
584 void getData(ndarray::Array<T, 2, 2> pts, ndarray::Array<T, 2, 2> fn,
585 ndarray::Array<int, 1, 1> indices)
const;
605 int numberOfNeighbors)
const;
621 void interpolate(ndarray::Array<T, 1, 1> mu, ndarray::Array<T, 1, 1> variance,
622 ndarray::Array<T, 1, 1>
const &vin,
int numberOfNeighbors)
const;
643 T selfInterpolate(ndarray::Array<T, 1, 1> variance,
int dex,
int numberOfNeighbors)
const;
659 void selfInterpolate(ndarray::Array<T, 1, 1> mu, ndarray::Array<T, 1, 1> variance,
int dex,
660 int numberOfNeighbors)
const;
683 void batchInterpolate(ndarray::Array<T, 1, 1> mu, ndarray::Array<T, 1, 1> variance,
684 ndarray::Array<T, 2, 2>
const &queries)
const;
703 void batchInterpolate(ndarray::Array<T, 1, 1> mu, ndarray::Array<T, 2, 2>
const &queries)
const;
709 void batchInterpolate(ndarray::Array<T, 2, 2> mu, ndarray::Array<T, 2, 2> variance,
710 ndarray::Array<T, 2, 2>
const &queries)
const;
716 void batchInterpolate(ndarray::Array<T, 2, 2> mu, ndarray::Array<T, 2, 2>
const &queries)
const;
734 void addPoint(ndarray::Array<T, 1, 1>
const &vin,
T f);
745 void addPoint(ndarray::Array<T, 1, 1>
const &vin, ndarray::Array<T, 1, 1>
const &f);
758 void removePoint(
int dex);
766 void setKrigingParameter(
T kk);
787 void setLambda(
T lambda);
806 int _npts, _useMaxMin, _dimensions, _room, _roomStep, _nFunctions;
808 T _krigingParameter, _lambda;
810 ndarray::Array<T, 1, 1> _max, _min;
811 ndarray::Array<T, 2, 2> _function;
822 #endif //#ifndef LSST_AFW_MATH_GAUSSIAN_PROCESS_H GaussianProcessTimer & operator=(GaussianProcessTimer const &)=default
Key< Flag > const & target
Stores values of a function sampled on an image and allows you to interpolate the function to unsampl...
void addToTotal(int i)
Adds time to _totalTime and adds counts to _interpolationCount.
This is a structure for keeping track of how long the interpolation methods spend on different parts ...
The parent class of covariogram functions for use in Gaussian Process interpolation.
void start()
Starts the timer for an individual call to an interpolation routine.
void addToVariance()
Adds time to _varianceTime.
~GaussianProcessTimer()=default
A base class for image defects.
void addToIteration()
Adds time to _iterationTime.
The data for GaussianProcess is stored in a KD tree to facilitate nearest-neighbor searches...
a Covariogram that recreates a neural network with one hidden layer and infinite units in that layer ...
Interface for DateTime class.
Covariogram()
construct a Covariogram assigning default values to the hyper parameters
void addToSearch()
Adds time to _searchTime.
void display()
Displays the current values of all times and _interpolationCount.
void reset()
Resets all of the data members of GaussianProcessTimer to zero.
table::Key< double > sigma1
void addToEigen()
Adds time to _eigenTime.