LSST Applications  21.0.0-147-g0e635eb1+1acddb5be5,22.0.0+052faf71bd,22.0.0+1ea9a8b2b2,22.0.0+6312710a6c,22.0.0+729191ecac,22.0.0+7589c3a021,22.0.0+9f079a9461,22.0.1-1-g7d6de66+b8044ec9de,22.0.1-1-g87000a6+536b1ee016,22.0.1-1-g8e32f31+6312710a6c,22.0.1-10-gd060f87+016f7cdc03,22.0.1-12-g9c3108e+df145f6f68,22.0.1-16-g314fa6d+c825727ab8,22.0.1-19-g93a5c75+d23f2fb6d8,22.0.1-19-gb93eaa13+aab3ef7709,22.0.1-2-g8ef0a89+b8044ec9de,22.0.1-2-g92698f7+9f079a9461,22.0.1-2-ga9b0f51+052faf71bd,22.0.1-2-gac51dbf+052faf71bd,22.0.1-2-gb66926d+6312710a6c,22.0.1-2-gcb770ba+09e3807989,22.0.1-20-g32debb5+b8044ec9de,22.0.1-23-gc2439a9a+fb0756638e,22.0.1-3-g496fd5d+09117f784f,22.0.1-3-g59f966b+1e6ba2c031,22.0.1-3-g849a1b8+f8b568069f,22.0.1-3-gaaec9c0+c5c846a8b1,22.0.1-32-g5ddfab5d3+60ce4897b0,22.0.1-4-g037fbe1+64e601228d,22.0.1-4-g8623105+b8044ec9de,22.0.1-5-g096abc9+d18c45d440,22.0.1-5-g15c806e+57f5c03693,22.0.1-7-gba73697+57f5c03693,master-g6e05de7fdc+c1283a92b8,master-g72cdda8301+729191ecac,w.2021.39
LSST Data Management Base Package
Public Member Functions | List of all members
lsst::meas::astrom::sip::LeastSqFitter1d< FittingFunc > Class Template Reference

Fit an lsst::afw::math::Function1 object to a set of data points in one dimension. More...

#include <LeastSqFitter1d.h>

Public Member Functions

 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) More...
 
Eigen::VectorXd getParams ()
 Return the best fit parameters as an Eigen::Matrix. More...
 
Eigen::VectorXd getErrors ()
 Return the 1 sigma uncertainties in the best fit parameters as an Eigen::Matrix. More...
 
FittingFunc getBestFitFunction ()
 Return the best fit polynomial as a lsst::afw::math::Function1 object. More...
 
double valueAt (double x)
 Calculate the value of the function at a given point. More...
 
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. More...
 
double getChiSq ()
 Return a measure of the goodness of fit. More...
 
double getReducedChiSq ()
 Return a measure of the goodness of fit. More...
 

Detailed Description

template<class FittingFunc>
class lsst::meas::astrom::sip::LeastSqFitter1d< FittingFunc >

Fit an lsst::afw::math::Function1 object to a set of data points in one dimension.

The class is templated over the kind of object to fit.

Input is a list of x ordinates for a set of points, the y coordinate, and the uncertainties, s. order is order of the polynomial to fit (e.g if the templated function is lsst::afw::math::PolynomialFunction1, then order=3 => fit a function of the form \(ax^2+bx+c\)

Template Parameters
FittingFuncThe 1d function to fit in both dimensions. Must inherit from lsst::afw::math::Function1
Parameters
xOrdinate of points to fit
yCo-ordinate of pionts to fit
s1 \(\sigma\) uncertainties in z
orderPolynomial order to fit
See also
LeastSqFitter1d

Definition at line 64 of file LeastSqFitter1d.h.

Constructor & Destructor Documentation

◆ LeastSqFitter1d()

template<class FittingFunc >
lsst::meas::astrom::sip::LeastSqFitter1d< FittingFunc >::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)

Template Parameters
FittingFuncThe type of function to fit. This function extends the base class of lsst::afw::math::Function1
Parameters
xvector of x positions of data
yvector of y positions of data
sVector of measured uncertainties in the values of z
orderOrder of 2d function to fit

Definition at line 104 of file LeastSqFitter1d.h.

106  : _x(x), _y(y), _s(s), _order(order) {
107  if (order == 0) {
108  throw LSST_EXCEPT(pex::exceptions::RuntimeError, "Fit order must be >= 1");
109  }
110 
111  _nData = _x.size();
112  if (_nData != static_cast<int>(_y.size())) {
113  throw LSST_EXCEPT(pex::exceptions::RuntimeError, "x and y vectors of different lengths");
114  }
115  if (_nData != static_cast<int>(_s.size())) {
116  throw LSST_EXCEPT(pex::exceptions::RuntimeError, "x and s vectors of different lengths");
117  }
118 
119  if (_nData < _order) {
120  throw LSST_EXCEPT(pex::exceptions::RuntimeError, "Fewer data points than parameters");
121  }
122 
123  initFunctions();
124 
125  Eigen::MatrixXd design(_nData, _order);
126  Eigen::VectorXd rhs(_nData);
127  for (int i = 0; i < _nData; ++i) {
128  rhs[i] = y[i] / _s[i];
129  for (int j = 0; j < _order; ++j) {
130  design(i, j) = func1d(_x[i], j) / _s[i];
131  }
132  }
133  _svd.compute(design, Eigen::ComputeThinU | Eigen::ComputeThinV);
134  _par = _svd.solve(rhs);
135 }
double x
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
int y
Definition: SpanSet.cc:48
T size(T... args)
table::Key< int > order

