LSSTApplications  10.0+286,10.0+36,10.0+46,10.0-2-g4f67435,10.1+152,10.1+37,11.0,11.0+1,11.0-1-g47edd16,11.0-1-g60db491,11.0-1-g7418c06,11.0-2-g04d2804,11.0-2-g68503cd,11.0-2-g818369d,11.0-2-gb8b8ce7
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 <vector>
31 
32 #include "boost/shared_ptr.hpp"
33 #include "Eigen/Core"
34 #include "Eigen/SVD"
35 
37 #include "lsst/pex/logging/Trace.h"
39 
40 namespace except = lsst::pex::exceptions;
41 namespace pexLogging = lsst::pex::logging;
42 
43 
44 namespace lsst {
45 namespace meas {
46 namespace astrom {
47 namespace sip {
48 
69 template <class FittingFunc>class LeastSqFitter1d {
70 public:
71  LeastSqFitter1d(const std::vector<double> &x, const std::vector<double> &y,
72  const std::vector<double> &s,
73  unsigned int order);
74 
75  Eigen::VectorXd getParams();
76  Eigen::VectorXd getErrors();
77  FittingFunc getBestFitFunction();
78  double valueAt(double x);
79  std::vector<double> residuals();
80 
81  double getChiSq();
82  double getReducedChiSq();
83 
84 
85 private:
86  void initFunctions();
87 
88  double func1d(double value, int exponent);
89 
90  std::vector<double> _x, _y, _s;
91  int _order; //Degree of polynomial to fit, e.g 4=> cubic
92  int _nData; //Number of data points, == _x.size()
93 
94  Eigen::JacobiSVD<Eigen::MatrixXd> _svd;
95  Eigen::VectorXd _par;
96 
97  std::vector<boost::shared_ptr<FittingFunc> > _funcArray;
98 };
99 
100 //The .cc part
101 
110 template<class FittingFunc> LeastSqFitter1d<FittingFunc>::LeastSqFitter1d(const std::vector<double> &x,
111  const std::vector<double> &y, const std::vector<double> &s, unsigned int order) :
112  _x(x), _y(y), _s(s), _order(order) {
113 
114  if (order == 0) {
115  throw LSST_EXCEPT(except::RuntimeError, "Fit order must be >= 1");
116  }
117 
118  _nData = _x.size();
119  if (_nData != static_cast<int>(_y.size())) {
120  throw LSST_EXCEPT(except::RuntimeError, "x and y vectors of different lengths");
121  }
122  if (_nData != static_cast<int>(_s.size())) {
123  throw LSST_EXCEPT(except::RuntimeError, "x and s vectors of different lengths");
124  }
125 
126  if (_nData < _order) {
127  throw LSST_EXCEPT(except::RuntimeError, "Fewer data points than parameters");
128  }
129 
130  initFunctions();
131 
132  Eigen::MatrixXd design(_nData, _order);
133  Eigen::VectorXd rhs(_nData);
134  for (int i = 0; i < _nData; ++i) {
135  rhs[i] = y[i] / _s[i];
136  for (int j = 0; j < _order; ++j) {
137  design(i, j) = func1d(_x[i], j) / _s[i];
138  }
139  }
140  _svd.compute(design, Eigen::ComputeThinU | Eigen::ComputeThinV);
141  _par = _svd.solve(rhs);
142 }
143 
144 
145 
147 template<class FittingFunc> Eigen::VectorXd LeastSqFitter1d<FittingFunc>::getParams() {
148 
149  Eigen::VectorXd vec = Eigen::VectorXd::Zero(_order);
150  for (int i = 0; i < _order; ++i) {
151  vec(i) = _par(i);
152  }
153  return vec;
154 }
155 
156 
158 template<class FittingFunc> Eigen::VectorXd LeastSqFitter1d<FittingFunc>::getErrors() {
159  Eigen::ArrayXd variance(_order);
160  for (int i = 0; i < _order; ++i) {
161  variance[i] = _svd.matrixV().row(i).dot(
162  (_svd.singularValues().array().inverse().square() * _svd.matrixV().col(i).array()).matrix()
163  );
164  }
165  return variance.sqrt().matrix();
166 }
167 
168 
170 template<class FittingFunc> FittingFunc LeastSqFitter1d<FittingFunc>::getBestFitFunction() {
171 
172  //FittingFunc and LeastSqFitter disagree on the definition of order of a function.
173  //LSF says that a linear function is order 2 (two coefficients), FF says only 1
174  FittingFunc func(_order - 1);
175 
176  for (int i = 0; i < _order; ++i) {
177  func.setParameter(i, _par(i));
178  }
179  return func;
180 }
181 
182 
184 template<class FittingFunc> double LeastSqFitter1d<FittingFunc>::valueAt(double x) {
185  FittingFunc f = getBestFitFunction();
186 
187  return f(x);
188 }
189 
190 
191 
194 template<class FittingFunc> std::vector<double> LeastSqFitter1d<FittingFunc>::residuals(){
195  std::vector<double> out;
196  out.reserve(_nData);
197 
198  FittingFunc f = getBestFitFunction();
199 
200  for (int i = 0; i < _nData; ++i) {
201  out.push_back(_y[i] - f(_x[i]));
202  }
203 
204  return out;
205 }
206 
207 
211 template<class FittingFunc> double LeastSqFitter1d<FittingFunc>::getChiSq() {
212  FittingFunc f = getBestFitFunction();
213 
214  double chisq = 0;
215  for (int i = 0; i < _nData; ++i) {
216  double val = _y[i] - f(_x[i]);
217  val /= _s[i];
218  chisq += pow(val, 2);
219  }
220 
221  return chisq;
222 }
223 
224 
230 template<class FittingFunc> double LeastSqFitter1d<FittingFunc>::getReducedChiSq() {
231  return getChiSq()/(double) (_nData - _order);
232 }
233 
234 
235 
238 template<class FittingFunc> void LeastSqFitter1d<FittingFunc>::initFunctions() {
239 
240  _funcArray.reserve(_order);
241 
242  std::vector<double> coeff;
243  coeff.reserve( _order);
244 
245  coeff.push_back(1.0);
246  for (int i = 0; i < _order; ++i) {
247  boost::shared_ptr<FittingFunc> p(new FittingFunc(coeff));
248  _funcArray.push_back(p);
249  coeff[i] = 0.0;
250  coeff.push_back(1.0); //coeff now looks like [0,0,...,0,1]
251  }
252 }
253 
254 template<class FittingFunc> double LeastSqFitter1d<FittingFunc>::func1d(double value, int exponent) {
255  return (*_funcArray[exponent])(value);
256 }
257 
258 }}}}
259 
260 #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.
definition of the Trace messaging facilities
double getReducedChiSq()
Return a measure of the goodness of fit. Where is the number of data points, and is the number of ...
Define a collection of useful Functions.
bool val
LeastSqFitter1d(const std::vector< double > &x, const std::vector< double > &y, const std::vector< double > &s, unsigned int order)
Eigen::JacobiSVD< Eigen::MatrixXd > _svd
double func1d(double value, int exponent)
double x
double chisq
#define LSST_EXCEPT(type,...)
Definition: Exception.h:46
double valueAt(double x)
Calculate the value of the function at a given point.
Eigen::JacobiSVD< Eigen::MatrixXd > _svd
std::vector< boost::shared_ptr< FittingFunc > > _funcArray
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.