LSST Applications g0da5cf3356+25b44625d0,g17e5ecfddb+50a5ac4092,g1c76d35bf8+585f0f68a2,g295839609d+8ef6456700,g2e2c1a68ba+cc1f6f037e,g38293774b4+62d12e78cb,g3b44f30a73+2891c76795,g48ccf36440+885b902d19,g4b2f1765b6+0c565e8f25,g5320a0a9f6+bd4bf1dc76,g56364267ca+403c24672b,g56b687f8c9+585f0f68a2,g5c4744a4d9+78cd207961,g5ffd174ac0+bd4bf1dc76,g6075d09f38+3075de592a,g667d525e37+cacede5508,g6f3e93b5a3+da81c812ee,g71f27ac40c+cacede5508,g7212e027e3+eb621d73aa,g774830318a+18d2b9fa6c,g7985c39107+62d12e78cb,g79ca90bc5c+fa2cc03294,g881bdbfe6c+cacede5508,g91fc1fa0cf+82a115f028,g961520b1fb+2534687f64,g96f01af41f+f2060f23b6,g9ca82378b8+cacede5508,g9d27549199+78cd207961,gb065e2a02a+ad48cbcda4,gb1df4690d6+585f0f68a2,gb35d6563ee+62d12e78cb,gbc3249ced9+bd4bf1dc76,gbec6a3398f+bd4bf1dc76,gd01420fc67+bd4bf1dc76,gd59336e7c4+c7bb92e648,gf46e8334de+81c9a61069,gfed783d017+bd4bf1dc76,v25.0.1.rc3
LSST Data Management Base Package
Loading...
Searching...
No Matches
_angle.cc
Go to the documentation of this file.
1/*
2 * LSST Data Management System
3 * See COPYRIGHT file at the top of the source tree.
4 *
5 * This product includes software developed by the
6 * LSST Project (http://www.lsst.org/).
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 LSST License Statement and
19 * the GNU General Public License along with this program. If not,
20 * see <https://www.lsstcorp.org/LegalNotices/>.
21 */
22#include "pybind11/pybind11.h"
23
24#include "lsst/sphgeom/python.h"
25
26#include "lsst/sphgeom/Angle.h"
28
29namespace py = pybind11;
30using namespace pybind11::literals;
31
32namespace lsst {
33namespace sphgeom {
34
35template <>
36void defineClass(py::class_<Angle> &cls) {
37 cls.def_static("nan", &Angle::nan);
38 cls.def_static("fromDegrees", &Angle::fromDegrees);
39 cls.def_static("fromRadians", &Angle::fromRadians);
40
41 cls.def(py::init<>());
42 cls.def(py::init<double>(), "radians"_a);
43 cls.def(py::init<Angle>(), "angle"_a);
44 // Construct an Angle from a NormalizedAngle, enabling implicit
45 // conversion from NormalizedAngle to Angle in python via
46 // py::implicitly_convertible
47 cls.def(py::init(
48 [](NormalizedAngle &a) {
49 return new Angle(a.asRadians());
50 }),
51 "normalizedAngle"_a);
52
53 cls.def("__eq__", &Angle::operator==, py::is_operator());
54 cls.def("__ne__", &Angle::operator!=, py::is_operator());
55 cls.def("__lt__", &Angle::operator<, py::is_operator());
56 cls.def("__gt__", &Angle::operator>, py::is_operator());
57 cls.def("__le__", &Angle::operator<=, py::is_operator());
58 cls.def("__ge__", &Angle::operator>=, py::is_operator());
59
60 cls.def("__neg__", (Angle(Angle::*)() const) & Angle::operator-);
61 cls.def("__add__", &Angle::operator+, py::is_operator());
62 cls.def("__sub__",
63 (Angle(Angle::*)(Angle const &) const) & Angle::operator-,
64 py::is_operator());
65 cls.def("__mul__", &Angle::operator*, py::is_operator());
66 cls.def("__rmul__", &Angle::operator*, py::is_operator());
67 cls.def("__truediv__", (Angle(Angle::*)(double) const) & Angle::operator/,
68 py::is_operator());
69 cls.def("__truediv__",
70 (double (Angle::*)(Angle const &) const) & Angle::operator/,
71 py::is_operator());
72
73 cls.def("__iadd__", &Angle::operator+=);
74 cls.def("__isub__", &Angle::operator-=);
75 cls.def("__imul__", &Angle::operator*=);
76 cls.def("__itruediv__", &Angle::operator/=);
77
78 cls.def("asDegrees", &Angle::asDegrees);
79 cls.def("asRadians", &Angle::asRadians);
80 cls.def("isNormalized", &Angle::isNormalized);
81 cls.def("isNan", &Angle::isNan);
82
83 cls.def("__str__", [](Angle const &self) {
84 return py::str("{!s}").format(self.asRadians());
85 });
86 cls.def("__repr__", [](Angle const &self) {
87 return py::str("Angle({!r})").format(self.asRadians());
88 });
89
90 cls.def("__reduce__", [cls](Angle const &self) {
91 return py::make_tuple(cls, py::make_tuple(self.asRadians()));
92 });
93}
94
95} // sphgeom
96} // lsst
This file declares a class for representing normalized angles.
table::Key< int > a
Angle represents an angle in radians.
Definition: Angle.h:43
static Angle nan()
Definition: Angle.h:45
static Angle fromDegrees(double a)
Definition: Angle.h:49
bool isNan() const
isNan returns true if the angle value is NaN.
Definition: Angle.h:91
bool isNormalized() const
isNormalized returns true if this angle lies in the range [0, 2π).
Definition: Angle.h:88
static Angle fromRadians(double a)
Definition: Angle.h:51
double asDegrees() const
asDegrees returns the value of this angle in units of degrees.
Definition: Angle.h:82
double asRadians() const
asRadians returns the value of this angle in units of radians.
Definition: Angle.h:85
NormalizedAngle is an angle that lies in the range [0, 2π), with one exception - a NormalizedAngle ca...
void defineClass(Pybind11Class &cls)
This file declares a class for representing angles.