Member Function Documentation

◆ getBestFitFunction()

template<class FittingFunc >
FittingFunc lsst::meas::astrom::sip::LeastSqFitter1d< FittingFunc >::getBestFitFunction

Return the best fit polynomial as a lsst::afw::math::Function1 object.

Definition at line 160 of file LeastSqFitter1d.h.

160  {
161  // FittingFunc and LeastSqFitter disagree on the definition of order of a function.
162  // LSF says that a linear function is order 2 (two coefficients), FF says only 1
163  FittingFunc func(_order - 1);
164 
165  for (int i = 0; i < _order; ++i) {
166  func.setParameter(i, _par(i));
167  }
168  return func;
169 }

◆ getChiSq()

template<class FittingFunc >
double lsst::meas::astrom::sip::LeastSqFitter1d< FittingFunc >::getChiSq

Return a measure of the goodness of fit.

\[ \chi_r^2 = \sum \left( \frac{y_i - f(x_i)}{s} \right)^2 \]

Definition at line 199 of file LeastSqFitter1d.h.

199  {
200  FittingFunc f = getBestFitFunction();
201 
202  double chisq = 0;
203  for (int i = 0; i < _nData; ++i) {
204  double val = _y[i] - f(_x[i]);
205  val /= _s[i];
206  chisq += pow(val, 2);
207  }
208 
209  return chisq;
210 }
FittingFunc getBestFitFunction()
Return the best fit polynomial as a lsst::afw::math::Function1 object.
T pow(T... args)
ImageT val
Definition: CR.cc:146

◆ getErrors()

template<class FittingFunc >
Eigen::VectorXd lsst::meas::astrom::sip::LeastSqFitter1d< FittingFunc >::getErrors

Return the 1 sigma uncertainties in the best fit parameters as an Eigen::Matrix.

Definition at line 149 of file LeastSqFitter1d.h.

149  {
150  Eigen::ArrayXd variance(_order);
151  for (int i = 0; i < _order; ++i) {
152  variance[i] = _svd.matrixV().row(i).dot(
153  (_svd.singularValues().array().inverse().square() * _svd.matrixV().col(i).array()).matrix());
154  }
155  return variance.sqrt().matrix();
156 }
afw::table::Key< afw::table::Array< VariancePixelT > > variance

◆ getParams()

template<class FittingFunc >
Eigen::VectorXd lsst::meas::astrom::sip::LeastSqFitter1d< FittingFunc >::getParams

Return the best fit parameters as an Eigen::Matrix.

Definition at line 139 of file LeastSqFitter1d.h.

139  {
140  Eigen::VectorXd vec = Eigen::VectorXd::Zero(_order);
141  for (int i = 0; i < _order; ++i) {
142  vec(i) = _par(i);
143  }
144  return vec;
145 }

◆ getReducedChiSq()

template<class FittingFunc >
double lsst::meas::astrom::sip::LeastSqFitter1d< FittingFunc >::getReducedChiSq

Return a measure of the goodness of fit.

\[ \chi_r^2 = \sum \left( \frac{y_i - f(x_i)}{s} \right)^2 \div (N-p) \]

Where \( N \) is the number of data points, and \( p \) is the number of parameters in the fit

Definition at line 218 of file LeastSqFitter1d.h.

218  {
219  return getChiSq() / (double)(_nData - _order);
220 }
double getChiSq()
Return a measure of the goodness of fit.

◆ residuals()

template<class FittingFunc >
std::vector< double > lsst::meas::astrom::sip::LeastSqFitter1d< FittingFunc >::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.

Definition at line 182 of file LeastSqFitter1d.h.

182  {
184  out.reserve(_nData);
185 
186  FittingFunc f = getBestFitFunction();
187 
188  for (int i = 0; i < _nData; ++i) {
189  out.push_back(_y[i] - f(_x[i]));
190  }
191 
192  return out;
193 }
T push_back(T... args)
T reserve(T... args)

◆ valueAt()

template<class FittingFunc >
double lsst::meas::astrom::sip::LeastSqFitter1d< FittingFunc >::valueAt ( double  x)

Calculate the value of the function at a given point.

Definition at line 173 of file LeastSqFitter1d.h.

173  {
174  FittingFunc f = getBestFitFunction();
175 
176  return f(x);
177 }

The documentation for this class was generated from the following file: