LSSTApplications  15.0+21,16.0+1,16.0+10,16.0+3,16.0+4,16.0-1-g2115a9e+4,16.0-1-g4515a79+8,16.0-1-g7bb14cc,16.0-1-g80120d7+6,16.0-1-g98efed3+6,16.0-1-gb7f560d+3,16.0-18-g7a076d417,16.0-2-g2ed7261+3,16.0-2-g311bfd2,16.0-2-g568a347+5,16.0-2-g7adb079,16.0-2-gd4c87cb+5,16.0-3-g099ede0,16.0-3-g150e024+5,16.0-3-g1f513a6+2,16.0-3-g958ce35,16.0-3-gc6a11d1,16.0-4-g84f75fb+7,16.0-4-gcfd1396+6,16.0-4-gde8cee2,16.0-5-g7bc0afb+5,16.0-5-g81851deb,16.0-5-g82b7855+1,16.0-5-gd32631f,16.0-5-gf14cb0b,16.0-6-g2dd73041+6,16.0-6-gcf12234+1,16.0-7-g95fb7bf+2,16.0-7-gc37dbc2+6,w.2018.28
LSSTDataManagementBasePackage
angle.cc
Go to the documentation of this file.
1 /*
2  * Developed for the LSST Data Management System.
3  * This product includes software developed by the LSST Project
4  * (https://www.lsst.org).
5  * See the COPYRIGHT file at the top-level directory of this distribution
6  * for details of code ownership.
7  *
8  * This program is free software: you can redistribute it and/or modify
9  * it under the terms of the GNU General Public License as published by
10  * the Free Software Foundation, either version 3 of the License, or
11  * (at your option) any later version.
12  *
13  * This program is distributed in the hope that it will be useful,
14  * but WITHOUT ANY WARRANTY; without even the implied warranty of
15  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  * GNU General Public License for more details.
17  *
18  * You should have received a copy of the GNU General Public License
19  * along with this program. If not, see <https://www.gnu.org/licenses/>.
20  */
21 
22 #include "pybind11/pybind11.h"
23 
24 #include "lsst/utils/python.h"
25 
26 #include "lsst/geom/Angle.h"
27 
28 namespace py = pybind11;
29 
30 namespace lsst {
31 namespace geom {
32 namespace {
33 
34 using PyAngle = py::class_<Angle>;
35 using PyAngleUnit = py::class_<AngleUnit>;
36 
37 template <typename OtherT>
38 void declareAngleComparisonOperators(PyAngle& cls) {
39  cls.def("__eq__", [](Angle const& self, OtherT const& other) { return self == other; },
40  py::is_operator());
41  cls.def("__ne__", [](Angle const& self, OtherT const& other) { return self != other; },
42  py::is_operator());
43  cls.def("__le__", [](Angle const& self, OtherT const& other) { return self <= other; },
44  py::is_operator());
45  cls.def("__ge__", [](Angle const& self, OtherT const& other) { return self >= other; },
46  py::is_operator());
47  cls.def("__lt__", [](Angle const& self, OtherT const& other) { return self < other; }, py::is_operator());
48  cls.def("__gt__", [](Angle const& self, OtherT const& other) { return self > other; }, py::is_operator());
49 }
50 
52  py::module mod("angle");
53 
54  /* AngleUnit */
55 
56  PyAngleUnit clsAngleUnit(mod, "AngleUnit");
57 
58  clsAngleUnit.def("__eq__", [](AngleUnit const& self, AngleUnit const& other) { return self == other; },
59  py::is_operator());
60  clsAngleUnit.def("__ne__", [](AngleUnit const& self, AngleUnit const& other) { return !(self == other); },
61  py::is_operator());
62  clsAngleUnit.def("_mul", [](AngleUnit const& self, double other) { return other * self; },
63  py::is_operator());
64  clsAngleUnit.def("_rmul", [](AngleUnit const& self, double other) { return other * self; },
65  py::is_operator());
66  mod.attr("radians") = py::cast(radians);
67  mod.attr("degrees") = py::cast(degrees);
68  mod.attr("hours") = py::cast(hours);
69  mod.attr("arcminutes") = py::cast(arcminutes);
70  mod.attr("arcseconds") = py::cast(arcseconds);
71  mod.attr("milliarcseconds") = py::cast(milliarcseconds);
72 
73  /* Angle */
74 
75  PyAngle clsAngle(mod, "Angle");
76 
77  clsAngle.def(py::init<double, AngleUnit>(), py::arg("val"), py::arg("units") = radians);
78  clsAngle.def(py::init<>());
79 
80  declareAngleComparisonOperators<Angle>(clsAngle);
81  declareAngleComparisonOperators<double>(clsAngle);
82  declareAngleComparisonOperators<int>(clsAngle);
83 
84  clsAngle.def("__mul__", [](Angle const& self, double other) { return self * other; }, py::is_operator());
85  clsAngle.def("__mul__", [](Angle const& self, int other) { return self * other; }, py::is_operator());
86  clsAngle.def("__rmul__", [](Angle const& self, double other) { return self * other; }, py::is_operator());
87  clsAngle.def("__rmul__", [](Angle const& self, int other) { return self * other; }, py::is_operator());
88  clsAngle.def("__imul__", [](Angle& self, double other) { return self *= other; });
89  clsAngle.def("__imul__", [](Angle& self, int other) { return self *= other; });
90  clsAngle.def("__add__", [](Angle const& self, Angle const& other) { return self + other; },
91  py::is_operator());
92  clsAngle.def("__sub__", [](Angle const& self, Angle const& other) { return self - other; },
93  py::is_operator());
94  clsAngle.def("__neg__", [](Angle const& self) { return -self; }, py::is_operator());
95  clsAngle.def("__iadd__", [](Angle& self, Angle const& other) { return self += other; });
96  clsAngle.def("__isub__", [](Angle& self, Angle const& other) { return self -= other; });
97  clsAngle.def("__truediv__", [](Angle const& self, double other) { return self / other; },
98  py::is_operator());
99  // Without an explicit wrapper, Python lets Angle / Angle -> Angle
100  clsAngle.def("__truediv__", [](Angle const& self, Angle const& other) {
101  throw py::type_error("unsupported operand type(s) for /: 'Angle' and 'Angle'");
102  });
103 
104  clsAngle.def("__float__", &Angle::operator double);
105  clsAngle.def("__abs__", [](Angle const& self) { return std::abs(self.asRadians()) * radians; });
106 
107  clsAngle.def("__reduce__", [clsAngle](Angle const& self) {
108  return py::make_tuple(clsAngle, py::make_tuple(py::cast(self.asRadians())));
109  });
110 
111  utils::python::addOutputOp(clsAngle, "__str__");
112  utils::python::addOutputOp(clsAngle, "__repr__");
113 
114  clsAngle.def("asAngularUnits", &Angle::asAngularUnits);
115  clsAngle.def("asRadians", &Angle::asRadians);
116  clsAngle.def("asDegrees", &Angle::asDegrees);
117  clsAngle.def("asHours", &Angle::asHours);
118  clsAngle.def("asArcminutes", &Angle::asArcminutes);
119  clsAngle.def("asArcseconds", &Angle::asArcseconds);
120  clsAngle.def("asMilliarcseconds", &Angle::asMilliarcseconds);
121  clsAngle.def("wrap", &Angle::wrap);
122  clsAngle.def("wrapCtr", &Angle::wrapCtr);
123  clsAngle.def("wrapNear", &Angle::wrapNear);
124  clsAngle.def("separation", &Angle::separation);
125 
126  /* Non-members */
127 
128  mod.attr("PI") = py::float_(PI);
129  mod.attr("TWOPI") = py::float_(TWOPI);
130  mod.attr("HALFPI") = py::float_(HALFPI);
131  mod.attr("ONE_OVER_PI") = py::float_(ONE_OVER_PI);
132  mod.attr("SQRTPI") = py::float_(SQRTPI);
133  mod.attr("INVSQRTPI") = py::float_(INVSQRTPI);
134  mod.attr("ROOT2") = py::float_(ROOT2);
135 
136  mod.def("degToRad", degToRad);
137  mod.def("radToDeg", radToDeg);
138  mod.def("radToArcsec", radToArcsec);
139  mod.def("radToMas", radToMas);
140  mod.def("arcsecToRad", arcsecToRad);
141  mod.def("masToRad", masToRad);
142  mod.def("isAngle", isAngle<Angle>);
143  mod.def("isAngle", isAngle<double>);
144 
145  return mod.ptr();
146 }
147 
148 } // namespace
149 } // namespace geom
150 } // namespace lsst
constexpr double asAngularUnits(AngleUnit const &units) const noexcept
Return an Angle&#39;s value in the specified units.
Definition: Angle.h:160
AngleUnit constexpr arcminutes
constant with units of arcminutes
Definition: Angle.h:108
Angle abs(Angle const &a)
Definition: Angle.h:106
void addOutputOp(PyClass &cls, std::string const &method)
Add __str__ or __repr__ method implemented by operator<<.
Definition: python.h:82
constexpr double asRadians() const noexcept
Return an Angle&#39;s value in radians.
Definition: Angle.h:163
constexpr double asDegrees() const noexcept
Return an Angle&#39;s value in degrees.
Definition: Angle.h:166
double constexpr TWOPI
Definition: Angle.h:40
AngleUnit constexpr hours
constant with units of hours
Definition: Angle.h:107
constexpr double radToDeg(double x) noexcept
Definition: Angle.h:52
Angle separation(Angle const &other) const noexcept
The signed difference between two Angles.
Definition: Angle.h:443
AngleUnit constexpr radians
constant with units of radians
Definition: Angle.h:105
Angle wrapCtr() const noexcept
Wrap this angle to the range [-π, π).
Definition: Angle.h:406
constexpr double radToArcsec(double x) noexcept
Definition: Angle.h:53
double constexpr ONE_OVER_PI
Definition: Angle.h:42
PYBIND11_PLUGIN(cameraSys)
Definition: cameraSys.cc:62
AngleUnit constexpr degrees
constant with units of degrees
Definition: Angle.h:106
A base class for image defects.
Definition: cameraGeom.dox:3
constexpr double asHours() const noexcept
Return an Angle&#39;s value in hours.
Definition: Angle.h:169
AngleUnit constexpr arcseconds
constant with units of arcseconds
Definition: Angle.h:109
constexpr double degToRad(double x) noexcept
Definition: Angle.h:51
constexpr double asArcseconds() const noexcept
Return an Angle&#39;s value in arcseconds.
Definition: Angle.h:175
Angle wrapNear(Angle const &refAng) const noexcept
Wrap this angle to a value x such that -π ≤ x - refAng ≤ π, approximately.
Definition: Angle.h:425
constexpr double asMilliarcseconds() const noexcept
Return an Angle&#39;s value in milliarcseconds.
Definition: Angle.h:178
double const INVSQRTPI
Definition: Angle.h:45
double constexpr HALFPI
Definition: Angle.h:41
constexpr double radToMas(double x) noexcept
Definition: Angle.h:54
double const SQRTPI
Definition: Angle.h:44
ItemVariant const * other
Definition: Schema.cc:55
constexpr double asArcminutes() const noexcept
Return an Angle&#39;s value in arcminutes.
Definition: Angle.h:172
lsst::geom::Angle Angle
Definition: misc.h:33
Angle wrap() const noexcept
Wrap this angle to the range [0, 2π).
Definition: Angle.h:397
AngleUnit constexpr milliarcseconds
constant with units of milliarcseconds
Definition: Angle.h:112
constexpr double masToRad(double x) noexcept
Definition: Angle.h:56
constexpr double arcsecToRad(double x) noexcept
Definition: Angle.h:55
double constexpr ROOT2
Definition: Angle.h:46
double constexpr PI
The ratio of a circle&#39;s circumference to diameter.
Definition: Angle.h:39