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
Public Member Functions | Private Member Functions | Private Attributes | List of all members
lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc > Class Template Reference

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

#include <LeastSqFitter2d.h>

Public Member Functions

 LeastSqFitter2d (const std::vector< double > &x, const std::vector< double > &y, const std::vector< double > &z, const std::vector< double > &s, int order)
 
Eigen::MatrixXd getParams ()
 
Eigen::MatrixXd getErrors ()
 Companion function to getParams(). Returns uncertainties in the parameters as a matrix. More...
 
double valueAt (double x, double y)
 Return the value of the best fit function at a given position (x,y) More...
 
std::vector< double > residuals ()
 
double getChiSq ()
 Return a measure of the goodness of fit.

\[ \chi^2 = \sum \left( \frac{z_i - f(x_i, y_i)}{s_i} \right)^2 \]

. More...

 
double getReducedChiSq ()
 Return a measure of the goodness of fit.

\[ \chi_r^2 = \sum \left( \frac{z_i - f(x_i, y_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. More...

 

Private Member Functions

void initFunctions ()
 
Eigen::MatrixXd expandParams (Eigen::VectorXd const &input) const
 Turn a flattened parameter-like vector into a triangular matrix. More...
 
double func2d (double x, double y, int term)
 
double func1d (double value, int exponent)
 

Private Attributes

std::vector< double > _x
 
std::vector< double > _y
 
std::vector< double > _z
 
std::vector< double > _s
 
int _order
 
int _nPar
 
int _nData
 
Eigen::JacobiSVD< Eigen::MatrixXd > _svd
 
Eigen::VectorXd _par
 
std::vector< boost::shared_ptr
< FittingFunc > > 
_funcArray
 

Detailed Description

template<class FittingFunc>
class lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >

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

The class is templated over the kind of object to fit. Note that we fit the 1d object in each dimension, not the 2d one.

Input is a list of x and y ordinates for a set of points, the z 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
yOrdinate of points to fit
zCo-ordinate of pionts to fit
s1 \(\sigma\) uncertainties in z
orderPolynomial order to fit
See Also
LeastSqFitter1d

Definition at line 71 of file LeastSqFitter2d.h.

Constructor & Destructor Documentation

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

Template Parameters
FittingFuncThe type of function to fit in each dimension. This function should extend the base class of lsst::afw::math::Function1. Note that this is a 1d function
Parameters
xvector of x positions of data
yvector of y positions of data
zValue of data for a given x,y. \(z = z_i = z_i(x_i, y_i)\)
sVector of measured uncertainties in the values of z
orderOrder of 2d function to fit

Definition at line 121 of file LeastSqFitter2d.h.

122  :
123  _x(x), _y(y), _z(z), _s(s), _order(order), _nPar(0), _par(1) {
124 
125  //_nPar, the number of terms to fix (x^2, xy, y^2 etc.) is \Sigma^(order+1) 1
126  _nPar = 0;
127  for (int i = 0; i < order; ++i) {
128  _nPar += i + 1;
129  }
130 
131  //Check input vectors are the same size
132  _nData = _x.size();
133  if (_nData != static_cast<int>(_y.size())) {
134  throw LSST_EXCEPT(except::RuntimeError, "x and y vectors of different lengths");
135  }
136  if (_nData != static_cast<int>(_s.size())) {
137  throw LSST_EXCEPT(except::RuntimeError, "x and s vectors of different lengths");
138  }
139  if (_nData != static_cast<int>(_z.size())) {
140  throw LSST_EXCEPT(except::RuntimeError, "x and z vectors of different lengths");
141  }
142 
143  for (int i = 0; i < _nData; ++i) {
144  if ( _s[i] == 0.0 ) {
145  std::string msg = "Illegal zero value for fit weight encountered.";
146  throw LSST_EXCEPT(except::RuntimeError, msg);
147  }
148  }
149 
150  if (_nData < _order) {
151  throw LSST_EXCEPT(except::RuntimeError, "Fewer data points than parameters");
152  }
153 
154  initFunctions();
155 
156  Eigen::MatrixXd design(_nData, _nPar);
157  Eigen::VectorXd rhs(_nData);
158  for (int i = 0; i < _nData; ++i) {
159  rhs[i] = z[i] / s[i];
160  for (int j = 0; j < _nPar; ++j) {
161  design(i, j) = func2d(_x[i], _y[i], j) / _s[i];
162  }
163  }
164  _svd.compute(design, Eigen::ComputeThinU | Eigen::ComputeThinV);
165  _par = _svd.solve(rhs);
166 }
double func2d(double x, double y, int term)
Eigen::JacobiSVD< Eigen::MatrixXd > _svd
#define LSST_EXCEPT(type,...)
Definition: Exception.h:46

Member Function Documentation

template<class FittingFunc >
Eigen::MatrixXd lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::expandParams ( Eigen::VectorXd const &  input) const
private

Turn a flattened parameter-like vector into a triangular matrix.

Definition at line 184 of file LeastSqFitter2d.h.

184  {
185  Eigen::MatrixXd out = Eigen::MatrixXd::Zero(_order, _order);
186  int count = 0;
187  for (int i = 0; i < _order; ++i) {
188  for (int j = 0; j < _order - i; ++j) {
189  out(i, j) = input[count++];
190  }
191  }
192  return out;
193 }
template<class FittingFunc >
double lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::func1d ( double  value,
int  exponent 
)
private

Definition at line 302 of file LeastSqFitter2d.h.

302  {
303  return (*_funcArray[exponent])(value);
304 }
std::vector< boost::shared_ptr< FittingFunc > > _funcArray
template<class FittingFunc >
double lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::func2d ( double  x,
double  y,
int  term 
)
private

Definition at line 280 of file LeastSqFitter2d.h.

280  {
281 
282  int yexp = 0; //y exponent
283  int xexp = 0; //x exponent
284 
285  for (int i = 0; i<term; ++i) {
286  yexp = (yexp + 1) % (_order - xexp);
287  if ( yexp == 0){
288  xexp++;
289  }
290  }
291 
292  double xcomp = func1d(x, xexp); //x component of polynomial
293  double ycomp = func1d(y, yexp); //y component
294 
295  #if 0 //A useful debugging printf statement
296  printf("The %i(th) function: x^%i * y^%i = %.1f * %.1f\n", term, xexp, yexp, xcomp, ycomp);
297  #endif
298 
299  return xcomp*ycomp;
300 }
int y
double func1d(double value, int exponent)
double x
template<class FittingFunc >
double lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::getChiSq ( )

Return a measure of the goodness of fit.

\[ \chi^2 = \sum \left( \frac{z_i - f(x_i, y_i)}{s_i} \right)^2 \]

.

Definition at line 197 of file LeastSqFitter2d.h.

197  {
198 
199  double chisq = 0;
200  for (int i = 0; i < _nData; ++i) {
201  double val = _z[i] - valueAt(_x[i], _y[i]);
202  val /= _s[i];
203  chisq += pow(val, 2);
204  }
205 
206  return chisq;
207 }
bool val
double valueAt(double x, double y)
Return the value of the best fit function at a given position (x,y)
double chisq
template<class FittingFunc >
Eigen::MatrixXd lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::getErrors ( )

Companion function to getParams(). Returns uncertainties in the parameters as a matrix.

Definition at line 249 of file LeastSqFitter2d.h.

249  {
250  Eigen::ArrayXd variance(_nPar);
251  for (int i = 0; i < _nPar; ++i) {
252  variance[i] = _svd.matrixV().row(i).dot(
253  (_svd.singularValues().array().inverse().square() * _svd.matrixV().col(i).array()).matrix()
254  );
255  }
256  return expandParams(variance.sqrt().matrix());
257 }
Eigen::JacobiSVD< Eigen::MatrixXd > _svd
Eigen::MatrixXd expandParams(Eigen::VectorXd const &input) const
Turn a flattened parameter-like vector into a triangular matrix.
template<class FittingFunc >
Eigen::MatrixXd lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::getParams ( )

Build up a triangular matrix of the parameters. The shape of the matrix is such that the values correspond to the coefficients of the following polynomials
1 y y^2 y^3
x xy xy^2 0
x^2 x^2y 0 0
x^3 0 0 0
(order==4)
where row x column < order

Definition at line 178 of file LeastSqFitter2d.h.

178  {
179  return expandParams(_par);
180 }
Eigen::MatrixXd expandParams(Eigen::VectorXd const &input) const
Turn a flattened parameter-like vector into a triangular matrix.
template<class FittingFunc >
double lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::getReducedChiSq ( )

Return a measure of the goodness of fit.

\[ \chi_r^2 = \sum \left( \frac{z_i - f(x_i, y_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 215 of file LeastSqFitter2d.h.

215  {
216  return getChiSq()/(double) (_nData - _nPar);
217 }
double getChiSq()
Return a measure of the goodness of fit. .
template<class FittingFunc >
void lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::initFunctions ( )
private

Definition at line 260 of file LeastSqFitter2d.h.

260  {
261  //Initialise the array of functions. _funcArray[i] is a object of type math::Function1 of order i
262  _funcArray.reserve(_order);
263 
264  std::vector<double> coeff;
265  coeff.reserve( _order);
266 
267  coeff.push_back(1);
268  for (int i = 0; i < _order; ++i) {
269  boost::shared_ptr<FittingFunc> p(new FittingFunc(coeff));
270  _funcArray.push_back(p);
271 
272  coeff[i] = 0;
273  coeff.push_back(1); //coeff now looks like [0,0,...,0,1]
274  }
275 }
std::vector< boost::shared_ptr< FittingFunc > > _funcArray
template<class FittingFunc >
std::vector< double > lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::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.

Definition at line 234 of file LeastSqFitter2d.h.

234  {
235 
236  std::vector<double> out;
237  out.reserve(_nData);
238 
239  for (int i = 0; i < _nData; ++i) {
240  out.push_back(_z[i] - valueAt(_x[i], _y[i]));
241  }
242 
243  return out;
244 }
double valueAt(double x, double y)
Return the value of the best fit function at a given position (x,y)
template<class FittingFunc >
double lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::valueAt ( double  x,
double  y 
)

Return the value of the best fit function at a given position (x,y)

Definition at line 222 of file LeastSqFitter2d.h.

222  {
223  double val = 0;
224 
225  //Sum the values of the different orders to get the value of the fitted function
226  for (int i = 0; i < _nPar; ++i) {
227  val += _par[i] * func2d(x, y, i);
228  }
229  return val;
230 }
int y
double func2d(double x, double y, int term)
bool val
double x

Member Data Documentation

template<class FittingFunc >
std::vector<boost::shared_ptr<FittingFunc> > lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::_funcArray
private

Definition at line 100 of file LeastSqFitter2d.h.

template<class FittingFunc >
int lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::_nData
private

Definition at line 95 of file LeastSqFitter2d.h.

template<class FittingFunc >
int lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::_nPar
private

Definition at line 94 of file LeastSqFitter2d.h.

template<class FittingFunc >
int lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::_order
private

Definition at line 93 of file LeastSqFitter2d.h.

template<class FittingFunc >
Eigen::VectorXd lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::_par
private

Definition at line 98 of file LeastSqFitter2d.h.

template<class FittingFunc >
std::vector<double> lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::_s
private

Definition at line 92 of file LeastSqFitter2d.h.

template<class FittingFunc >
Eigen::JacobiSVD<Eigen::MatrixXd> lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::_svd
private

Definition at line 97 of file LeastSqFitter2d.h.

template<class FittingFunc >
std::vector<double> lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::_x
private

Definition at line 92 of file LeastSqFitter2d.h.

template<class FittingFunc >
std::vector<double> lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::_y
private

Definition at line 92 of file LeastSqFitter2d.h.

template<class FittingFunc >
std::vector<double> lsst::meas.astrom.sip::LeastSqFitter2d< FittingFunc >::_z
private

Definition at line 92 of file LeastSqFitter2d.h.


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