LSST Applications g0f08755f38+82efc23009,g12f32b3c4e+e7bdf1200e,g1653933729+a8ce1bb630,g1a0ca8cf93+50eff2b06f,g28da252d5a+52db39f6a5,g2bbee38e9b+37c5a29d61,g2bc492864f+37c5a29d61,g2cdde0e794+c05ff076ad,g3156d2b45e+41e33cbcdc,g347aa1857d+37c5a29d61,g35bb328faa+a8ce1bb630,g3a166c0a6a+37c5a29d61,g3e281a1b8c+fb992f5633,g414038480c+7f03dfc1b0,g41af890bb2+11b950c980,g5fbc88fb19+17cd334064,g6b1c1869cb+12dd639c9a,g781aacb6e4+a8ce1bb630,g80478fca09+72e9651da0,g82479be7b0+04c31367b4,g858d7b2824+82efc23009,g9125e01d80+a8ce1bb630,g9726552aa6+8047e3811d,ga5288a1d22+e532dc0a0b,gae0086650b+a8ce1bb630,gb58c049af0+d64f4d3760,gc28159a63d+37c5a29d61,gcf0d15dbbd+2acd6d4d48,gd7358e8bfb+778a810b6e,gda3e153d99+82efc23009,gda6a2b7d83+2acd6d4d48,gdaeeff99f8+1711a396fd,ge2409df99d+6b12de1076,ge79ae78c31+37c5a29d61,gf0baf85859+d0a5978c5a,gf3967379c6+4954f8c433,gfb92a5be7c+82efc23009,gfec2e1e490+2aaed99252,w.2024.46
LSST Data Management Base Package
Loading...
Searching...
No Matches
PolynomialTransform.h
Go to the documentation of this file.
1// -*- LSST-C++ -*-
2
3/*
4 * LSST Data Management System
5 * Copyright 2016 LSST/AURA
6 *
7 * This product includes software developed by the
8 * LSST Project (http://www.lsst.org/).
9 *
10 * This program is free software: you can redistribute it and/or modify
11 * it under the terms of the GNU General Public License as published by
12 * the Free Software Foundation, either version 3 of the License, or
13 * (at your option) any later version.
14 *
15 * This program is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the LSST License Statement and
21 * the GNU General Public License along with this program. If not,
22 * see <http://www.lsstcorp.org/LegalNotices/>.
23 */
24#ifndef LSST_MEAS_ASTROM_PolynomialTransform_INCLUDED
25#define LSST_MEAS_ASTROM_PolynomialTransform_INCLUDED
26
27#include "ndarray/eigen.h"
28#include "lsst/geom/Point.h"
30
31namespace lsst {
32namespace meas {
33namespace astrom {
34
35class SipForwardTransform;
36class SipReverseTransform;
37class ScaledPolynomialTransform;
38
46public:
51
56
61
72 PolynomialTransform(ndarray::Array<double const, 2, 0> const& xCoeffs,
73 ndarray::Array<double const, 2, 0> const& yCoeffs);
74
81
88
95
102
104 void swap(PolynomialTransform& other);
105
107 int getOrder() const { return _xCoeffs.getSize<0>() - 1; }
108
115 ndarray::Array<double const, 2, 2> getXCoeffs() const { return _xCoeffs.shallow(); }
116
123 ndarray::Array<double const, 2, 2> getYCoeffs() const { return _yCoeffs.shallow(); }
124
129
133 geom::Point2D operator()(geom::Point2D const& in) const;
134
135private:
137
144
145 ndarray::Array<double, 2, 2> _xCoeffs;
146 ndarray::Array<double, 2, 2> _yCoeffs;
147 mutable Eigen::VectorXd _u; // workspace for operator() and linearize
148 mutable Eigen::VectorXd _v;
149};
150
158public:
166
174 static ScaledPolynomialTransform convert(SipForwardTransform const& sipForward);
175
183 static ScaledPolynomialTransform convert(SipReverseTransform const& sipReverse);
184
196 geom::AffineTransform const& outputScalingInverse);
197
199
201
203
205
206 void swap(ScaledPolynomialTransform& other);
207
209 PolynomialTransform const& getPoly() const { return _poly; }
210
212 geom::AffineTransform const& getInputScaling() const { return _inputScaling; }
213
215 geom::AffineTransform const& getOutputScalingInverse() const { return _outputScalingInverse; }
216
221
225 geom::Point2D operator()(geom::Point2D const& in) const;
226
227private:
230 geom::AffineTransform _inputScaling;
231 geom::AffineTransform _outputScalingInverse;
232};
233
241
249
250} // namespace astrom
251} // namespace meas
252} // namespace lsst
253
254#endif // !LSST_MEAS_ASTROM_PolynomialTransform_INCLUDED
An affine coordinate transformation consisting of a linear transformation and an offset.
A 2-d coordinate transform represented by a pair of standard polynomials (one for each coordinate).
geom::Point2D operator()(geom::Point2D const &in) const
Apply the transform to a point.
PolynomialTransform & operator=(PolynomialTransform const &other)
Copy assignment.
friend PolynomialTransform compose(geom::AffineTransform const &t1, PolynomialTransform const &t2)
Return a PolynomialTransform that is equivalent to the composition t1(t2())
static PolynomialTransform convert(ScaledPolynomialTransform const &other)
Convert a ScaledPolynomialTransform to an equivalent PolynomialTransform.
geom::AffineTransform linearize(geom::Point2D const &in) const
Return an approximate affine transform at the given point.
PolynomialTransform(ndarray::Array< double const, 2, 0 > const &xCoeffs, ndarray::Array< double const, 2, 0 > const &yCoeffs)
Construct a new transform from existing coefficient arrays.
ndarray::Array< double const, 2, 2 > getXCoeffs() const
2-D polynomial coefficients that compute the output x coordinate.
ndarray::Array< double const, 2, 2 > getYCoeffs() const
2-D polynomial coefficients that compute the output x coordinate.
void swap(PolynomialTransform &other)
Lightweight swap.
int getOrder() const
Return the order of the polynomials.
A fitter class for scaled polynomial transforms.
A 2-d coordinate transform represented by a lazy composition of an AffineTransform,...
ScaledPolynomialTransform(ScaledPolynomialTransform &&other)=default
ScaledPolynomialTransform(ScaledPolynomialTransform const &other)=default
ScaledPolynomialTransform & operator=(ScaledPolynomialTransform const &other)=default
geom::AffineTransform const & getOutputScalingInverse() const
Return the affine transform applied to points after the polynomial transform.
ScaledPolynomialTransform & operator=(ScaledPolynomialTransform &&other)=default
geom::AffineTransform linearize(geom::Point2D const &in) const
Return an approximate affine transform at the given point.
geom::Point2D operator()(geom::Point2D const &in) const
Apply the transform to a point.
void swap(ScaledPolynomialTransform &other)
static ScaledPolynomialTransform convert(PolynomialTransform const &poly)
Convert a PolynomialTransform to an equivalent ScaledPolynomialTransform.
PolynomialTransform const & getPoly() const
Return the polynomial transform applied after the input scaling.
geom::AffineTransform const & getInputScaling() const
Return the first affine transform applied to input points.
ScaledPolynomialTransform(PolynomialTransform const &poly, geom::AffineTransform const &inputScaling, geom::AffineTransform const &outputScalingInverse)
Construct a new ScaledPolynomialTransform from its constituents.
A transform that maps pixel coordinates to intermediate world coordinates according to the SIP conven...
A transform that maps intermediate world coordinates to pixel coordinates according to the SIP conven...
Low-level polynomials (including special polynomials) in C++.
PolynomialTransform compose(geom::AffineTransform const &t1, PolynomialTransform const &t2)
Return a PolynomialTransform that is equivalent to the composition t1(t2())
table::Key< int > order