File LeastSqFitter1d.h¶
-
namespace
lsst
Class for a simple mapping implementing a generic AstrometryTransform.
Remove all non-astronomical counts from the Chunk Exposure’s pixels.
Forward declarations for lsst::utils::Cache
For details on the Cache class, see the Cache.h file.
It uses a template rather than a pointer so that the derived classes can use the specifics of the transform. The class simplePolyMapping overloads a few routines.
A base class for image defects
Numeric constants used by the Integrate.h integrator routines.
Compute Image Statistics
- Note
Gauss-Kronrod-Patterson quadrature coefficients for use in quadpack routine qng. These coefficients were calculated with 101 decimal digit arithmetic by L. W. Fullerton, Bell Labs, Nov 1981.
- Note
The Statistics class itself can only handle lsst::afw::image::MaskedImage() types. The philosophy has been to handle other types by making them look like lsst::afw::image::MaskedImage() and reusing that code. Users should have no need to instantiate a Statistics object directly, but should use the overloaded makeStatistics() factory functions.
-
namespace
meas
-
namespace
astrom
-
namespace
sip
-
template<class
FittingFunc
>
classLeastSqFitter1d
- #include <LeastSqFitter1d.h>
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\)
- See
- Template Parameters
FittingFunc
: The 1d function to fit in both dimensions. Must inherit from lsst::afw::math::Function1
- Parameters
x
: Ordinate of points to fity
: Co-ordinate of pionts to fits
: 1 \(\sigma\) uncertainties in zorder
: Polynomial order to fit
Public 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)
- Template Parameters
FittingFunc
: The type of function to fit. This function extends the base class of lsst::afw::math::Function1
- Parameters
x
: vector of x positions of datay
: vector of y positions of datas
: Vector of measured uncertainties in the values of zorder
: Order of 2d function to fit
-
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.
-
FittingFunc
getBestFitFunction
() Return the best fit polynomial as a lsst::afw::math::Function1 object.
-
double
valueAt
(double x) Calculate the value of the function at a given point.
-
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.
-
double
getChiSq
() Return a measure of the goodness of fit.
\[ \chi_r^2 = \sum \left( \frac{y_i - f(x_i)}{s} \right)^2 \].
-
double
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.
-
template<class
-
namespace
-
namespace