34 #include "boost/format.hpp"
35 #include "boost/shared_ptr.hpp"
43 namespace ex = pex::exceptions;
53 _style(style), _orderX(orderX), _orderY(orderY < 0 ? orderX : orderY), _weighting(weighting) {
55 throw LSST_EXCEPT(lsst::pex::exceptions::InvalidParameterError,
57 "due to a limitation in math::Chebyshev1Function2")
68 template<
typename PixelT>
69 class ApproximateChebyshev :
public Approximate<PixelT> {
72 math::
makeApproximate(std::vector<
double> const &xVec, std::vector<
double> const &yVec,
73 image::MaskedImage<T> const& im, geom::Box2I const& bbox,
76 virtual ~ApproximateChebyshev();
78 math::Chebyshev1Function2<
double>
_poly;
80 ApproximateChebyshev(std::vector<
double> const &xVec, std::vector<
double> const &yVec,
81 image::MaskedImage<
PixelT> const& im, geom::Box2I const& bbox,
82 ApproximateControl const& ctrl);
84 doGetImage(
int orderX,
int orderY) const;
85 virtual
PTR(
image::MaskedImage<typename Approximate<
PixelT>::OutPixelT>)
86 doGetMaskedImage(
int orderX,
int orderY) const;
95 solveMatrix_Eigen(Eigen::MatrixXd &a,
97 Eigen::Map<Eigen::VectorXd> &c
99 Eigen::PartialPivLU<Eigen::MatrixXd> lu(a);
107 template<
typename PixelT>
108 ApproximateChebyshev<PixelT>::ApproximateChebyshev(
109 std::vector<double>
const &xVec,
110 std::vector<double>
const &yVec,
113 ApproximateControl
const& ctrl
115 : Approximate<
PixelT>(xVec, yVec, bbox, ctrl),
116 _poly(math::Chebyshev1Function2<double>(ctrl.getOrderX(), geom::Box2D(bbox)))
120 std::vector<double>
const& coeffs =
_poly.getParameters();
121 assert(std::accumulate(coeffs.begin(), coeffs.end(), 0.0) == 0.0);
124 int const nTerm =
_poly.getNParameters();
136 std::vector<std::vector<double> > termCoeffs(nTerm);
138 for (
int i = 0; i != nTerm; ++i) {
139 termCoeffs[i].reserve(nData);
142 for (
int iy = 0; iy != im.
getHeight(); ++iy) {
143 double const y = yVec[iy];
145 for (
int ix = 0; ix != im.
getWidth(); ++ix) {
146 double const x = xVec[ix];
148 for (
int i = 0; i != nTerm; ++i) {
149 _poly.setParameter(i, 1.0);
150 termCoeffs[i].push_back(
_poly(x, y));
151 _poly.setParameter(i, 0.0);
156 Eigen::MatrixXd A; A.setZero(nTerm, nTerm);
157 Eigen::VectorXd
b; b.setZero(nTerm);
162 for (
int iy = 0; iy != im.
getHeight(); ++iy) {
165 double const val = ptr.image();
166 double const ivar = ctrl.getWeighting() ? 1/ptr.variance() : 1.0;
171 for (
int i = 0; i != nTerm; ++i) {
172 double const c_i = termCoeffs[i][
alpha];
173 double const tmp = c_i*ivar;
176 for (
int j = 0; j < i; ++j) {
177 double const c_j = termCoeffs[j][
alpha];
184 if (ctrl.getWeighting()) {
186 "No valid points to fit. Variance is likely zero. Try weighting=False");
190 "No valid points to fit (even though weighting is False). "
191 "Check that binSize & approxOrderX settings are appropriate for image size.");
195 for (
int j = 0; j != nTerm; ++j) {
196 for (
int i = j + 1; i != nTerm; ++i) {
203 std::vector<double> cvec(nTerm);
204 Eigen::Map<Eigen::VectorXd> c(&cvec[0], nTerm);
206 solveMatrix_Eigen(A, b, c);
208 _poly.setParameters(cvec);
212 template<
typename PixelT>
213 ApproximateChebyshev<PixelT>::~ApproximateChebyshev() {
223 template<
typename PixelT>
229 if (orderX < 0) orderX = Approximate<PixelT>::_ctrl.getOrderX();
230 if (orderY < 0) orderY = Approximate<PixelT>::_ctrl.getOrderY();
233 (orderX == Approximate<PixelT>::_ctrl.getOrderX() &&
234 orderY == Approximate<PixelT>::_ctrl.getOrderY()) ?
_poly :
_poly.truncate(orderX);
239 for (
int iy = 0; iy != im->getHeight(); ++iy) {
243 for (
typename ImageT::x_iterator ptr = im->
row_begin(iy),
244 end = im->
row_end(iy); ptr != end; ++ptr, ++ix) {
261 template<
typename PixelT>
263 ApproximateChebyshev<
PixelT>::doGetMaskedImage(
270 PTR(MImageT) mi(new MImageT(Approximate<
PixelT>::_bbox));
271 PTR(typename MImageT::Image) im = mi->getImage();
273 for (
int iy = 0; iy != im->getHeight(); ++iy) {
277 for (
typename MImageT::Image::x_iterator ptr = im->
row_begin(iy),
278 end = im->
row_end(iy); ptr != end; ++ptr, ++ix) {
293 template<
typename PixelT>
294 PTR(Approximate<PixelT>)
296 std::vector<
double> const &y,
298 geom::Box2I const& bbox,
302 switch (ctrl.getStyle()) {
303 case ApproximateControl::CHEBYSHEV:
306 throw LSST_EXCEPT(lsst::pex::exceptions::InvalidParameterError,
307 str(
boost::format(
"Unknown ApproximationStyle: %d") % ctrl.getStyle()));
315 #define INSTANTIATE(PIXEL_T) \
317 PTR(Approximate<PIXEL_T>) makeApproximate( \
318 std::vector<double> const &x, std::vector<double> const &y, \
319 image::MaskedImage<PIXEL_T> const& im, \
320 geom::Box2I const& bbox, \
321 ApproximateControl const& ctrl)
x_iterator row_begin(int y) const
Return an x_iterator to the start of the image.
2-dimensional weighted sum of Chebyshev polynomials of the first kind.
for(FootprintList::const_iterator ptr=feet->begin(), end=feet->end();ptr!=end;++ptr)
Include files required for standard LSST Exception handling.
int getHeight() const
Return the number of rows in the image.
#define INSTANTIATE(MATCH)
x_iterator row_end(int y) const
Return an x_iterator to the end of the image.
Define a collection of useful Functions.
Approximate values for a MaskedImage.
An integer coordinate rectangle.
table::Key< table::Array< Kernel::Pixel > > image
math::Chebyshev1Function2< double > _poly
Control how to make an approximation.
geom::Box2I const & _bbox
A class to manipulate images, masks, and variance as a single object.
Style
Choose the type of approximation to use.
ApproximateControl(Style style, int orderX, int orderY=-1, bool weighting=true)
ctor
#define LSST_EXCEPT(type,...)
boost::shared_ptr< Approximate< PixelT > > makeApproximate(std::vector< double > const &x, std::vector< double > const &y, image::MaskedImage< PixelT > const &im, geom::Box2I const &bbox, ApproximateControl const &ctrl)
Construct a new Approximate object, inferring the type from the type of the given MaskedImage...
int getWidth() const
Return the number of columns in the image.
afw::table::Key< double > b
Implementation of the Class MaskedImage.
A class to represent a 2-dimensional array of pixels.