LSST Applications  21.0.0-172-gfb10e10a+18fedfabac,22.0.0+297cba6710,22.0.0+80564b0ff1,22.0.0+8d77f4f51a,22.0.0+a28f4c53b1,22.0.0+dcf3732eb2,22.0.1-1-g7d6de66+2a20fdde0d,22.0.1-1-g8e32f31+297cba6710,22.0.1-1-geca5380+7fa3b7d9b6,22.0.1-12-g44dc1dc+2a20fdde0d,22.0.1-15-g6a90155+515f58c32b,22.0.1-16-g9282f48+790f5f2caa,22.0.1-2-g92698f7+dcf3732eb2,22.0.1-2-ga9b0f51+7fa3b7d9b6,22.0.1-2-gd1925c9+bf4f0e694f,22.0.1-24-g1ad7a390+a9625a72a8,22.0.1-25-g5bf6245+3ad8ecd50b,22.0.1-25-gb120d7b+8b5510f75f,22.0.1-27-g97737f7+2a20fdde0d,22.0.1-32-gf62ce7b1+aa4237961e,22.0.1-4-g0b3f228+2a20fdde0d,22.0.1-4-g243d05b+871c1b8305,22.0.1-4-g3a563be+32dcf1063f,22.0.1-4-g44f2e3d+9e4ab0f4fa,22.0.1-42-gca6935d93+ba5e5ca3eb,22.0.1-5-g15c806e+85460ae5f3,22.0.1-5-g58711c4+611d128589,22.0.1-5-g75bb458+99c117b92f,22.0.1-6-g1c63a23+7fa3b7d9b6,22.0.1-6-g50866e6+84ff5a128b,22.0.1-6-g8d3140d+720564cf76,22.0.1-6-gd805d02+cc5644f571,22.0.1-8-ge5750ce+85460ae5f3,master-g6e05de7fdc+babf819c66,master-g99da0e417a+8d77f4f51a,w.2021.48
LSST Data Management Base Package
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 
31 namespace lsst {
32 namespace meas {
33 namespace astrom {
34 
35 class SipForwardTransform;
36 class SipReverseTransform;
37 class ScaledPolynomialTransform;
38 
46 public:
51 
55  static PolynomialTransform convert(SipForwardTransform const& other);
56 
60  static PolynomialTransform convert(SipReverseTransform const& other);
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 
135 private:
137 
141  friend class SipForwardTransform;
142  friend class SipReverseTransform;
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 
158 public:
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 
227 private:
229  PolynomialTransform _poly;
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 > getYCoeffs() const
2-D polynomial coefficients that compute the output x coordinate.
void swap(PolynomialTransform &other)
Lightweight swap.
ndarray::Array< double const, 2, 2 > getXCoeffs() const
2-D polynomial coefficients that compute the output x coordinate.
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 & operator=(ScaledPolynomialTransform &&other)=default
ScaledPolynomialTransform & operator=(ScaledPolynomialTransform const &other)=default
ScaledPolynomialTransform(ScaledPolynomialTransform const &other)=default
geom::AffineTransform linearize(geom::Point2D const &in) const
Return an approximate affine transform at the given point.
geom::AffineTransform const & getOutputScalingInverse() const
Return the affine transform applied to points after the polynomial transform.
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.
geom::AffineTransform const & getInputScaling() const
Return the first affine transform applied to input points.
PolynomialTransform const & getPoly() const
Return the polynomial transform applied after the input scaling.
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...
Definition: SipTransform.h:136
A transform that maps intermediate world coordinates to pixel coordinates according to the SIP conven...
Definition: SipTransform.h:246
Low-level polynomials (including special polynomials) in C++.
Definition: Basis1d.h:26
PolynomialTransform compose(geom::AffineTransform const &t1, PolynomialTransform const &t2)
Return a PolynomialTransform that is equivalent to the composition t1(t2())
A base class for image defects.
table::Key< int > order