lsst.meas.astrom  14.0-7-g0d69b06+3
LeastSqFitter2d.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 LEAST_SQ_FITTER_2D
27 #define LEAST_SQ_FITTER_2D
28 
29 
30 #include <cstdio>
31 #include <memory>
32 #include <vector>
33 
34 #include "Eigen/Core"
35 #include "Eigen/SVD"
36 
39 
40 namespace except = lsst::pex::exceptions;
41 
42 
43 namespace lsst {
44 namespace meas {
45 namespace astrom {
46 namespace sip {
47 
69 template <class FittingFunc>class LeastSqFitter2d {
70 public:
72  const std::vector<double> &s, int order);
73 
74  Eigen::MatrixXd getParams();
75  Eigen::MatrixXd getErrors();
76  double valueAt(double x, double y);
78 
79  double getChiSq();
80  double getReducedChiSq();
81 
82 private:
83  void initFunctions();
84 
85  Eigen::MatrixXd expandParams(Eigen::VectorXd const & input) const;
86 
87  double func2d(double x, double y, int term);
88  double func1d(double value, int exponent);
89 
90  std::vector<double> _x, _y, _z, _s;
91  int _order; //Degree of polynomial to fit, e.g 4=> cubic
92  int _nPar; //Number of parameters in fitting eqn, e.g x^2, xy, y^2, x^3,
93  int _nData; //Number of data points, == _x.size()
94 
95  Eigen::JacobiSVD<Eigen::MatrixXd> _svd;
96  Eigen::VectorXd _par;
97 
99 
100 
101 };
102 
103 
104 
105 
106 //The .cc part
107 namespace sip = lsst::meas::astrom::sip;
108 namespace math = lsst::afw::math;
109 
120  const std::vector<double> &y, const std::vector<double> &z, const std::vector<double> &s, int order) :
121  _x(x), _y(y), _z(z), _s(s), _order(order), _nPar(0), _par(1) {
122 
123  //_nPar, the number of terms to fix (x^2, xy, y^2 etc.) is \Sigma^(order+1) 1
124  _nPar = 0;
125  for (int i = 0; i < order; ++i) {
126  _nPar += i + 1;
127  }
128 
129  //Check input vectors are the same size
130  _nData = _x.size();
131  if (_nData != static_cast<int>(_y.size())) {
132  throw LSST_EXCEPT(except::RuntimeError, "x and y vectors of different lengths");
133  }
134  if (_nData != static_cast<int>(_s.size())) {
135  throw LSST_EXCEPT(except::RuntimeError, "x and s vectors of different lengths");
136  }
137  if (_nData != static_cast<int>(_z.size())) {
138  throw LSST_EXCEPT(except::RuntimeError, "x and z vectors of different lengths");
139  }
140 
141  for (int i = 0; i < _nData; ++i) {
142  if ( _s[i] == 0.0 ) {
143  std::string msg = "Illegal zero value for fit weight encountered.";
144  throw LSST_EXCEPT(except::RuntimeError, msg);
145  }
146  }
147 
148  if (_nData < _order) {
149  throw LSST_EXCEPT(except::RuntimeError, "Fewer data points than parameters");
150  }
151 
152  initFunctions();
153 
154  Eigen::MatrixXd design(_nData, _nPar);
155  Eigen::VectorXd rhs(_nData);
156  for (int i = 0; i < _nData; ++i) {
157  rhs[i] = z[i] / s[i];
158  for (int j = 0; j < _nPar; ++j) {
159  design(i, j) = func2d(_x[i], _y[i], j) / _s[i];
160  }
161  }
162  _svd.compute(design, Eigen::ComputeThinU | Eigen::ComputeThinV);
163  _par = _svd.solve(rhs);
164 }
165 
166 
176 template<class FittingFunc> Eigen::MatrixXd LeastSqFitter2d<FittingFunc>::getParams() {
177  return expandParams(_par);
178 }
179 
181 template<class FittingFunc>
182 Eigen::MatrixXd LeastSqFitter2d<FittingFunc>::expandParams(Eigen::VectorXd const & input) const {
183  Eigen::MatrixXd out = Eigen::MatrixXd::Zero(_order, _order);
184  int count = 0;
185  for (int i = 0; i < _order; ++i) {
186  for (int j = 0; j < _order - i; ++j) {
187  out(i, j) = input[count++];
188  }
189  }
190  return out;
191 }
192 
195 template<class FittingFunc> double LeastSqFitter2d<FittingFunc>::getChiSq() {
196 
197  double chisq = 0;
198  for (int i = 0; i < _nData; ++i) {
199  double val = _z[i] - valueAt(_x[i], _y[i]);
200  val /= _s[i];
201  chisq += pow(val, 2);
202  }
203 
204  return chisq;
205 }
206 
207 
213 template<class FittingFunc> double LeastSqFitter2d<FittingFunc>::getReducedChiSq() {
214  return getChiSq()/(double) (_nData - _nPar);
215 }
216 
217 
218 
220 template<class FittingFunc> double LeastSqFitter2d<FittingFunc>::valueAt(double x, double y){
221  double val = 0;
222 
223  //Sum the values of the different orders to get the value of the fitted function
224  for (int i = 0; i < _nPar; ++i) {
225  val += _par[i] * func2d(x, y, i);
226  }
227  return val;
228 }
229 
233 
235  out.reserve(_nData);
236 
237  for (int i = 0; i < _nData; ++i) {
238  out.push_back(_z[i] - valueAt(_x[i], _y[i]));
239  }
240 
241  return out;
242 }
243 
244 
246 template<class FittingFunc>
248  Eigen::ArrayXd variance(_nPar);
249  for (int i = 0; i < _nPar; ++i) {
250  variance[i] = _svd.matrixV().row(i).dot(
251  (_svd.singularValues().array().inverse().square() * _svd.matrixV().col(i).array()).matrix()
252  );
253  }
254  return expandParams(variance.sqrt().matrix());
255 }
256 
257 
258 template<class FittingFunc> void LeastSqFitter2d<FittingFunc>::initFunctions() {
259  //Initialise the array of functions. _funcArray[i] is a object of type math::Function1 of order i
260  _funcArray.reserve(_order);
261 
262  std::vector<double> coeff;
263  coeff.reserve( _order);
264 
265  coeff.push_back(1);
266  for (int i = 0; i < _order; ++i) {
267  std::shared_ptr<FittingFunc> p(new FittingFunc(coeff));
268  _funcArray.push_back(p);
269 
270  coeff[i] = 0;
271  coeff.push_back(1); //coeff now looks like [0,0,...,0,1]
272  }
273 }
274 
275 
276 //The ith term in the fitting polynomial is of the form x^a * y^b. This function figures
277 //out the value of a and b, then calculates the value of the ith term at the given x and y
278 template<class FittingFunc> double LeastSqFitter2d<FittingFunc>::func2d(double x, double y, int term) {
279 
280  int yexp = 0; //y exponent
281  int xexp = 0; //x exponent
282 
283  for (int i = 0; i<term; ++i) {
284  yexp = (yexp + 1) % (_order - xexp);
285  if ( yexp == 0){
286  xexp++;
287  }
288  }
289 
290  double xcomp = func1d(x, xexp); //x component of polynomial
291  double ycomp = func1d(y, yexp); //y component
292 
293  #if 0 //A useful debugging printf statement
294  printf("The %i(th) function: x^%i * y^%i = %.1f * %.1f\n", term, xexp, yexp, xcomp, ycomp);
295  #endif
296 
297  return xcomp*ycomp;
298 }
299 
300 template<class FittingFunc> double LeastSqFitter2d<FittingFunc>::func1d(double value, int exponent) {
301  return (*_funcArray[exponent])(value);
302 }
303 
304 
305 
306 }}}}
307 
308 #endif
Fit an lsst::afw::math::Function1 object to a set of data points in two dimensions.
afw::table::Key< afw::table::Array< VariancePixelT > > variance
double getReducedChiSq()
Return a measure of the goodness of fit.
std::vector< double > residuals()
Return a vector of residuals of the fit (i.e the difference between the input z values, and the value of the fitting function at that point.
double getChiSq()
Return a measure of the goodness of fit.
double valueAt(double x, double y)
Return the value of the best fit function at a given position (x,y)
STL class.
T push_back(T... args)
Eigen::MatrixXd getParams()
Build up a triangular matrix of the parameters.
double x
T size(T... args)
#define LSST_EXCEPT(type,...)
Eigen::MatrixXd getErrors()
Companion function to getParams(). Returns uncertainties in the parameters as a matrix.
LeastSqFitter2d(const std::vector< double > &x, const std::vector< double > &y, const std::vector< double > &z, const std::vector< double > &s, int order)
Fit a 2d polynomial to a set of data points z(x, y)
T reserve(T... args)