LSSTApplications  11.0-13-gbb96280,12.1.rc1,12.1.rc1+1,12.1.rc1+2,12.1.rc1+5,12.1.rc1+8,12.1.rc1-1-g06d7636+1,12.1.rc1-1-g253890b+5,12.1.rc1-1-g3d31b68+7,12.1.rc1-1-g3db6b75+1,12.1.rc1-1-g5c1385a+3,12.1.rc1-1-g83b2247,12.1.rc1-1-g90cb4cf+6,12.1.rc1-1-g91da24b+3,12.1.rc1-2-g3521f8a,12.1.rc1-2-g39433dd+4,12.1.rc1-2-g486411b+2,12.1.rc1-2-g4c2be76,12.1.rc1-2-gc9c0491,12.1.rc1-2-gda2cd4f+6,12.1.rc1-3-g3391c73+2,12.1.rc1-3-g8c1bd6c+1,12.1.rc1-3-gcf4b6cb+2,12.1.rc1-4-g057223e+1,12.1.rc1-4-g19ed13b+2,12.1.rc1-4-g30492a7
LSSTDataManagementBasePackage
GaussianProcess.h
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 
3 /*
4  * LSST Data Management System
5  * Copyright 2008, 2009, 2010 LSST Corporation.
6  *
7  * This product includes software developed by the
8  * LSST Project (http://www.lsst.org/).
9  *
10  * This program is free software: you can redistribute it and/or modify
11  * it under the terms of the GNU General Public License as published by
12  * the Free Software Foundation, either version 3 of the License, or
13  * (at your option) any later version.
14  *
15  * This program is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the LSST License Statement and
21  * the GNU General Public License along with this program. If not,
22  * see <http://www.lsstcorp.org/LegalNotices/>.
23  */
24 
25 
26 #ifndef LSST_AFW_MATH_GAUSSIAN_PROCESS_H
27 #define LSST_AFW_MATH_GAUSSIAN_PROCESS_H
28 
29 #include <Eigen/Dense>
30 #include <stdexcept>
31 
32 #include "ndarray/eigen.h"
33 #include <memory>
34 
35 #include "lsst/daf/base/Citizen.h"
36 #include "lsst/daf/base/DateTime.h"
37 #include "lsst/pex/exceptions.h"
38 #include "lsst/pex/policy.h"
39 #include "lsst/pex/policy/Policy.h"
40 
41 namespace lsst {
42 namespace afw {
43 namespace math {
44 
67 
68 public:
69 
71 
76  void reset();
77 
81  void start();
82 
86  void addToEigen();
87 
91  void addToVariance();
92 
96  void addToSearch();
97 
101  void addToIteration();
102 
108  void addToTotal(int i);
109 
113  void display();
114 
115 private:
119 };
120 
121 
130 template <typename T>
132 public:
133  virtual ~Covariogram();
134 
135 #ifndef SWIG
136  // No copying
137  Covariogram (const Covariogram&) = delete;
138  Covariogram& operator=(const Covariogram&) = delete;
139 
140  // No moving
141  Covariogram (Covariogram&&) = delete;
142  Covariogram& operator=(Covariogram&&) = delete;
143 #endif
144 
148  explicit Covariogram():lsst::daf::base::Citizen(typeid(this)){};
149 
150 
158  virtual T operator() (ndarray::Array<const T,1,1> const &p1,
159  ndarray::Array<const T,1,1> const &p2
160  ) const;
161 
162 };
163 
173 template <typename T>
175 
176 public:
177  virtual ~SquaredExpCovariogram();
178 
179  explicit SquaredExpCovariogram();
180 
185  void setEllSquared(double ellSquared);
186 
187  virtual T operator() (ndarray::Array<const T,1,1> const &,
188  ndarray::Array<const T,1,1> const &
189  ) const;
190 
191 private:
192  double _ellSquared;
193 
194 };
195 
208 template <typename T>
210 
211 public:
212  virtual ~NeuralNetCovariogram();
213 
214  explicit NeuralNetCovariogram();
215 
219  void setSigma0(double sigma0);
220 
224  void setSigma1(double sigma1);
225 
226  virtual T operator() (ndarray::Array<const T,1,1> const &,
227  ndarray::Array<const T,1,1> const &
228  ) const;
229 
230 private:
231  double _sigma0,_sigma1;
232 
233 
234 };
235 
246 template <typename T>
247 class KdTree {
248 public:
249 
250 #ifndef SWIG
251  // No copying
252  KdTree (const KdTree&) = delete;
253  KdTree& operator=(const KdTree&) = delete;
254 
255  // No moving
256  KdTree (KdTree&&) = delete;
257  KdTree& operator=(KdTree&&) = delete;
258 
259  // Add default constructor (which is no longer generated)
260  KdTree() = default;
261 #endif
262 
271  void Initialize(ndarray::Array<T,2,2> const &dt);
272 
288  void findNeighbors(ndarray::Array<int,1,1> neighdex,
289  ndarray::Array<double,1,1> dd,
290  ndarray::Array<const T,1,1> const &v,
291  int n_nn) const;
292 
293 
301  T getData(int ipt, int idim) const;
302 
303 
325  ndarray::Array<T,1,1> getData(int ipt) const;
326 
334  void addPoint(ndarray::Array<const T,1,1> const &v);
335 
344  void removePoint(int dex);
345 
349  int getPoints() const;
350 
358  void getTreeNode(ndarray::Array<int,1,1> const &v,int dex) const;
359 
360 
361 private:
362  ndarray::Array<int,2,2> _tree;
363  ndarray::Array<int,1,1> _inn;
364  ndarray::Array<T,2,2> _data;
365 
367 
368  //_tree stores the relationships between data points
369  //_tree[i][DIMENSION] is the dimension on which the ith node segregates its daughters
370  //
371  //_tree[i][LT] is the branch of the tree down which the daughters' DIMENSIONth component
372  //is less than the parent's
373  //
374  //_tree[i][GEQ] is the branch of the tree down which the daughters' DIMENSIONth component is
375  //greather than or equal to the parent's
376  //
377  //_tree[i][PARENT] is the parent node of the ith node
378 
379  //_data actually stores the data points
380 
381 
384 
385  //_room denotes the capacity of _data and _tree. It will usually be larger
386  //than _pts so that we do not have to reallocate
387  //_tree and _data every time we add a new point to the tree
388 
389  mutable ndarray::Array<double,1,1> _neighborDistances;
390  mutable ndarray::Array<int,1,1> _neighborCandidates;
391 
404  void _organize(ndarray::Array<int,1,1> const &use,
405  int ct,
406  int parent,
407  int dir);
408 
414  int _findNode(ndarray::Array<const T,1,1> const &v) const;
415 
432  void _lookForNeighbors(ndarray::Array<const T,1,1> const &v,
433  int consider,
434  int from) const;
435 
439  int _testTree() const;
440 
454  int _walkUpTree(int target,
455  int dir,
456  int root) const;
457 
466  void _count(int where,int *ct) const;
467 
473  void _descend(int root);
474 
480  void _reassign(int target);
481 
485  double _distance(ndarray::Array<const T,1,1> const &p1,
486  ndarray::Array<const T,1,1> const &p2) const;
487 };
488 
489 
513 template <typename T>
515 
516 public:
517 
518 #ifndef SWIG
519  // No copying
520  GaussianProcess (const GaussianProcess&) = delete;
521  GaussianProcess& operator=(const GaussianProcess&) = delete;
522 
523  // No moving
524  GaussianProcess (GaussianProcess&&) = delete;
526 #endif
527 
540  GaussianProcess(ndarray::Array<T,2,2> const &dataIn,
541  ndarray::Array<T,1,1> const &ff,
542  std::shared_ptr< Covariogram<T> > const &covarIn);
543 
563  GaussianProcess(ndarray::Array<T,2,2> const &dataIn,
564  ndarray::Array<T,1,1> const &mn,
565  ndarray::Array<T,1,1> const &mx,
566  ndarray::Array<T,1,1> const &ff,
567  std::shared_ptr< Covariogram<T> > const &covarIn);
580  GaussianProcess(ndarray::Array<T,2,2> const &dataIn,
581  ndarray::Array<T,2,2> const &ff,
582  std::shared_ptr< Covariogram<T> > const &covarIn);
599  GaussianProcess(ndarray::Array<T,2,2> const &dataIn,
600  ndarray::Array<T,1,1> const &mn,
601  ndarray::Array<T,1,1> const &mx,
602  ndarray::Array<T,2,2> const &ff,
603  std::shared_ptr< Covariogram<T> > const &covarIn);
604 
622  T interpolate(ndarray::Array<T,1,1> variance,
623  ndarray::Array<T,1,1> const &vin,
624  int numberOfNeighbors) const;
640  void interpolate(ndarray::Array<T,1,1> mu,
641  ndarray::Array<T,1,1> variance,
642  ndarray::Array<T,1,1> const &vin,
643  int numberOfNeighbors) const;
644 
664  T selfInterpolate(ndarray::Array<T,1,1> variance,
665  int dex,
666  int numberOfNeighbors) const;
667 
682  void selfInterpolate(ndarray::Array<T,1,1> mu,
683  ndarray::Array<T,1,1> variance,
684  int dex,
685  int numberOfNeighbors) const;
686 
708  void batchInterpolate(ndarray::Array<T,1,1> mu,
709  ndarray::Array<T,1,1> variance,
710  ndarray::Array<T,2,2> const &queries) const;
711 
729  void batchInterpolate(ndarray::Array<T,1,1> mu,
730  ndarray::Array<T,2,2> const &queries) const;
731 
736  void batchInterpolate(ndarray::Array<T,2,2> mu,
737  ndarray::Array<T,2,2> variance,
738  ndarray::Array<T,2,2> const &queries) const;
739 
744  void batchInterpolate(ndarray::Array<T,2,2> mu,
745  ndarray::Array<T,2,2> const &queries) const;
746 
747 
764  void addPoint(ndarray::Array<T,1,1> const &vin,T f);
765 
775  void addPoint(ndarray::Array<T,1,1> const &vin,ndarray::Array<T,1,1> const &f);
776 
788  void removePoint(int dex);
789 
796  void setKrigingParameter(T kk);
797 
804  void setCovariogram(std::shared_ptr< Covariogram<T> > const &covar);
805 
817  void setLambda(T lambda);
818 
819 
835 
836  private:
838 
840 
841  ndarray::Array<T,1,1> _max,_min;
842  ndarray::Array<T,2,2> _function;
843 
845 
846  std::shared_ptr< Covariogram<T> > _covariogram;
848 
849 };
850 
851 }}}
852 
853 
854 
855 #endif //#ifndef LSST_AFW_MATH_GAUSSIAN_PROCESS_H
void setCovariogram(std::shared_ptr< Covariogram< T > > const &covar)
Assign a different covariogram to this GaussianProcess.
GaussianProcess(const GaussianProcess &)=delete
ndarray::Array< double, 1, 1 > _neighborDistances
int _findNode(ndarray::Array< const T, 1, 1 > const &v) const
Find the point already in the tree that would be the parent of a point not in the tree...
void setSigma1(double sigma1)
set the _sigma1 hyper parameter
void addPoint(ndarray::Array< const T, 1, 1 > const &v)
Add a point to the tree. Allot more space in _tree and data if needed.
ndarray::Array< int, 1, 1 > _neighborCandidates
T getData(int ipt, int idim) const
Return one element of one node on the tree.
void setSigma0(double sigma0)
set the _sigma0 hyper parameter
void _organize(ndarray::Array< int, 1, 1 > const &use, int ct, int parent, int dir)
Find the daughter point of a node in the tree and segregate the points around it. ...
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 ...
ndarray::Array< T, 1, 1 > _min
void _lookForNeighbors(ndarray::Array< const T, 1, 1 > const &v, int consider, int from) const
This method actually looks for the neighbors, determining whether or not to descend branches of the t...
void findNeighbors(ndarray::Array< int, 1, 1 > neighdex, ndarray::Array< double, 1, 1 > dd, ndarray::Array< const T, 1, 1 > const &v, int n_nn) const
Find the nearest neighbors of a point.
The parent class of covariogram functions for use in Gaussian Process interpolation.
ndarray::Array< int, 2, 2 > _tree
void setKrigingParameter(T kk)
Assign a value to the Kriging paramter.
double _distance(ndarray::Array< const T, 1, 1 > const &p1, ndarray::Array< const T, 1, 1 > const &p2) const
calculate the Euclidean distance between the points p1 and p2
void start()
Starts the timer for an individual call to an interpolation routine.
void _count(int where, int *ct) const
A method which counts the number of nodes descended from a given node (used by removePoint(int)) ...
void addToVariance()
Adds time to _varianceTime.
void setLambda(T lambda)
set the value of the hyperparameter _lambda
ndarray::Array< T, 2, 2 > _function
Covariogram & operator=(const Covariogram &)=delete
void Initialize(ndarray::Array< T, 2, 2 > const &dt)
Build a KD Tree to store the data for GaussianProcess.
T selfInterpolate(ndarray::Array< T, 1, 1 > variance, int dex, int numberOfNeighbors) const
This method will interpolate the function on a data point for purposes of optimizing hyper parameters...
void batchInterpolate(ndarray::Array< T, 1, 1 > mu, ndarray::Array< T, 1, 1 > variance, ndarray::Array< T, 2, 2 > const &queries) const
Interpolate a list of query points using all of the input data (rather than nearest neighbors) ...
ndarray::Array< T, 1, 1 > _max
void addToIteration()
Adds time to _iterationTime.
std::shared_ptr< Covariogram< T > > _covariogram
The data for GaussianProcess is stored in a KD tree to facilitate nearest-neighbor searches...
a Covariogram that falls off as the negative exponent of the square of the distance between the point...
void _descend(int root)
Descend the tree from a node which has been removed, reassigning severed nodes as you go...
a Covariogram that recreates a neural network with one hidden layer and infinite units in that layer ...
GaussianProcessTimer & getTimes() const
Give the user acces to _timer, an object keeping track of the time spent on various processes within ...
Interface for DateTime class.
Covariogram()
construct a Covariogram assigning default values to the hyper parameters
Citizen(const std::type_info &)
Definition: Citizen.cc:174
afw::table::Key< double > sigma1
void removePoint(int dex)
Remove a point from the tree. Reorganize what remains so that the tree remains self-consistent.
void getTreeNode(ndarray::Array< int, 1, 1 > const &v, int dex) const
Return the _tree information for a given data point.
ndarray::Array< int, 1, 1 > _inn
virtual T operator()(ndarray::Array< const T, 1, 1 > const &, ndarray::Array< const T, 1, 1 > const &) const
Actually evaluate the covariogram function relating two points you want to interpolate from...
void addToSearch()
Adds time to _searchTime.
int _walkUpTree(int target, int dir, int root) const
A method to make sure that every data point in the tree is in the correct position relative to its pa...
void removePoint(int dex)
This will remove a point from the data set.
int _testTree() const
Make sure that the tree is properly constructed. Returns 1 of it is. Return zero if not...
ndarray::Array< T, 2, 2 > _data
virtual T operator()(ndarray::Array< const T, 1, 1 > const &, ndarray::Array< const T, 1, 1 > const &) const
Actually evaluate the covariogram function relating two points you want to interpolate from...
KdTree & operator=(const KdTree &)=delete
virtual T operator()(ndarray::Array< const T, 1, 1 > const &p1, ndarray::Array< const T, 1, 1 > const &p2) const
Actually evaluate the covariogram function relating two points you want to interpolate from...
void addPoint(ndarray::Array< T, 1, 1 > const &vin, T f)
Add a point to the pool of data used by GaussianProcess for interpolation.
GaussianProcess & operator=(const GaussianProcess &)=delete
Citizen is a class that should be among all LSST classes base classes, and handles basic memory manag...
Definition: Citizen.h:53
Include files required for standard LSST Exception handling.
int getPoints() const
return the number of data points stored in the tree
T interpolate(ndarray::Array< T, 1, 1 > variance, ndarray::Array< T, 1, 1 > const &vin, int numberOfNeighbors) const
Interpolate the function value at one point using a specified number of nearest neighbors.
void setEllSquared(double ellSquared)
set the _ellSquared hyper parameter (the square of the characteristic length scale of the covariogram...
void display()
Displays the current values of all times and _interpolationCount.
void _reassign(int target)
Reassign nodes to the tree that were severed by a call to removePoint(int)
void reset()
Resets all of the data members of GaussianProcessTimer to zero.
void addToEigen()
Adds time to _eigenTime.