LSSTApplications  11.0-13-gbb96280,12.1+18,12.1+7,12.1-1-g14f38d3+72,12.1-1-g16c0db7+5,12.1-1-g5961e7a+84,12.1-1-ge22e12b+23,12.1-11-g06625e2+4,12.1-11-g0d7f63b+4,12.1-19-gd507bfc,12.1-2-g7dda0ab+38,12.1-2-gc0bc6ab+81,12.1-21-g6ffe579+2,12.1-21-gbdb6c2a+4,12.1-24-g941c398+5,12.1-3-g57f6835+7,12.1-3-gf0736f3,12.1-37-g3ddd237,12.1-4-gf46015e+5,12.1-5-g06c326c+20,12.1-5-g648ee80+3,12.1-5-gc2189d7+4,12.1-6-ga608fc0+1,12.1-7-g3349e2a+5,12.1-7-gfd75620+9,12.1-9-g577b946+5,12.1-9-gc4df26a+10
LSSTDataManagementBasePackage
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:
69  LeastSqFitter1d(const std::vector<double> &x, const std::vector<double> &y,
70  const std::vector<double> &s,
71  unsigned int order);
72 
73  Eigen::VectorXd getParams();
74  Eigen::VectorXd getErrors();
75  FittingFunc getBestFitFunction();
76  double valueAt(double x);
77  std::vector<double> residuals();
78 
79  double getChiSq();
80  double getReducedChiSq();
81 
82 
83 private:
84  void initFunctions();
85 
86  double func1d(double value, int exponent);
87 
88  std::vector<double> _x, _y, _s;
89  int _order; //Degree of polynomial to fit, e.g 4=> cubic
90  int _nData; //Number of data points, == _x.size()
91 
92  Eigen::JacobiSVD<Eigen::MatrixXd> _svd;
93  Eigen::VectorXd _par;
94 
95  std::vector<std::shared_ptr<FittingFunc> > _funcArray;
96 };
97 
98 //The .cc part
99 
108 template<class FittingFunc> LeastSqFitter1d<FittingFunc>::LeastSqFitter1d(const std::vector<double> &x,
109  const std::vector<double> &y, const std::vector<double> &s, unsigned int order) :
110  _x(x), _y(y), _s(s), _order(order) {
111 
112  if (order == 0) {
113  throw LSST_EXCEPT(except::RuntimeError, "Fit order must be >= 1");
114  }
115 
116  _nData = _x.size();
117  if (_nData != static_cast<int>(_y.size())) {
118  throw LSST_EXCEPT(except::RuntimeError, "x and y vectors of different lengths");
119  }
120  if (_nData != static_cast<int>(_s.size())) {
121  throw LSST_EXCEPT(except::RuntimeError, "x and s vectors of different lengths");
122  }
123 
124  if (_nData < _order) {
125  throw LSST_EXCEPT(except::RuntimeError, "Fewer data points than parameters");
126  }
127 
128  initFunctions();
129 
130  Eigen::MatrixXd design(_nData, _order);
131  Eigen::VectorXd rhs(_nData);
132  for (int i = 0; i < _nData; ++i) {
133  rhs[i] = y[i] / _s[i];
134  for (int j = 0; j < _order; ++j) {
135  design(i, j) = func1d(_x[i], j) / _s[i];
136  }
137  }
138  _svd.compute(design, Eigen::ComputeThinU | Eigen::ComputeThinV);
139  _par = _svd.solve(rhs);
140 }
141 
142 
143 
145 template<class FittingFunc> Eigen::VectorXd LeastSqFitter1d<FittingFunc>::getParams() {
146 
147  Eigen::VectorXd vec = Eigen::VectorXd::Zero(_order);
148  for (int i = 0; i < _order; ++i) {
149  vec(i) = _par(i);
150  }
151  return vec;
152 }
153 
154 
156 template<class FittingFunc> Eigen::VectorXd LeastSqFitter1d<FittingFunc>::getErrors() {
157  Eigen::ArrayXd variance(_order);
158  for (int i = 0; i < _order; ++i) {
159  variance[i] = _svd.matrixV().row(i).dot(
160  (_svd.singularValues().array().inverse().square() * _svd.matrixV().col(i).array()).matrix()
161  );
162  }
163  return variance.sqrt().matrix();
164 }
165 
166 
168 template<class FittingFunc> FittingFunc LeastSqFitter1d<FittingFunc>::getBestFitFunction() {
169 
170  //FittingFunc and LeastSqFitter disagree on the definition of order of a function.
171  //LSF says that a linear function is order 2 (two coefficients), FF says only 1
172  FittingFunc func(_order - 1);
173 
174  for (int i = 0; i < _order; ++i) {
175  func.setParameter(i, _par(i));
176  }
177  return func;
178 }
179 
180 
182 template<class FittingFunc> double LeastSqFitter1d<FittingFunc>::valueAt(double x) {
183  FittingFunc f = getBestFitFunction();
184 
185  return f(x);
186 }
187 
188 
189 
192 template<class FittingFunc> std::vector<double> LeastSqFitter1d<FittingFunc>::residuals(){
193  std::vector<double> out;
194  out.reserve(_nData);
195 
196  FittingFunc f = getBestFitFunction();
197 
198  for (int i = 0; i < _nData; ++i) {
199  out.push_back(_y[i] - f(_x[i]));
200  }
201 
202  return out;
203 }
204 
205 
209 template<class FittingFunc> double LeastSqFitter1d<FittingFunc>::getChiSq() {
210  FittingFunc f = getBestFitFunction();
211 
212  double chisq = 0;
213  for (int i = 0; i < _nData; ++i) {
214  double val = _y[i] - f(_x[i]);
215  val /= _s[i];
216  chisq += pow(val, 2);
217  }
218 
219  return chisq;
220 }
221 
222 
228 template<class FittingFunc> double LeastSqFitter1d<FittingFunc>::getReducedChiSq() {
229  return getChiSq()/(double) (_nData - _order);
230 }
231 
232 
233 
236 template<class FittingFunc> void LeastSqFitter1d<FittingFunc>::initFunctions() {
237 
238  _funcArray.reserve(_order);
239 
240  std::vector<double> coeff;
241  coeff.reserve( _order);
242 
243  coeff.push_back(1.0);
244  for (int i = 0; i < _order; ++i) {
245  std::shared_ptr<FittingFunc> p(new FittingFunc(coeff));
246  _funcArray.push_back(p);
247  coeff[i] = 0.0;
248  coeff.push_back(1.0); //coeff now looks like [0,0,...,0,1]
249  }
250 }
251 
252 template<class FittingFunc> double LeastSqFitter1d<FittingFunc>::func1d(double value, int exponent) {
253  return (*_funcArray[exponent])(value);
254 }
255 
256 }}}}
257 
258 #endif
int y
Eigen::VectorXd getErrors()
Return the 1 sigma uncertainties in the best fit parameters as an Eigen::Matrix.
FittingFunc getBestFitFunction()
Return the best fit polynomial as a lsst::afw::math::Function1 object.
double getReducedChiSq()
Return a measure of the goodness of fit.
Define a collection of useful Functions.
LeastSqFitter1d(const std::vector< double > &x, const std::vector< double > &y, const std::vector< double > &s, unsigned int order)
Fit a 1d polynomial to a set of data points z(x, y)
Eigen::JacobiSVD< Eigen::MatrixXd > _svd
double func1d(double value, int exponent)
double x
double chisq
#define LSST_EXCEPT(type,...)
Create an exception with a given type and message and optionally other arguments (dependent on the ty...
Definition: Exception.h:46
void initFunctions()
Initialise the array of functions.
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.
std::vector< std::shared_ptr< FittingFunc > > _funcArray
double valueAt(double x)
Calculate the value of the function at a given point.
ImageT val
Definition: CR.cc:159
Eigen::JacobiSVD< Eigen::MatrixXd > _svd
Eigen::VectorXd getParams()
Return the best fit parameters as an Eigen::Matrix.
double getChiSq()
Return a measure of the goodness of fit.
Fit an lsst::afw::math::Function1 object to a set of data points in one dimension.