angle-axis.hpp
Go to the documentation of this file.
1 /*
2  * Copyright 2014-2020 CNRS INRIA
3  */
4 
5 #ifndef __eigenpy_angle_axis_hpp__
6 #define __eigenpy_angle_axis_hpp__
7 
8 #include "eigenpy/eigenpy.hpp"
9 
10 #include <Eigen/Core>
11 #include <Eigen/Geometry>
12 
13 namespace eigenpy
14 {
15 
16  namespace bp = boost::python;
17 
18  template<typename AngleAxis> class AngleAxisVisitor;
19 
20  template<typename Scalar>
21  struct call< Eigen::AngleAxis<Scalar> >
22  {
23  typedef Eigen::AngleAxis<Scalar> AngleAxis;
24 
25  static inline void expose()
26  {
28  }
29 
30  static inline bool isApprox(const AngleAxis & self, const AngleAxis & other,
31  const Scalar & prec = Eigen::NumTraits<Scalar>::dummy_precision())
32  {
33  return self.isApprox(other,prec);
34  }
35  };
36 
37  template<typename AngleAxis>
38  class AngleAxisVisitor
39  : public bp::def_visitor< AngleAxisVisitor<AngleAxis> >
40  {
41 
42  typedef typename AngleAxis::Scalar Scalar;
43  typedef typename AngleAxis::Vector3 Vector3;
44  typedef typename AngleAxis::Matrix3 Matrix3;
45 
46  typedef typename Eigen::Quaternion<Scalar,0> Quaternion;
47  typedef Eigen::RotationBase<AngleAxis,3> RotationBase;
48 
49  BOOST_PYTHON_FUNCTION_OVERLOADS(isApproxAngleAxis_overload,call<AngleAxis>::isApprox,2,3)
50 
51  public:
52 
53  template<class PyClass>
54  void visit(PyClass& cl) const
55  {
56  cl
57  .def(bp::init<>(bp::arg("self"),"Default constructor"))
58  .def(bp::init<Scalar,Vector3>
59  ((bp::arg("self"),bp::arg("angle"),bp::arg("axis")),
60  "Initialize from angle and axis."))
61  .def(bp::init<Matrix3>
62  ((bp::arg("self"),bp::arg("rotation matrix")),
63  "Initialize from a rotation matrix"))
64  .def(bp::init<Quaternion>((bp::arg("self"),bp::arg("quaternion")),
65  "Initialize from a quaternion."))
66  .def(bp::init<AngleAxis>((bp::arg("self"),bp::arg("copy")),
67  "Copy constructor."))
68 
69  /* --- Properties --- */
70  .add_property("axis",
71  bp::make_function((Vector3 & (AngleAxis::*)())&AngleAxis::axis,
72  bp::return_internal_reference<>()),
73  &AngleAxisVisitor::setAxis,"The rotation axis.")
74  .add_property("angle",
75  (Scalar (AngleAxis::*)()const)&AngleAxis::angle,
76  &AngleAxisVisitor::setAngle,"The rotation angle.")
77 
78  /* --- Methods --- */
79  .def("inverse",&AngleAxis::inverse,
80  bp::arg("self"),
81  "Return the inverse rotation.")
82  .def("fromRotationMatrix",&AngleAxis::template fromRotationMatrix<Matrix3>,
83  (bp::arg("self"),bp::arg("rotation matrix")),
84  "Sets *this from a 3x3 rotation matrix",
85  bp::return_self<>())
86  .def("toRotationMatrix",
87 // bp::arg("self"),
88  &AngleAxis::toRotationMatrix,
89  "Constructs and returns an equivalent 3x3 rotation matrix.")
90  .def("matrix",&AngleAxis::matrix,
91  bp::arg("self"),
92  "Returns an equivalent rotation matrix.")
93 
94  .def("isApprox",
96  isApproxAngleAxis_overload(bp::args("self","other","prec"),
97  "Returns true if *this is approximately equal to other, within the precision determined by prec."))
98 
99  /* --- Operators --- */
100  .def(bp::self * bp::other<Vector3>())
101  .def(bp::self * bp::other<Quaternion>())
102  .def(bp::self * bp::self)
103  .def("__eq__",&AngleAxisVisitor::__eq__)
104  .def("__ne__",&AngleAxisVisitor::__ne__)
105 
106  .def("__str__",&print)
107  .def("__repr__",&print)
108  ;
109  }
110 
111  private:
112 
113  static void setAxis(AngleAxis & self, const Vector3 & axis)
114  { self.axis() = axis; }
115 
116  static void setAngle( AngleAxis & self, const Scalar & angle)
117  { self.angle() = angle; }
118 
119  static bool __eq__(const AngleAxis & u, const AngleAxis & v)
120  { return u.axis() == v.axis() && v.angle() == u.angle(); }
121 
122  static bool __ne__(const AngleAxis & u, const AngleAxis & v)
123  { return !__eq__(u,v); }
124 
125  static std::string print(const AngleAxis & self)
126  {
127  std::stringstream ss;
128  ss << "angle: " << self.angle() << std::endl;
129  ss << "axis: " << self.axis().transpose() << std::endl;
130 
131  return ss.str();
132  }
133 
134  public:
135 
136  static void expose()
137  {
138  bp::class_<AngleAxis>("AngleAxis",
139  "AngleAxis representation of a rotation.\n\n",
140  bp::no_init)
142 
143  // Cast to Eigen::RotationBase
144  bp::implicitly_convertible<AngleAxis,RotationBase>();
145  }
146 
147  };
148 
149 } // namespace eigenpy
150 
151 #endif //ifndef __eigenpy_angle_axis_hpp__
Eigen::Quaternion< Scalar, 0 > Quaternion
Definition: angle-axis.hpp:46
static void setAngle(AngleAxis &self, const Scalar &angle)
Definition: angle-axis.hpp:116
AngleAxis::Scalar Scalar
Definition: angle-axis.hpp:42
boost::python::object matrix()
Definition: bnpy.cpp:20
Eigen::AngleAxis< Scalar > AngleAxis
Definition: angle-axis.hpp:23
Definition: complex.cpp:7
static std::string print(const AngleAxis &self)
Definition: angle-axis.hpp:125
Allows a template specialization.
Definition: expose.hpp:16
static bool isApprox(const AngleAxis &self, const AngleAxis &other, const Scalar &prec=Eigen::NumTraits< Scalar >::dummy_precision())
Definition: angle-axis.hpp:30
static bool __ne__(const AngleAxis &u, const AngleAxis &v)
Definition: angle-axis.hpp:122
Eigen::RotationBase< AngleAxis, 3 > RotationBase
Definition: angle-axis.hpp:47
static bool __eq__(const AngleAxis &u, const AngleAxis &v)
Definition: angle-axis.hpp:119
AngleAxis::Matrix3 Matrix3
Definition: angle-axis.hpp:44
AngleAxis::Vector3 Vector3
Definition: angle-axis.hpp:43
static void setAxis(AngleAxis &self, const Vector3 &axis)
Definition: angle-axis.hpp:113


eigenpy
Author(s): Justin Carpentier, Nicolas Mansard
autogenerated on Sat Apr 17 2021 02:37:59