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
astrometryTransform.cc
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 /*
3  * This file is part of jointcal.
4  *
5  * Developed for the LSST Data Management System.
6  * This product includes software developed by the LSST Project
7  * (https://www.lsst.org).
8  * See the COPYRIGHT file at the top-level directory of this distribution
9  * for details of code ownership.
10  *
11  * This program is free software: you can redistribute it and/or modify
12  * it under the terms of the GNU General Public License as published by
13  * the Free Software Foundation, either version 3 of the License, or
14  * (at your option) any later version.
15  *
16  * This program is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19  * GNU General Public License for more details.
20  *
21  * You should have received a copy of the GNU General Public License
22  * along with this program. If not, see <https://www.gnu.org/licenses/>.
23  */
24 
25 #include "astshim.h"
26 #include "pybind11/pybind11.h"
27 #include "pybind11/eigen.h"
28 #include "ndarray/pybind11.h"
29 #include "ndarray/eigen.h"
30 #include "Eigen/Core"
31 
32 #include "lsst/utils/python.h"
33 
34 #include "lsst/jointcal/Frame.h"
36 #include "lsst/jointcal/Point.h"
37 
38 namespace py = pybind11;
39 using namespace pybind11::literals;
40 
41 namespace lsst {
42 namespace jointcal {
43 namespace {
44 
45 void declareAstrometryTransform(py::module &mod) {
46  py::class_<AstrometryTransform, std::shared_ptr<AstrometryTransform>> cls(mod, "AstrometryTransform");
47 
48  cls.def("apply",
49  (jointcal::Point(AstrometryTransform::*)(const jointcal::Point &) const) &
50  AstrometryTransform::apply,
51  "inPos"_a);
52  cls.def("apply",
53  (jointcal::Frame(AstrometryTransform::*)(Frame const &, bool) const) & AstrometryTransform::apply,
54  "inputframe"_a, "inscribed"_a);
55  cls.def("getNpar", &AstrometryTransform::getNpar);
56  cls.def("offsetParams", &AstrometryTransform::offsetParams);
57  cls.def("linearApproximation", &AstrometryTransform::linearApproximation);
58  cls.def("computeDerivative", [](AstrometryTransform const &self, Point const &where, double const step) {
59  AstrometryTransformLinear derivative;
60  self.computeDerivative(where, derivative, step);
61  return derivative;
62  });
63 
64  utils::python::addOutputOp(cls, "__str__");
65 }
66 
67 void declareAstrometryTransformIdentity(py::module &mod) {
68  py::class_<AstrometryTransformIdentity, std::shared_ptr<AstrometryTransformIdentity>, AstrometryTransform>
69  cls(mod, "AstrometryTransformIdentity");
70 }
71 
72 void declareAstrometryTransformPolynomial(py::module &mod) {
73  py::class_<AstrometryTransformPolynomial, std::shared_ptr<AstrometryTransformPolynomial>,
74  AstrometryTransform>
75  cls(mod, "AstrometryTransformPolynomial");
76 
77  cls.def(py::init<const unsigned>(), "order"_a);
78  cls.def("getOrder", &AstrometryTransformPolynomial::getOrder);
79  cls.def("getCoefficient", py::overload_cast<std::size_t, std::size_t, std::size_t>(
80  &AstrometryTransformPolynomial::getCoefficient, py::const_));
81  cls.def("setCoefficient", [](AstrometryTransformPolynomial &self, std::size_t powX, std::size_t powY,
82  std::size_t whichCoord,
83  double value) { self.getCoefficient(powX, powY, whichCoord) = value; });
84  cls.def("determinant", &AstrometryTransformPolynomial::determinant);
85  cls.def("getNpar", &AstrometryTransformPolynomial::getNpar);
86  cls.def("toAstMap", &AstrometryTransformPolynomial::toAstMap);
87  cls.def("write", [](AstrometryTransformPolynomial const &self) {
89  self.write(result);
90  return result.str();
91  });
92  cls.def("read", [](AstrometryTransformPolynomial &self, std::string const &str) {
93  std::istringstream istr(str);
94  self.read(istr);
95  });
96 }
97 
98 void declareAstrometryTransformLinear(py::module &mod) {
99  py::class_<AstrometryTransformLinear, std::shared_ptr<AstrometryTransformLinear>,
100  AstrometryTransformPolynomial>
101  cls(mod, "AstrometryTransformLinear");
102 
103  cls.def(py::init<>());
104 }
105 
106 void declareAstrometryTransformLinearShift(py::module &mod) {
107  py::class_<AstrometryTransformLinearShift, std::shared_ptr<AstrometryTransformLinearShift>,
108  AstrometryTransformLinear>
109  cls(mod, "AstrometryTransformLinearShift");
110 }
111 
112 void declareAstrometryTransformLinearRot(py::module &mod) {
113  py::class_<AstrometryTransformLinearRot, std::shared_ptr<AstrometryTransformLinearRot>,
114  AstrometryTransformLinear>
115  cls(mod, "AstrometryTransformLinearRot");
116 }
117 
118 void declareAstrometryTransformLinearScale(py::module &mod) {
119  py::class_<AstrometryTransformLinearScale, std::shared_ptr<AstrometryTransformLinearScale>,
120  AstrometryTransformLinear>
121  cls(mod, "AstrometryTransformLinearScale");
122 }
123 
124 void declareAstrometryTransformSkyWcs(py::module &mod) {
125  py::class_<AstrometryTransformSkyWcs, std::shared_ptr<AstrometryTransformSkyWcs>, AstrometryTransform>
126  cls(mod, "AstrometryTransformSkyWcs");
127  cls.def("getSkyWcs", &AstrometryTransformSkyWcs::getSkyWcs);
128 }
129 
130 void declareBaseTanWcs(py::module &mod) {
131  py::class_<BaseTanWcs, std::shared_ptr<BaseTanWcs>, AstrometryTransform> cls(mod, "BaseTanWcs");
132 }
133 
134 void declareTanPixelToRaDec(py::module &mod) {
135  py::class_<TanPixelToRaDec, std::shared_ptr<TanPixelToRaDec>, AstrometryTransform> cls(mod,
136  "TanPixelToRaDec");
137 }
138 
139 void declareTanRaDecToPixel(py::module &mod) {
140  py::class_<TanRaDecToPixel, std::shared_ptr<TanRaDecToPixel>, AstrometryTransform> cls(mod,
141  "TanRaDecToPixel");
142 }
143 
144 void declareTanSipPixelToRaDec(py::module &mod) {
145  py::class_<TanSipPixelToRaDec, std::shared_ptr<TanSipPixelToRaDec>, BaseTanWcs> cls(mod,
146  "TanSipPixelToRaDec");
147 }
148 
149 PYBIND11_MODULE(astrometryTransform, mod) {
150  py::module::import("astshim.mapping");
151  py::module::import("lsst.jointcal.frame");
152  py::module::import("lsst.jointcal.star");
153  declareAstrometryTransform(mod);
154  declareAstrometryTransformIdentity(mod);
155  declareAstrometryTransformPolynomial(mod);
156  declareAstrometryTransformLinear(mod);
157  declareAstrometryTransformLinearShift(mod);
158  declareAstrometryTransformLinearRot(mod);
159  declareAstrometryTransformLinearScale(mod);
160  declareAstrometryTransformSkyWcs(mod);
161  declareBaseTanWcs(mod);
162  declareTanPixelToRaDec(mod);
163  declareTanRaDecToPixel(mod);
164  declareTanSipPixelToRaDec(mod);
165 
166  // utility functions
167  mod.def("inversePolyTransform", &inversePolyTransform, "forward"_a, "domain"_a, "precision"_a,
168  "maxOrder"_a = 9, "nSteps"_a = 50);
169 }
170 } // namespace
171 } // namespace jointcal
172 } // namespace lsst
py::object result
Definition: _schema.cc:429
int const step
rectangle with sides parallel to axes.
Definition: Frame.h:38
A point in a plane.
Definition: Point.h:37
PYBIND11_MODULE(_cameraGeom, mod)
Definition: _cameraGeom.cc:38
void addOutputOp(PyClass &cls, std::string const &method)
Add __str__ or __repr__ method implemented by operator<<.
Definition: python.h:87
std::shared_ptr< AstrometryTransformPolynomial > inversePolyTransform(AstrometryTransform const &forward, Frame const &domain, double const precision, std::size_t maxOrder=9, std::size_t nSteps=50)
Approximate the inverse by a polynomial, to some precision.
A base class for image defects.