/// @file /// Similarity group Sim(2) - scaling, rotation and translation in 2d. #ifndef SOPHUS_SIM2_HPP #define SOPHUS_SIM2_HPP #include "rxso2.hpp" #include "sim_details.hpp" namespace Sophus { template class Sim2; using Sim2d = Sim2; using Sim2f = Sim2; } // namespace Sophus namespace Eigen { namespace internal { template struct traits> { using Scalar = Scalar_; using TranslationType = Sophus::Vector2; using RxSO2Type = Sophus::RxSO2; }; template struct traits, Options>> : traits> { using Scalar = Scalar_; using TranslationType = Map, Options>; using RxSO2Type = Map, Options>; }; template struct traits const, Options>> : traits const> { using Scalar = Scalar_; using TranslationType = Map const, Options>; using RxSO2Type = Map const, Options>; }; } // namespace internal } // namespace Eigen namespace Sophus { /// Sim2 base type - implements Sim2 class but is storage agnostic. /// /// Sim(2) is the group of rotations and translation and scaling in 2d. It is /// the semi-direct product of R+xSO(2) and the 2d Euclidean vector space. The /// class is represented using a composition of RxSO2 for scaling plus /// rotation and a 2-vector for translation. /// /// Sim(2) is neither compact, nor a commutative group. /// /// See RxSO2 for more details of the scaling + rotation representation in 2d. /// template class Sim2Base { public: using Scalar = typename Eigen::internal::traits::Scalar; using TranslationType = typename Eigen::internal::traits::TranslationType; using RxSO2Type = typename Eigen::internal::traits::RxSO2Type; /// Degrees of freedom of manifold, number of dimensions in tangent space /// (two for translation, one for rotation and one for scaling). static int constexpr DoF = 4; /// Number of internal parameters used (2-tuple for complex number, two for /// translation). static int constexpr num_parameters = 4; /// Group transformations are 3x3 matrices. static int constexpr N = 3; using Transformation = Matrix; using Point = Vector2; using HomogeneousPoint = Vector3; using Line = ParametrizedLine2; using Tangent = Vector; using Adjoint = Matrix; /// For binary operations the return type is determined with the /// ScalarBinaryOpTraits feature of Eigen. This allows mixing concrete and Map /// types, as well as other compatible scalar types such as Ceres::Jet and /// double scalars with SIM2 operations. template using ReturnScalar = typename Eigen::ScalarBinaryOpTraits< Scalar, typename OtherDerived::Scalar>::ReturnType; template using Sim2Product = Sim2>; template using PointProduct = Vector2>; template using HomogeneousPointProduct = Vector3>; /// Adjoint transformation /// /// This function return the adjoint transformation ``Ad`` of the group /// element ``A`` such that for all ``x`` it holds that /// ``hat(Ad_A * x) = A * hat(x) A^{-1}``. See hat-operator below. /// SOPHUS_FUNC Adjoint Adj() const { Adjoint res; res.setZero(); res.template block<2, 2>(0, 0) = rxso2().matrix(); res(0, 2) = translation()[1]; res(1, 2) = -translation()[0]; res.template block<2, 1>(0, 3) = -translation(); res(2, 2) = Scalar(1); res(3, 3) = Scalar(1); return res; } /// Returns copy of instance casted to NewScalarType. /// template SOPHUS_FUNC Sim2 cast() const { return Sim2(rxso2().template cast(), translation().template cast()); } /// Returns group inverse. /// SOPHUS_FUNC Sim2 inverse() const { RxSO2 invR = rxso2().inverse(); return Sim2(invR, invR * (translation() * Scalar(-1))); } /// Logarithmic map /// /// Computes the logarithm, the inverse of the group exponential which maps /// element of the group (rigid body transformations) to elements of the /// tangent space (twist). /// /// To be specific, this function computes ``vee(logmat(.))`` with /// ``logmat(.)`` being the matrix logarithm and ``vee(.)`` the vee-operator /// of Sim(2). /// SOPHUS_FUNC Tangent log() const { /// The derivation of the closed-form Sim(2) logarithm for is done /// analogously to the closed-form solution of the SE(2) logarithm, see /// J. Gallier, D. Xu, "Computing exponentials of skew symmetric matrices /// and logarithms of orthogonal matrices", IJRA 2002. /// https:///pdfs.semanticscholar.org/cfe3/e4b39de63c8cabd89bf3feff7f5449fc981d.pdf /// (Sec. 6., pp. 8) Tangent res; Vector2 const theta_sigma = rxso2().log(); Scalar const theta = theta_sigma[0]; Scalar const sigma = theta_sigma[1]; Matrix2 const Omega = SO2::hat(theta); Matrix2 const W_inv = details::calcWInv(Omega, theta, sigma, scale()); res.segment(0, 2) = W_inv * translation(); res[2] = theta; res[3] = sigma; return res; } /// Returns 3x3 matrix representation of the instance. /// /// It has the following form: /// /// | s*R t | /// | o 1 | /// /// where ``R`` is a 2x2 rotation matrix, ``s`` a scale factor, ``t`` a /// translation 2-vector and ``o`` a 2-column vector of zeros. /// SOPHUS_FUNC Transformation matrix() const { Transformation homogenious_matrix; homogenious_matrix.template topLeftCorner<2, 3>() = matrix2x3(); homogenious_matrix.row(2) = Matrix(Scalar(0), Scalar(0), Scalar(1)); return homogenious_matrix; } /// Returns the significant first two rows of the matrix above. /// SOPHUS_FUNC Matrix matrix2x3() const { Matrix matrix; matrix.template topLeftCorner<2, 2>() = rxso2().matrix(); matrix.col(2) = translation(); return matrix; } /// Assignment-like operator from OtherDerived. /// template SOPHUS_FUNC Sim2Base& operator=( Sim2Base const& other) { rxso2() = other.rxso2(); translation() = other.translation(); return *this; } /// Group multiplication, which is rotation plus scaling concatenation. /// /// Note: That scaling is calculated with saturation. See RxSO2 for /// details. /// template SOPHUS_FUNC Sim2Product operator*( Sim2Base const& other) const { return Sim2Product( rxso2() * other.rxso2(), translation() + rxso2() * other.translation()); } /// Group action on 2-points. /// /// This function rotates, scales and translates a two dimensional point /// ``p`` by the Sim(2) element ``(bar_sR_foo, t_bar)`` (= similarity /// transformation): /// /// ``p_bar = bar_sR_foo * p_foo + t_bar``. /// template ::value>::type> SOPHUS_FUNC PointProduct operator*( Eigen::MatrixBase const& p) const { return rxso2() * p + translation(); } /// Group action on homogeneous 2-points. See above for more details. /// template ::value>::type> SOPHUS_FUNC HomogeneousPointProduct operator*( Eigen::MatrixBase const& p) const { const PointProduct tp = rxso2() * p.template head<2>() + p(2) * translation(); return HomogeneousPointProduct(tp(0), tp(1), p(2)); } /// Group action on lines. /// /// This function rotates, scales and translates a parametrized line /// ``l(t) = o + t * d`` by the Sim(2) element: /// /// Origin ``o`` is rotated, scaled and translated /// Direction ``d`` is rotated /// SOPHUS_FUNC Line operator*(Line const& l) const { Line rotatedLine = rxso2() * l; return Line(rotatedLine.origin() + translation(), rotatedLine.direction()); } /// Returns internal parameters of Sim(2). /// /// It returns (c[0], c[1], t[0], t[1]), /// with c being the complex number, t the translation 3-vector. /// SOPHUS_FUNC Sophus::Vector params() const { Sophus::Vector p; p << rxso2().params(), translation(); return p; } /// In-place group multiplication. This method is only valid if the return /// type of the multiplication is compatible with this SO2's Scalar type. /// template >::value>::type> SOPHUS_FUNC Sim2Base& operator*=( Sim2Base const& other) { *static_cast(this) = *this * other; return *this; } /// Setter of non-zero complex number. /// /// Precondition: ``z`` must not be close to zero. /// SOPHUS_FUNC void setComplex(Vector2 const& z) { rxso2().setComplex(z); } /// Accessor of complex number. /// SOPHUS_FUNC typename Eigen::internal::traits::RxSO2Type::ComplexType const& complex() const { return rxso2().complex(); } /// Returns Rotation matrix /// SOPHUS_FUNC Matrix2 rotationMatrix() const { return rxso2().rotationMatrix(); } /// Mutator of SO2 group. /// SOPHUS_FUNC RxSO2Type& rxso2() { return static_cast(this)->rxso2(); } /// Accessor of SO2 group. /// SOPHUS_FUNC RxSO2Type const& rxso2() const { return static_cast(this)->rxso2(); } /// Returns scale. /// SOPHUS_FUNC Scalar scale() const { return rxso2().scale(); } /// Setter of complex number using rotation matrix ``R``, leaves scale as is. /// SOPHUS_FUNC void setRotationMatrix(Matrix2& R) { rxso2().setRotationMatrix(R); } /// Sets scale and leaves rotation as is. /// /// Note: This function as a significant computational cost, since it has to /// call the square root twice. /// SOPHUS_FUNC void setScale(Scalar const& scale) { rxso2().setScale(scale); } /// Setter of complexnumber using scaled rotation matrix ``sR``. /// /// Precondition: The 2x2 matrix must be "scaled orthogonal" /// and have a positive determinant. /// SOPHUS_FUNC void setScaledRotationMatrix(Matrix2 const& sR) { rxso2().setScaledRotationMatrix(sR); } /// Mutator of translation vector /// SOPHUS_FUNC TranslationType& translation() { return static_cast(this)->translation(); } /// Accessor of translation vector /// SOPHUS_FUNC TranslationType const& translation() const { return static_cast(this)->translation(); } }; /// Sim2 using default storage; derived from Sim2Base. template class Sim2 : public Sim2Base> { public: using Base = Sim2Base>; using Scalar = Scalar_; using Transformation = typename Base::Transformation; using Point = typename Base::Point; using HomogeneousPoint = typename Base::HomogeneousPoint; using Tangent = typename Base::Tangent; using Adjoint = typename Base::Adjoint; using RxSo2Member = RxSO2; using TranslationMember = Vector2; using Base::operator=; EIGEN_MAKE_ALIGNED_OPERATOR_NEW /// Default constructor initializes similarity transform to the identity. /// SOPHUS_FUNC Sim2(); /// Copy constructor /// SOPHUS_FUNC Sim2(Sim2 const& other) = default; /// Copy-like constructor from OtherDerived. /// template SOPHUS_FUNC Sim2(Sim2Base const& other) : rxso2_(other.rxso2()), translation_(other.translation()) { static_assert(std::is_same::value, "must be same Scalar type"); } /// Constructor from RxSO2 and translation vector /// template SOPHUS_FUNC Sim2(RxSO2Base const& rxso2, Eigen::MatrixBase const& translation) : rxso2_(rxso2), translation_(translation) { static_assert(std::is_same::value, "must be same Scalar type"); static_assert(std::is_same::value, "must be same Scalar type"); } /// Constructor from complex number and translation vector. /// /// Precondition: complex number must not be close to zero. /// template SOPHUS_FUNC Sim2(Vector2 const& complex_number, Eigen::MatrixBase const& translation) : rxso2_(complex_number), translation_(translation) { static_assert(std::is_same::value, "must be same Scalar type"); } /// Constructor from 3x3 matrix /// /// Precondition: Top-left 2x2 matrix needs to be "scaled-orthogonal" with /// positive determinant. The last row must be ``(0, 0, 1)``. /// SOPHUS_FUNC explicit Sim2(Matrix const& T) : rxso2_((T.template topLeftCorner<2, 2>()).eval()), translation_(T.template block<2, 1>(0, 2)) {} /// This provides unsafe read/write access to internal data. Sim(2) is /// represented by a complex number (two parameters) and a 2-vector. When /// using direct write access, the user needs to take care of that the /// complex number is not set close to zero. /// SOPHUS_FUNC Scalar* data() { // rxso2_ and translation_ are laid out sequentially with no padding return rxso2_.data(); } /// Const version of data() above. /// SOPHUS_FUNC Scalar const* data() const { // rxso2_ and translation_ are laid out sequentially with no padding return rxso2_.data(); } /// Accessor of RxSO2 /// SOPHUS_FUNC RxSo2Member& rxso2() { return rxso2_; } /// Mutator of RxSO2 /// SOPHUS_FUNC RxSo2Member const& rxso2() const { return rxso2_; } /// Mutator of translation vector /// SOPHUS_FUNC TranslationMember& translation() { return translation_; } /// Accessor of translation vector /// SOPHUS_FUNC TranslationMember const& translation() const { return translation_; } /// Returns derivative of exp(x).matrix() wrt. ``x_i at x=0``. /// SOPHUS_FUNC static Transformation Dxi_exp_x_matrix_at_0(int i) { return generator(i); } /// Derivative of Lie bracket with respect to first element. /// /// This function returns ``D_a [a, b]`` with ``D_a`` being the /// differential operator with respect to ``a``, ``[a, b]`` being the lie /// bracket of the Lie algebra sim(2). /// See ``lieBracket()`` below. /// /// Group exponential /// /// This functions takes in an element of tangent space and returns the /// corresponding element of the group Sim(2). /// /// The first two components of ``a`` represent the translational part /// ``upsilon`` in the tangent space of Sim(2), the following two components /// of ``a`` represents the rotation ``theta`` and the final component /// represents the logarithm of the scaling factor ``sigma``. /// To be more specific, this function computes ``expmat(hat(a))`` with /// ``expmat(.)`` being the matrix exponential and ``hat(.)`` the hat-operator /// of Sim(2), see below. /// SOPHUS_FUNC static Sim2 exp(Tangent const& a) { // For the derivation of the exponential map of Sim(N) see // H. Strasdat, "Local Accuracy and Global Consistency for Efficient Visual // SLAM", PhD thesis, 2012. // http:///hauke.strasdat.net/files/strasdat_thesis_2012.pdf (A.5, pp. 186) Vector2 const upsilon = a.segment(0, 2); Scalar const theta = a[2]; Scalar const sigma = a[3]; RxSO2 rxso2 = RxSO2::exp(a.template tail<2>()); Matrix2 const Omega = SO2::hat(theta); Matrix2 const W = details::calcW(Omega, theta, sigma); return Sim2(rxso2, W * upsilon); } /// Returns the ith infinitesimal generators of Sim(2). /// /// The infinitesimal generators of Sim(2) are: /// /// ``` /// | 0 0 1 | /// G_0 = | 0 0 0 | /// | 0 0 0 | /// /// | 0 0 0 | /// G_1 = | 0 0 1 | /// | 0 0 0 | /// /// | 0 -1 0 | /// G_2 = | 1 0 0 | /// | 0 0 0 | /// /// | 1 0 0 | /// G_3 = | 0 1 0 | /// | 0 0 0 | /// ``` /// /// Precondition: ``i`` must be in [0, 3]. /// SOPHUS_FUNC static Transformation generator(int i) { SOPHUS_ENSURE(i >= 0 || i <= 3, "i should be in range [0,3]."); Tangent e; e.setZero(); e[i] = Scalar(1); return hat(e); } /// hat-operator /// /// It takes in the 4-vector representation and returns the corresponding /// matrix representation of Lie algebra element. /// /// Formally, the hat()-operator of Sim(2) is defined as /// /// ``hat(.): R^4 -> R^{3x3}, hat(a) = sum_i a_i * G_i`` (for i=0,...,6) /// /// with ``G_i`` being the ith infinitesimal generator of Sim(2). /// /// The corresponding inverse is the vee()-operator, see below. /// SOPHUS_FUNC static Transformation hat(Tangent const& a) { Transformation Omega; Omega.template topLeftCorner<2, 2>() = RxSO2::hat(a.template tail<2>()); Omega.col(2).template head<2>() = a.template head<2>(); Omega.row(2).setZero(); return Omega; } /// Lie bracket /// /// It computes the Lie bracket of Sim(2). To be more specific, it computes /// /// ``[omega_1, omega_2]_sim2 := vee([hat(omega_1), hat(omega_2)])`` /// /// with ``[A,B] := AB-BA`` being the matrix commutator, ``hat(.)`` the /// hat()-operator and ``vee(.)`` the vee()-operator of Sim(2). /// SOPHUS_FUNC static Tangent lieBracket(Tangent const& a, Tangent const& b) { Vector2 const upsilon1 = a.template head<2>(); Vector2 const upsilon2 = b.template head<2>(); Scalar const theta1 = a[2]; Scalar const theta2 = b[2]; Scalar const sigma1 = a[3]; Scalar const sigma2 = b[3]; Tangent res; res[0] = -theta1 * upsilon2[1] + theta2 * upsilon1[1] + sigma1 * upsilon2[0] - sigma2 * upsilon1[0]; res[1] = theta1 * upsilon2[0] - theta2 * upsilon1[0] + sigma1 * upsilon2[1] - sigma2 * upsilon1[1]; res[2] = Scalar(0); res[3] = Scalar(0); return res; } /// Draw uniform sample from Sim(2) manifold. /// /// Translations are drawn component-wise from the range [-1, 1]. /// The scale factor is drawn uniformly in log2-space from [-1, 1], /// hence the scale is in [0.5, 2]. /// template static Sim2 sampleUniform(UniformRandomBitGenerator& generator) { std::uniform_real_distribution uniform(Scalar(-1), Scalar(1)); return Sim2(RxSO2::sampleUniform(generator), Vector2(uniform(generator), uniform(generator))); } /// vee-operator /// /// It takes the 3x3-matrix representation ``Omega`` and maps it to the /// corresponding 4-vector representation of Lie algebra. /// /// This is the inverse of the hat()-operator, see above. /// /// Precondition: ``Omega`` must have the following structure: /// /// | d -c a | /// | c d b | /// | 0 0 0 | /// SOPHUS_FUNC static Tangent vee(Transformation const& Omega) { Tangent upsilon_omega_sigma; upsilon_omega_sigma.template head<2>() = Omega.col(2).template head<2>(); upsilon_omega_sigma.template tail<2>() = RxSO2::vee(Omega.template topLeftCorner<2, 2>()); return upsilon_omega_sigma; } protected: RxSo2Member rxso2_; TranslationMember translation_; }; template Sim2::Sim2() : translation_(TranslationMember::Zero()) { static_assert(std::is_standard_layout::value, "Assume standard layout for the use of offsetof check below."); static_assert( offsetof(Sim2, rxso2_) + sizeof(Scalar) * RxSO2::num_parameters == offsetof(Sim2, translation_), "This class assumes packed storage and hence will only work " "correctly depending on the compiler (options) - in " "particular when using [this->data(), this-data() + " "num_parameters] to access the raw data in a contiguous fashion."); } } // namespace Sophus namespace Eigen { /// Specialization of Eigen::Map for ``Sim2``; derived from Sim2Base. /// /// Allows us to wrap Sim2 objects around POD array. template class Map, Options> : public Sophus::Sim2Base, Options>> { public: using Base = Sophus::Sim2Base, Options>>; using Scalar = Scalar_; using Transformation = typename Base::Transformation; using Point = typename Base::Point; using HomogeneousPoint = typename Base::HomogeneousPoint; using Tangent = typename Base::Tangent; using Adjoint = typename Base::Adjoint; using Base::operator=; using Base::operator*=; using Base::operator*; SOPHUS_FUNC Map(Scalar* coeffs) : rxso2_(coeffs), translation_(coeffs + Sophus::RxSO2::num_parameters) {} /// Mutator of RxSO2 /// SOPHUS_FUNC Map, Options>& rxso2() { return rxso2_; } /// Accessor of RxSO2 /// SOPHUS_FUNC Map, Options> const& rxso2() const { return rxso2_; } /// Mutator of translation vector /// SOPHUS_FUNC Map, Options>& translation() { return translation_; } /// Accessor of translation vector SOPHUS_FUNC Map, Options> const& translation() const { return translation_; } protected: Map, Options> rxso2_; Map, Options> translation_; }; /// Specialization of Eigen::Map for ``Sim2 const``; derived from Sim2Base. /// /// Allows us to wrap RxSO2 objects around POD array. template class Map const, Options> : public Sophus::Sim2Base const, Options>> { public: using Base = Sophus::Sim2Base const, Options>>; using Scalar = Scalar_; using Transformation = typename Base::Transformation; using Point = typename Base::Point; using HomogeneousPoint = typename Base::HomogeneousPoint; using Tangent = typename Base::Tangent; using Adjoint = typename Base::Adjoint; using Base::operator*=; using Base::operator*; SOPHUS_FUNC Map(Scalar const* coeffs) : rxso2_(coeffs), translation_(coeffs + Sophus::RxSO2::num_parameters) {} /// Accessor of RxSO2 /// SOPHUS_FUNC Map const, Options> const& rxso2() const { return rxso2_; } /// Accessor of translation vector /// SOPHUS_FUNC Map const, Options> const& translation() const { return translation_; } protected: Map const, Options> const rxso2_; Map const, Options> const translation_; }; } // namespace Eigen #endif