26 #ifndef LEAST_SQ_FITTER_2D
27 #define LEAST_SQ_FITTER_2D
33 #include "boost/shared_ptr.hpp"
41 namespace except = lsst::pex::exceptions;
42 namespace pexLogging = lsst::pex::logging;
73 LeastSqFitter2d(
const std::vector<double> &
x,
const std::vector<double> &
y,
const std::vector<double> &z,
74 const std::vector<double> &s,
int order);
78 double valueAt(
double x,
double y);
87 Eigen::MatrixXd
expandParams(Eigen::VectorXd
const & input)
const;
89 double func2d(
double x,
double y,
int term);
90 double func1d(
double value,
int exponent);
97 Eigen::JacobiSVD<Eigen::MatrixXd>
_svd;
109 namespace sip = lsst::meas::astrom::sip;
110 namespace math = lsst::afw::math;
122 const std::vector<double> &
y,
const std::vector<double> &z,
const std::vector<double> &s,
int order) :
123 _x(x), _y(y), _z(z), _s(s), _order(order), _nPar(0), _par(1) {
127 for (
int i = 0; i < order; ++i) {
133 if (
_nData != static_cast<int>(
_y.size())) {
134 throw LSST_EXCEPT(except::RuntimeError,
"x and y vectors of different lengths");
136 if (
_nData != static_cast<int>(
_s.size())) {
137 throw LSST_EXCEPT(except::RuntimeError,
"x and s vectors of different lengths");
139 if (
_nData != static_cast<int>(
_z.size())) {
140 throw LSST_EXCEPT(except::RuntimeError,
"x and z vectors of different lengths");
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.";
151 throw LSST_EXCEPT(except::RuntimeError,
"Fewer data points than parameters");
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) {
164 _svd.compute(design, Eigen::ComputeThinU | Eigen::ComputeThinV);
179 return expandParams(_par);
183 template<
class FittingFunc>
185 Eigen::MatrixXd out = Eigen::MatrixXd::Zero(_order, _order);
187 for (
int i = 0; i < _order; ++i) {
188 for (
int j = 0; j < _order - i; ++j) {
189 out(i, j) = input[count++];
200 for (
int i = 0; i < _nData; ++i) {
201 double val = _z[i] - valueAt(_x[i], _y[i]);
203 chisq += pow(val, 2);
216 return getChiSq()/(double) (_nData - _nPar);
226 for (
int i = 0; i < _nPar; ++i) {
227 val += _par[i] * func2d(x, y, i);
236 std::vector<double> out;
239 for (
int i = 0; i < _nData; ++i) {
240 out.push_back(_z[i] - valueAt(_x[i], _y[i]));
248 template<
class FittingFunc>
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()
256 return expandParams(variance.sqrt().matrix());
262 _funcArray.reserve(_order);
264 std::vector<double> coeff;
265 coeff.reserve( _order);
268 for (
int i = 0; i < _order; ++i) {
269 boost::shared_ptr<FittingFunc> p(
new FittingFunc(coeff));
270 _funcArray.push_back(p);
285 for (
int i = 0; i<term; ++i) {
286 yexp = (yexp + 1) % (_order - xexp);
292 double xcomp = func1d(x, xexp);
293 double ycomp = func1d(y, yexp);
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);
303 return (*_funcArray[exponent])(value);
Fit an lsst::afw::math::Function1 object to a set of data points in two dimensions.
Eigen::MatrixXd expandParams(Eigen::VectorXd const &input) const
Turn a flattened parameter-like vector into a triangular matrix.
double func2d(double x, double y, int term)
double getReducedChiSq()
Return a measure of the goodness of fit. Where is the number of data points, and is the number of ...
definition of the Trace messaging facilities
std::vector< double > residuals()
Eigen::JacobiSVD< Eigen::MatrixXd > _svd
Define a collection of useful Functions.
double getChiSq()
Return a measure of the goodness of fit. .
double valueAt(double x, double y)
Return the value of the best fit function at a given position (x,y)
Eigen::MatrixXd getParams()
#define LSST_EXCEPT(type,...)
double func1d(double value, int exponent)
Eigen::MatrixXd getErrors()
Companion function to getParams(). Returns uncertainties in the parameters as a matrix.
LeastSqFitter2d(const std::vector< double > &x, const std::vector< double > &y, const std::vector< double > &z, const std::vector< double > &s, int order)
Eigen::JacobiSVD< Eigen::MatrixXd > _svd
std::vector< boost::shared_ptr< FittingFunc > > _funcArray