26 #ifndef LEAST_SQ_FITTER_1D
27 #define LEAST_SQ_FITTER_1D
39 namespace except = lsst::pex::exceptions;
70 const std::vector<double> &s,
86 double func1d(
double value,
int exponent);
92 Eigen::JacobiSVD<Eigen::MatrixXd>
_svd;
109 const std::vector<double> &
y,
const std::vector<double> &s,
unsigned int order) :
110 _x(x), _y(y), _s(s), _order(order) {
113 throw LSST_EXCEPT(except::RuntimeError,
"Fit order must be >= 1");
117 if (
_nData != static_cast<int>(
_y.size())) {
118 throw LSST_EXCEPT(except::RuntimeError,
"x and y vectors of different lengths");
120 if (
_nData != static_cast<int>(
_s.size())) {
121 throw LSST_EXCEPT(except::RuntimeError,
"x and s vectors of different lengths");
125 throw LSST_EXCEPT(except::RuntimeError,
"Fewer data points than parameters");
131 Eigen::VectorXd rhs(
_nData);
132 for (
int i = 0; i <
_nData; ++i) {
133 rhs[i] = y[i] /
_s[i];
134 for (
int j = 0; j <
_order; ++j) {
138 _svd.compute(design, Eigen::ComputeThinU | Eigen::ComputeThinV);
147 Eigen::VectorXd vec = Eigen::VectorXd::Zero(_order);
148 for (
int i = 0; i < _order; ++i) {
157 Eigen::ArrayXd variance(_order);
158 for (
int i = 0; i < _order; ++i) {
159 variance[i] =
_svd.matrixV().row(i).dot(
160 (
_svd.singularValues().array().inverse().square() *
_svd.matrixV().col(i).array()).matrix()
163 return variance.sqrt().matrix();
172 FittingFunc func(_order - 1);
174 for (
int i = 0; i < _order; ++i) {
175 func.setParameter(i, _par(i));
183 FittingFunc f = getBestFitFunction();
193 std::vector<double> out;
196 FittingFunc f = getBestFitFunction();
198 for (
int i = 0; i < _nData; ++i) {
199 out.push_back(_y[i] - f(_x[i]));
210 FittingFunc f = getBestFitFunction();
213 for (
int i = 0; i < _nData; ++i) {
214 double val = _y[i] - f(_x[i]);
216 chisq += pow(val, 2);
229 return getChiSq()/(double) (_nData - _order);
238 _funcArray.reserve(_order);
240 std::vector<double> coeff;
241 coeff.reserve( _order);
243 coeff.push_back(1.0);
244 for (
int i = 0; i < _order; ++i) {
245 std::shared_ptr<FittingFunc> p(
new FittingFunc(coeff));
246 _funcArray.push_back(p);
248 coeff.push_back(1.0);
253 return (*_funcArray[exponent])(value);
double getReducedChiSq()
Return a measure of the goodness of fit. Where is the number of data points, and is the number of ...
Eigen::JacobiSVD< Eigen::MatrixXd > _svd
std::vector< std::shared_ptr< FittingFunc > > _funcArray
std::vector< double > residuals()
Define a collection of useful Functions.
Eigen::VectorXd getParams()
Return the best fit parameters as an Eigen::Matrix.
LeastSqFitter1d(const std::vector< double > &x, const std::vector< double > &y, const std::vector< double > &s, unsigned int order)
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.
#define LSST_EXCEPT(type,...)
double getChiSq()
Return a measure of the goodness of fit. .
double func1d(double value, int exponent)
Eigen::JacobiSVD< Eigen::MatrixXd > _svd
Fit an lsst::afw::math::Function1 object to a set of data points in one dimension.