Class sym::Rot3#

template<typename ScalarType>
class Rot3#

Autogenerated C++ implementation of symforce.geo.rot3.Rot3.

Group of three-dimensional orthogonal matrices with determinant +1, representing rotations in 3D space. Backed by a quaternion with (x, y, z, w) storage.

Public Types

using Scalar = ScalarType#
using Self = Rot3<Scalar>#
using DataVec = Eigen::Matrix<Scalar, 4, 1>#
using TangentVec = Eigen::Matrix<Scalar, 3, 1>#
using SelfJacobian = Eigen::Matrix<Scalar, 3, 3>#
using Vector3 = Eigen::Matrix<Scalar, 3, 1>#

Public Functions

inline explicit Rot3(const DataVec &data, const bool normalize = true)#

Construct from data vec

Parameters:

normalize – Project to the manifold on construction. This ensures numerical stability as this constructor is called after each codegen operation. Constructing from a normalized vector may be faster, e.g. with FromStorage.

inline Rot3()#
inline const DataVec &Data() const#
inline Eigen::Quaternion<Scalar> Quaternion() const#
inline explicit Rot3(const Eigen::Quaternion<Scalar> &q)#
inline Eigen::AngleAxis<Scalar> AngleAxis() const#
inline explicit Rot3(const Eigen::AngleAxis<Scalar> &aa)#
inline Rot3 ToPositiveReal() const#
const Vector3 ComposeWithPoint(const Vector3 &right) const#

Left-multiplication. Either rotation concatenation or point transform.

const Scalar ToTangentNorm(const Scalar epsilon) const#

Returns the norm of the tangent vector corresponding to this rotation

This is equal to the angle that should be rotated through to get this Rot3, in radians. Using this function directly is usually more efficient than computing the norm of the tangent vector, both in symbolic and generated code; by default, symbolic APIs will not automatically simplify to this

const Eigen::Matrix<Scalar, 3, 3> ToRotationMatrix() const#

Converts to a rotation matrix

const Vector3 ToYawPitchRoll() const#

Compute the yaw, pitch, and roll Euler angles in radians of this rotation

Euler angles are subject to gimbal lock: https://en.wikipedia.org/wiki/Gimbal_lock

This means that when the pitch is close to +/- pi/2, the yaw and roll angles are not uniquely defined, so the returned values are not unique in this case.

Returns: Scalar: Yaw angle [radians] Scalar: Pitch angle [radians] Scalar: Roll angle [radians]

inline void ToStorage(Scalar *const vec) const#
inline Self Inverse() const#
inline Self Compose(const Self &b) const#
inline Vector3 Compose(const Vector3 &point) const#
inline Self Between(const Self &b) const#
inline Self InverseWithJacobian(SelfJacobian *const res_D_a = nullptr) const#
inline Self ComposeWithJacobians(const Self &b, SelfJacobian *const res_D_a = nullptr, SelfJacobian *const res_D_b = nullptr) const#
inline Self BetweenWithJacobians(const Self &b, SelfJacobian *const res_D_a = nullptr, SelfJacobian *const res_D_b = nullptr) const#
template<typename Other>
inline auto operator*(const Other &b) const -> decltype(Compose(b))#
inline TangentVec ToTangent(const Scalar epsilon = kDefaultEpsilon<Scalar>) const#
inline Self Retract(const TangentVec &vec, const Scalar epsilon = kDefaultEpsilon<Scalar>) const#
inline TangentVec LocalCoordinates(const Self &b, const Scalar epsilon = kDefaultEpsilon<Scalar>) const#
inline Self Interpolate(const Self b, const Scalar alpha, const Scalar epsilon = kDefaultEpsilon<Scalar>) const#
inline bool IsApprox(const Self &b, const Scalar tol) const#
template<typename ToScalar>
inline Rot3<ToScalar> Cast() const#
inline bool operator==(const Rot3 &rhs) const#

Public Static Functions

static inline Rot3 FromQuaternion(const Eigen::Quaternion<Scalar> &q)#
static inline Rot3 FromAngleAxis(const Eigen::AngleAxis<Scalar> &aa)#
static inline Rot3 FromAngleAxis(const Scalar angle, const Vector3 &axis)#
static inline Rot3 FromRotationMatrix(const Eigen::Matrix<Scalar, 3, 3> &mat)#
template<typename Generator>
static inline Rot3 Random(Generator &gen)#
static const sym::Rot3<Scalar> RandomFromUniformSamples(const Scalar u1, const Scalar u2, const Scalar u3)#

Generate a random element of SO3 from three variables uniformly sampled in [0, 1].

static const sym::Rot3<Scalar> FromYawPitchRoll(const Scalar yaw, const Scalar pitch, const Scalar roll)#

Construct from yaw, pitch, and roll Euler angles in radians

static const sym::Rot3<Scalar> FromYawPitchRoll(const Vector3 &ypr)#

This function was autogenerated from a symbolic function. Do not modify by hand.

Symbolic function: <lambda>

Args: ypr: Matrix31

Outputs: res: Rot3

static const sym::Rot3<Scalar> FromTwoUnitVectors(const Vector3 &a, const Vector3 &b, const Scalar epsilon)#

Return a rotation that transforms a to b. Both inputs are three-vectors that are expected to be normalized.

Reference: http://lolengine.net/blog/2013/09/18/beautiful-maths-quaternion-from-vectors

static inline constexpr int32_t StorageDim()#
static inline Rot3 FromStorage(const Scalar *const vec)#
static inline Self Identity()#
static inline constexpr int32_t TangentDim()#
static inline Self FromTangent(const TangentVec &vec, const Scalar epsilon = kDefaultEpsilon<Scalar>)#