isce3 0.25.0
Loading...
Searching...
No Matches
EulerAngles.h
1#pragma once
2#define EIGEN_MPL2_ONLY
3
4#include <Eigen/Geometry>
5
6#include "DenseMatrix.h"
7#include "Quaternion.h"
8#include "Vector.h"
9
10namespace isce3 { namespace core {
11
14public:
21 EulerAngles(double y, double p, double r) : _yaw(y), _pitch(p), _roll(r) {}
22
30 explicit EulerAngles(const Mat3& rotmat);
31
36 explicit EulerAngles(const Quaternion& quat);
37
42 Mat3 toRotationMatrix() const
43 {
44 return (Eigen::AngleAxisd(yaw(), Eigen::Vector3d::UnitZ()) *
45 Eigen::AngleAxisd(pitch(), Eigen::Vector3d::UnitY()) *
46 Eigen::AngleAxisd(roll(), Eigen::Vector3d::UnitX()))
47 .toRotationMatrix();
48 }
49
51 double yaw() const { return _yaw; }
52
54 double pitch() const { return _pitch; }
55
57 double roll() const { return _roll; }
58
64
74 bool isApprox(const EulerAngles& other, double prec = 1e-7) const;
75
81 Vec3 rotate(const Eigen::Vector3d& vec) const;
82
89 {
90 _yaw += rhs._yaw;
91 _pitch += rhs._pitch;
92 _roll += rhs._roll;
93 return *this;
94 }
95
102 {
103 _yaw -= rhs._yaw;
104 _pitch -= rhs._pitch;
105 _roll -= rhs._roll;
106 return *this;
107 }
108
116 {
117 *this = EulerAngles(this->toRotationMatrix() * rhs.toRotationMatrix());
118 return *this;
119 }
120
121private:
122 double _yaw, _pitch, _roll;
123};
124
125// Regular/non-member functions
126
133inline EulerAngles operator+(const EulerAngles& lhs, const EulerAngles& rhs)
134{
135 auto sum {lhs};
136 sum += rhs;
137 return sum;
138}
139
146inline EulerAngles operator-(const EulerAngles& lhs, const EulerAngles& rhs)
147{
148 auto sub {lhs};
149 sub -= rhs;
150 return sub;
151}
152
159inline EulerAngles operator*(const EulerAngles& lhs, const EulerAngles& rhs)
160{
161 auto mul {lhs};
162 mul *= rhs;
163 return mul;
164}
165
166}} // namespace isce3::core
Representation of 3-2-1 Euler angle sequence of rotations.
Definition EulerAngles.h:13
double yaw() const
Get yaw in radians.
Definition EulerAngles.h:51
Mat3 toRotationMatrix() const
Convert to rotation matrix.
Definition EulerAngles.h:42
EulerAngles & operator-=(const EulerAngles &rhs)
Overloaded in-place operator -= on this.
Definition EulerAngles.h:101
double roll() const
Get roll in radians.
Definition EulerAngles.h:57
Vec3 rotate(const Eigen::Vector3d &vec) const
Rotate a 3-D vector by Euler object in YPR order.
Definition EulerAngles.cpp:48
EulerAngles & operator+=(const EulerAngles &rhs)
Overloaded in-place operator += on this.
Definition EulerAngles.h:88
Quaternion toQuaternion() const
Convert to isce3 Quaternion object.
Definition EulerAngles.cpp:46
EulerAngles(double y, double p, double r)
Construct from yaw, pitch, and roll angles.
Definition EulerAngles.h:21
EulerAngles & operator*=(const EulerAngles &rhs)
Overloaded in-place operator *= on this.
Definition EulerAngles.h:115
bool isApprox(const EulerAngles &other, double prec=1e-7) const
Check if *this is approximately equals to other within a desired precision.
Definition EulerAngles.cpp:33
double pitch() const
Get pitch in radians.
Definition EulerAngles.h:54
Quaternion representation of rotations, based on double precision Eigen::Quaterniond.
Definition Quaternion.h:25
base interpolator is an abstract base class
Definition BinarySearchFunc.cpp:5

Generated for ISCE3.0 by doxygen 1.13.2.