lsst.meas.astrom  14.0-7-g0d69b06+3
LeastSqFitter1d.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_1D
27 #define LEAST_SQ_FITTER_1D
28 
29 #include <cstdio>
30 #include <memory>
31 #include <vector>
32 
33 #include "Eigen/Core"
34 #include "Eigen/SVD"
35 
38 
39 namespace except = lsst::pex::exceptions;
40 
41 
42 namespace lsst {
43 namespace meas {
44 namespace astrom {
45 namespace sip {
46 
67 template <class FittingFunc>class LeastSqFitter1d {
68 public:
70  const std::vector<double> &s, int order);
71 
72  Eigen::VectorXd getParams();
73  Eigen::VectorXd getErrors();
74  FittingFunc getBestFitFunction();
75  double valueAt(double x);
77 
78  double getChiSq();
79  double getReducedChiSq();
80 
81 
82 private:
83  void initFunctions();
84 
85  double func1d(double value, int exponent);
86 
87  std::vector<double> _x, _y, _s;
88  int _order; //Degree of polynomial to fit, e.g 4=> cubic
89  int _nData; //Number of data points, == _x.size()
90 
91  Eigen::JacobiSVD<Eigen::MatrixXd> _svd;
92  Eigen::VectorXd _par;
93 
95 };
96 
97 //The .cc part
98 
108  const std::vector<double> &y, const std::vector<double> &s, int order) :
109  _x(x), _y(y), _s(s), _order(order) {
110 
111  if (order == 0) {
112  throw LSST_EXCEPT(except::RuntimeError, "Fit order must be >= 1");
113  }
114 
115  _nData = _x.size();
116  if (_nData != static_cast<int>(_y.size())) {
117  throw LSST_EXCEPT(except::RuntimeError, "x and y vectors of different lengths");
118  }
119  if (_nData != static_cast<int>(_s.size())) {
120  throw LSST_EXCEPT(except::RuntimeError, "x and s vectors of different lengths");
121  }
122 
123  if (_nData < _order) {
124  throw LSST_EXCEPT(except::RuntimeError, "Fewer data points than parameters");
125  }
126 
127  initFunctions();
128 
129  Eigen::MatrixXd design(_nData, _order);
130  Eigen::VectorXd rhs(_nData);
131  for (int i = 0; i < _nData; ++i) {
132  rhs[i] = y[i] / _s[i];
133  for (int j = 0; j < _order; ++j) {
134  design(i, j) = func1d(_x[i], j) / _s[i];
135  }
136  }
137  _svd.compute(design, Eigen::ComputeThinU | Eigen::ComputeThinV);
138  _par = _svd.solve(rhs);
139 }
140 
141 
142 
144 template<class FittingFunc> Eigen::VectorXd LeastSqFitter1d<FittingFunc>::getParams() {
145 
146  Eigen::VectorXd vec = Eigen::VectorXd::Zero(_order);
147  for (int i = 0; i < _order; ++i) {
148  vec(i) = _par(i);
149  }
150  return vec;
151 }
152 
153 
155 template<class FittingFunc> Eigen::VectorXd LeastSqFitter1d<FittingFunc>::getErrors() {
156  Eigen::ArrayXd variance(_order);
157  for (int i = 0; i < _order; ++i) {
158  variance[i] = _svd.matrixV().row(i).dot(
159  (_svd.singularValues().array().inverse().square() * _svd.matrixV().col(i).array()).matrix()
160  );
161  }
162  return variance.sqrt().matrix();
163 }
164 
165 
167 template<class FittingFunc> FittingFunc LeastSqFitter1d<FittingFunc>::getBestFitFunction() {
168 
169  //FittingFunc and LeastSqFitter disagree on the definition of order of a function.
170  //LSF says that a linear function is order 2 (two coefficients), FF says only 1
171  FittingFunc func(_order - 1);
172 
173  for (int i = 0; i < _order; ++i) {
174  func.setParameter(i, _par(i));
175  }
176  return func;
177 }
178 
179 
181 template<class FittingFunc> double LeastSqFitter1d<FittingFunc>::valueAt(double x) {
182  FittingFunc f = getBestFitFunction();
183 
184  return f(x);
185 }
186 
187 
188 
193  out.reserve(_nData);
194 
195  FittingFunc f = getBestFitFunction();
196 
197  for (int i = 0; i < _nData; ++i) {
198  out.push_back(_y[i] - f(_x[i]));
199  }
200 
201  return out;
202 }
203 
204 
208 template<class FittingFunc> double LeastSqFitter1d<FittingFunc>::getChiSq() {
209  FittingFunc f = getBestFitFunction();
210 
211  double chisq = 0;
212  for (int i = 0; i < _nData; ++i) {
213  double val = _y[i] - f(_x[i]);
214  val /= _s[i];
215  chisq += pow(val, 2);
216  }
217 
218  return chisq;
219 }
220 
221 
227 template<class FittingFunc> double LeastSqFitter1d<FittingFunc>::getReducedChiSq() {
228  return getChiSq()/(double) (_nData - _order);
229 }
230 
231 
232 
235 template<class FittingFunc> void LeastSqFitter1d<FittingFunc>::initFunctions() {
236 
237  _funcArray.reserve(_order);
238 
239  std::vector<double> coeff;
240  coeff.reserve( _order);
241 
242  coeff.push_back(1.0);
243  for (int i = 0; i < _order; ++i) {
244  std::shared_ptr<FittingFunc> p(new FittingFunc(coeff));
245  _funcArray.push_back(p);
246  coeff[i] = 0.0;
247  coeff.push_back(1.0); //coeff now looks like [0,0,...,0,1]
248  }
249 }
250 
251 template<class FittingFunc> double LeastSqFitter1d<FittingFunc>::func1d(double value, int exponent) {
252  return (*_funcArray[exponent])(value);
253 }
254 
255 }}}}
256 
257 #endif
double getReducedChiSq()
Return a measure of the goodness of fit.
afw::table::Key< afw::table::Array< VariancePixelT > > variance
std::vector< double > residuals()
Return a vector of residuals of the fit (i.e the difference between the input y values, and the value of the fitting function at that point.
Eigen::VectorXd getParams()
Return the best fit parameters as an Eigen::Matrix.
T push_back(T... args)
Eigen::VectorXd getErrors()
Return the 1 sigma uncertainties in the best fit parameters as an Eigen::Matrix.
double valueAt(double x)
Calculate the value of the function at a given point.
FittingFunc getBestFitFunction()
Return the best fit polynomial as a lsst::afw::math::Function1 object.
LeastSqFitter1d(const std::vector< double > &x, const std::vector< double > &y, const std::vector< double > &s, int order)
Fit a 1d polynomial to a set of data points z(x, y)
double x
T size(T... args)
#define LSST_EXCEPT(type,...)
double getChiSq()
Return a measure of the goodness of fit.
T reserve(T... args)
Fit an lsst::afw::math::Function1 object to a set of data points in one dimension.