26 #ifndef LEAST_SQ_FITTER_1D 27 #define LEAST_SQ_FITTER_1D 36 #include "lsst/pex/exceptions/Runtime.h" 37 #include "lsst/afw/math/FunctionLibrary.h" 39 namespace except = lsst::pex::exceptions;
69 LeastSqFitter1d(
const std::vector<double> &x,
const std::vector<double> &y,
70 const std::vector<double> &s,
int order);
85 double func1d(
double value,
int exponent);
87 std::vector<double> _x, _y, _s;
91 Eigen::JacobiSVD<Eigen::MatrixXd> _svd;
94 std::vector<std::shared_ptr<FittingFunc> > _funcArray;
108 const std::vector<double> &y,
const std::vector<double> &s,
int order) :
109 _x(x), _y(y), _s(s), _order(order) {
112 throw LSST_EXCEPT(except::RuntimeError,
"Fit order must be >= 1");
116 if (_nData != static_cast<int>(_y.size())) {
117 throw LSST_EXCEPT(except::RuntimeError,
"x and y vectors of different lengths");
119 if (_nData != static_cast<int>(_s.size())) {
120 throw LSST_EXCEPT(except::RuntimeError,
"x and s vectors of different lengths");
123 if (_nData < _order) {
124 throw LSST_EXCEPT(except::RuntimeError,
"Fewer data points than parameters");
129 Eigen::MatrixXd design(_nData, _order);
130 Eigen::VectorXd rhs(_nData);
131 for (
int i = 0; i < _nData; ++i) {
132 rhs[i] = y[i] / _s[i];
133 for (
int j = 0; j < _order; ++j) {
134 design(i, j) = func1d(_x[i], j) / _s[i];
137 _svd.compute(design, Eigen::ComputeThinU | Eigen::ComputeThinV);
138 _par = _svd.solve(rhs);
146 Eigen::VectorXd vec = Eigen::VectorXd::Zero(_order);
147 for (
int i = 0; i < _order; ++i) {
156 Eigen::ArrayXd variance(_order);
157 for (
int i = 0; i < _order; ++i) {
158 variance[i] = _svd.matrixV().row(i).dot(
159 (_svd.singularValues().array().inverse().square() * _svd.matrixV().col(i).array()).matrix()
162 return variance.sqrt().matrix();
171 FittingFunc func(_order - 1);
173 for (
int i = 0; i < _order; ++i) {
174 func.setParameter(i, _par(i));
192 std::vector<double> out;
197 for (
int i = 0; i < _nData; ++i) {
198 out.push_back(_y[i] - f(_x[i]));
212 for (
int i = 0; i < _nData; ++i) {
213 double val = _y[i] - f(_x[i]);
215 chisq += pow(val, 2);
228 return getChiSq()/(double) (_nData - _order);
237 _funcArray.reserve(_order);
239 std::vector<double> coeff;
240 coeff.reserve( _order);
242 coeff.push_back(1.0);
243 for (
int i = 0; i < _order; ++i) {
244 std::shared_ptr<FittingFunc> p(
new FittingFunc(coeff));
245 _funcArray.push_back(p);
247 coeff.push_back(1.0);
252 return (*_funcArray[exponent])(value);
double getReducedChiSq()
Return a measure of the goodness of fit.
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.
Eigen::VectorXd getParams()
Return the best fit parameters as an Eigen::Matrix.
Eigen::VectorXd getErrors()
Return the 1 sigma uncertainties in the best fit parameters as an Eigen::Matrix.
double valueAt(double x)
Calculate the value of the function at a given point.
FittingFunc getBestFitFunction()
Return the best fit polynomial as a lsst::afw::math::Function1 object.
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)
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.