Eigne Matrix Class Library
Dependents: MPC_current_control HydraulicControlBoard_SW AHRS Test_ekf ... more
src/Geometry/OrthoMethods.h@1:3b8049da21b8, 2019-09-24 (annotated)
- Committer:
- jsoh91
- Date:
- Tue Sep 24 00:18:23 2019 +0000
- Revision:
- 1:3b8049da21b8
- Parent:
- 0:13a5d365ba16
ignore and revise some of error parts
Who changed what in which revision?
User | Revision | Line number | New contents of line |
---|---|---|---|
ykuroda | 0:13a5d365ba16 | 1 | // This file is part of Eigen, a lightweight C++ template library |
ykuroda | 0:13a5d365ba16 | 2 | // for linear algebra. |
ykuroda | 0:13a5d365ba16 | 3 | // |
ykuroda | 0:13a5d365ba16 | 4 | // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr> |
ykuroda | 0:13a5d365ba16 | 5 | // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com> |
ykuroda | 0:13a5d365ba16 | 6 | // |
ykuroda | 0:13a5d365ba16 | 7 | // This Source Code Form is subject to the terms of the Mozilla |
ykuroda | 0:13a5d365ba16 | 8 | // Public License v. 2.0. If a copy of the MPL was not distributed |
ykuroda | 0:13a5d365ba16 | 9 | // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. |
ykuroda | 0:13a5d365ba16 | 10 | |
ykuroda | 0:13a5d365ba16 | 11 | #ifndef EIGEN_ORTHOMETHODS_H |
ykuroda | 0:13a5d365ba16 | 12 | #define EIGEN_ORTHOMETHODS_H |
ykuroda | 0:13a5d365ba16 | 13 | |
ykuroda | 0:13a5d365ba16 | 14 | namespace Eigen { |
ykuroda | 0:13a5d365ba16 | 15 | |
ykuroda | 0:13a5d365ba16 | 16 | /** \geometry_module |
ykuroda | 0:13a5d365ba16 | 17 | * |
ykuroda | 0:13a5d365ba16 | 18 | * \returns the cross product of \c *this and \a other |
ykuroda | 0:13a5d365ba16 | 19 | * |
ykuroda | 0:13a5d365ba16 | 20 | * Here is a very good explanation of cross-product: http://xkcd.com/199/ |
ykuroda | 0:13a5d365ba16 | 21 | * \sa MatrixBase::cross3() |
ykuroda | 0:13a5d365ba16 | 22 | */ |
ykuroda | 0:13a5d365ba16 | 23 | template<typename Derived> |
ykuroda | 0:13a5d365ba16 | 24 | template<typename OtherDerived> |
ykuroda | 0:13a5d365ba16 | 25 | inline typename MatrixBase<Derived>::template cross_product_return_type<OtherDerived>::type |
ykuroda | 0:13a5d365ba16 | 26 | MatrixBase<Derived>::cross(const MatrixBase<OtherDerived>& other) const |
ykuroda | 0:13a5d365ba16 | 27 | { |
ykuroda | 0:13a5d365ba16 | 28 | EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,3) |
ykuroda | 0:13a5d365ba16 | 29 | EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3) |
ykuroda | 0:13a5d365ba16 | 30 | |
ykuroda | 0:13a5d365ba16 | 31 | // Note that there is no need for an expression here since the compiler |
ykuroda | 0:13a5d365ba16 | 32 | // optimize such a small temporary very well (even within a complex expression) |
ykuroda | 0:13a5d365ba16 | 33 | typename internal::nested<Derived,2>::type lhs(derived()); |
ykuroda | 0:13a5d365ba16 | 34 | typename internal::nested<OtherDerived,2>::type rhs(other.derived()); |
ykuroda | 0:13a5d365ba16 | 35 | return typename cross_product_return_type<OtherDerived>::type( |
ykuroda | 0:13a5d365ba16 | 36 | numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)), |
ykuroda | 0:13a5d365ba16 | 37 | numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)), |
ykuroda | 0:13a5d365ba16 | 38 | numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0)) |
ykuroda | 0:13a5d365ba16 | 39 | ); |
ykuroda | 0:13a5d365ba16 | 40 | } |
ykuroda | 0:13a5d365ba16 | 41 | |
ykuroda | 0:13a5d365ba16 | 42 | namespace internal { |
ykuroda | 0:13a5d365ba16 | 43 | |
ykuroda | 0:13a5d365ba16 | 44 | template< int Arch,typename VectorLhs,typename VectorRhs, |
ykuroda | 0:13a5d365ba16 | 45 | typename Scalar = typename VectorLhs::Scalar, |
ykuroda | 0:13a5d365ba16 | 46 | bool Vectorizable = bool((VectorLhs::Flags&VectorRhs::Flags)&PacketAccessBit)> |
ykuroda | 0:13a5d365ba16 | 47 | struct cross3_impl { |
ykuroda | 0:13a5d365ba16 | 48 | static inline typename internal::plain_matrix_type<VectorLhs>::type |
ykuroda | 0:13a5d365ba16 | 49 | run(const VectorLhs& lhs, const VectorRhs& rhs) |
ykuroda | 0:13a5d365ba16 | 50 | { |
ykuroda | 0:13a5d365ba16 | 51 | return typename internal::plain_matrix_type<VectorLhs>::type( |
ykuroda | 0:13a5d365ba16 | 52 | numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)), |
ykuroda | 0:13a5d365ba16 | 53 | numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)), |
ykuroda | 0:13a5d365ba16 | 54 | numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0)), |
ykuroda | 0:13a5d365ba16 | 55 | 0 |
ykuroda | 0:13a5d365ba16 | 56 | ); |
ykuroda | 0:13a5d365ba16 | 57 | } |
ykuroda | 0:13a5d365ba16 | 58 | }; |
ykuroda | 0:13a5d365ba16 | 59 | |
ykuroda | 0:13a5d365ba16 | 60 | } |
ykuroda | 0:13a5d365ba16 | 61 | |
ykuroda | 0:13a5d365ba16 | 62 | /** \geometry_module |
ykuroda | 0:13a5d365ba16 | 63 | * |
ykuroda | 0:13a5d365ba16 | 64 | * \returns the cross product of \c *this and \a other using only the x, y, and z coefficients |
ykuroda | 0:13a5d365ba16 | 65 | * |
ykuroda | 0:13a5d365ba16 | 66 | * The size of \c *this and \a other must be four. This function is especially useful |
ykuroda | 0:13a5d365ba16 | 67 | * when using 4D vectors instead of 3D ones to get advantage of SSE/AltiVec vectorization. |
ykuroda | 0:13a5d365ba16 | 68 | * |
ykuroda | 0:13a5d365ba16 | 69 | * \sa MatrixBase::cross() |
ykuroda | 0:13a5d365ba16 | 70 | */ |
ykuroda | 0:13a5d365ba16 | 71 | template<typename Derived> |
ykuroda | 0:13a5d365ba16 | 72 | template<typename OtherDerived> |
ykuroda | 0:13a5d365ba16 | 73 | inline typename MatrixBase<Derived>::PlainObject |
ykuroda | 0:13a5d365ba16 | 74 | MatrixBase<Derived>::cross3(const MatrixBase<OtherDerived>& other) const |
ykuroda | 0:13a5d365ba16 | 75 | { |
ykuroda | 0:13a5d365ba16 | 76 | EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,4) |
ykuroda | 0:13a5d365ba16 | 77 | EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,4) |
ykuroda | 0:13a5d365ba16 | 78 | |
ykuroda | 0:13a5d365ba16 | 79 | typedef typename internal::nested<Derived,2>::type DerivedNested; |
ykuroda | 0:13a5d365ba16 | 80 | typedef typename internal::nested<OtherDerived,2>::type OtherDerivedNested; |
ykuroda | 0:13a5d365ba16 | 81 | DerivedNested lhs(derived()); |
ykuroda | 0:13a5d365ba16 | 82 | OtherDerivedNested rhs(other.derived()); |
ykuroda | 0:13a5d365ba16 | 83 | |
ykuroda | 0:13a5d365ba16 | 84 | return internal::cross3_impl<Architecture::Target, |
ykuroda | 0:13a5d365ba16 | 85 | typename internal::remove_all<DerivedNested>::type, |
ykuroda | 0:13a5d365ba16 | 86 | typename internal::remove_all<OtherDerivedNested>::type>::run(lhs,rhs); |
ykuroda | 0:13a5d365ba16 | 87 | } |
ykuroda | 0:13a5d365ba16 | 88 | |
ykuroda | 0:13a5d365ba16 | 89 | /** \returns a matrix expression of the cross product of each column or row |
ykuroda | 0:13a5d365ba16 | 90 | * of the referenced expression with the \a other vector. |
ykuroda | 0:13a5d365ba16 | 91 | * |
ykuroda | 0:13a5d365ba16 | 92 | * The referenced matrix must have one dimension equal to 3. |
ykuroda | 0:13a5d365ba16 | 93 | * The result matrix has the same dimensions than the referenced one. |
ykuroda | 0:13a5d365ba16 | 94 | * |
ykuroda | 0:13a5d365ba16 | 95 | * \geometry_module |
ykuroda | 0:13a5d365ba16 | 96 | * |
ykuroda | 0:13a5d365ba16 | 97 | * \sa MatrixBase::cross() */ |
ykuroda | 0:13a5d365ba16 | 98 | template<typename ExpressionType, int Direction> |
ykuroda | 0:13a5d365ba16 | 99 | template<typename OtherDerived> |
ykuroda | 0:13a5d365ba16 | 100 | const typename VectorwiseOp<ExpressionType,Direction>::CrossReturnType |
ykuroda | 0:13a5d365ba16 | 101 | VectorwiseOp<ExpressionType,Direction>::cross(const MatrixBase<OtherDerived>& other) const |
ykuroda | 0:13a5d365ba16 | 102 | { |
ykuroda | 0:13a5d365ba16 | 103 | EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3) |
ykuroda | 0:13a5d365ba16 | 104 | EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename OtherDerived::Scalar>::value), |
ykuroda | 0:13a5d365ba16 | 105 | YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY) |
ykuroda | 0:13a5d365ba16 | 106 | |
ykuroda | 0:13a5d365ba16 | 107 | CrossReturnType res(_expression().rows(),_expression().cols()); |
ykuroda | 0:13a5d365ba16 | 108 | if(Direction==Vertical) |
ykuroda | 0:13a5d365ba16 | 109 | { |
ykuroda | 0:13a5d365ba16 | 110 | eigen_assert(CrossReturnType::RowsAtCompileTime==3 && "the matrix must have exactly 3 rows"); |
ykuroda | 0:13a5d365ba16 | 111 | res.row(0) = (_expression().row(1) * other.coeff(2) - _expression().row(2) * other.coeff(1)).conjugate(); |
ykuroda | 0:13a5d365ba16 | 112 | res.row(1) = (_expression().row(2) * other.coeff(0) - _expression().row(0) * other.coeff(2)).conjugate(); |
ykuroda | 0:13a5d365ba16 | 113 | res.row(2) = (_expression().row(0) * other.coeff(1) - _expression().row(1) * other.coeff(0)).conjugate(); |
ykuroda | 0:13a5d365ba16 | 114 | } |
ykuroda | 0:13a5d365ba16 | 115 | else |
ykuroda | 0:13a5d365ba16 | 116 | { |
ykuroda | 0:13a5d365ba16 | 117 | eigen_assert(CrossReturnType::ColsAtCompileTime==3 && "the matrix must have exactly 3 columns"); |
ykuroda | 0:13a5d365ba16 | 118 | res.col(0) = (_expression().col(1) * other.coeff(2) - _expression().col(2) * other.coeff(1)).conjugate(); |
ykuroda | 0:13a5d365ba16 | 119 | res.col(1) = (_expression().col(2) * other.coeff(0) - _expression().col(0) * other.coeff(2)).conjugate(); |
ykuroda | 0:13a5d365ba16 | 120 | res.col(2) = (_expression().col(0) * other.coeff(1) - _expression().col(1) * other.coeff(0)).conjugate(); |
ykuroda | 0:13a5d365ba16 | 121 | } |
ykuroda | 0:13a5d365ba16 | 122 | return res; |
ykuroda | 0:13a5d365ba16 | 123 | } |
ykuroda | 0:13a5d365ba16 | 124 | |
ykuroda | 0:13a5d365ba16 | 125 | namespace internal { |
ykuroda | 0:13a5d365ba16 | 126 | |
ykuroda | 0:13a5d365ba16 | 127 | template<typename Derived, int Size = Derived::SizeAtCompileTime> |
ykuroda | 0:13a5d365ba16 | 128 | struct unitOrthogonal_selector |
ykuroda | 0:13a5d365ba16 | 129 | { |
ykuroda | 0:13a5d365ba16 | 130 | typedef typename plain_matrix_type<Derived>::type VectorType; |
ykuroda | 0:13a5d365ba16 | 131 | typedef typename traits<Derived>::Scalar Scalar; |
ykuroda | 0:13a5d365ba16 | 132 | typedef typename NumTraits<Scalar>::Real RealScalar; |
ykuroda | 0:13a5d365ba16 | 133 | typedef typename Derived::Index Index; |
ykuroda | 0:13a5d365ba16 | 134 | typedef Matrix<Scalar,2,1> Vector2; |
ykuroda | 0:13a5d365ba16 | 135 | static inline VectorType run(const Derived& src) |
ykuroda | 0:13a5d365ba16 | 136 | { |
ykuroda | 0:13a5d365ba16 | 137 | VectorType perp = VectorType::Zero(src.size()); |
ykuroda | 0:13a5d365ba16 | 138 | Index maxi = 0; |
ykuroda | 0:13a5d365ba16 | 139 | Index sndi = 0; |
ykuroda | 0:13a5d365ba16 | 140 | src.cwiseAbs().maxCoeff(&maxi); |
ykuroda | 0:13a5d365ba16 | 141 | if (maxi==0) |
ykuroda | 0:13a5d365ba16 | 142 | sndi = 1; |
ykuroda | 0:13a5d365ba16 | 143 | RealScalar invnm = RealScalar(1)/(Vector2() << src.coeff(sndi),src.coeff(maxi)).finished().norm(); |
ykuroda | 0:13a5d365ba16 | 144 | perp.coeffRef(maxi) = -numext::conj(src.coeff(sndi)) * invnm; |
ykuroda | 0:13a5d365ba16 | 145 | perp.coeffRef(sndi) = numext::conj(src.coeff(maxi)) * invnm; |
ykuroda | 0:13a5d365ba16 | 146 | |
ykuroda | 0:13a5d365ba16 | 147 | return perp; |
ykuroda | 0:13a5d365ba16 | 148 | } |
ykuroda | 0:13a5d365ba16 | 149 | }; |
ykuroda | 0:13a5d365ba16 | 150 | |
ykuroda | 0:13a5d365ba16 | 151 | template<typename Derived> |
ykuroda | 0:13a5d365ba16 | 152 | struct unitOrthogonal_selector<Derived,3> |
ykuroda | 0:13a5d365ba16 | 153 | { |
ykuroda | 0:13a5d365ba16 | 154 | typedef typename plain_matrix_type<Derived>::type VectorType; |
ykuroda | 0:13a5d365ba16 | 155 | typedef typename traits<Derived>::Scalar Scalar; |
ykuroda | 0:13a5d365ba16 | 156 | typedef typename NumTraits<Scalar>::Real RealScalar; |
ykuroda | 0:13a5d365ba16 | 157 | static inline VectorType run(const Derived& src) |
ykuroda | 0:13a5d365ba16 | 158 | { |
ykuroda | 0:13a5d365ba16 | 159 | VectorType perp; |
ykuroda | 0:13a5d365ba16 | 160 | /* Let us compute the crossed product of *this with a vector |
ykuroda | 0:13a5d365ba16 | 161 | * that is not too close to being colinear to *this. |
ykuroda | 0:13a5d365ba16 | 162 | */ |
ykuroda | 0:13a5d365ba16 | 163 | |
ykuroda | 0:13a5d365ba16 | 164 | /* unless the x and y coords are both close to zero, we can |
ykuroda | 0:13a5d365ba16 | 165 | * simply take ( -y, x, 0 ) and normalize it. |
ykuroda | 0:13a5d365ba16 | 166 | */ |
ykuroda | 0:13a5d365ba16 | 167 | if((!isMuchSmallerThan(src.x(), src.z())) |
ykuroda | 0:13a5d365ba16 | 168 | || (!isMuchSmallerThan(src.y(), src.z()))) |
ykuroda | 0:13a5d365ba16 | 169 | { |
ykuroda | 0:13a5d365ba16 | 170 | RealScalar invnm = RealScalar(1)/src.template head<2>().norm(); |
ykuroda | 0:13a5d365ba16 | 171 | perp.coeffRef(0) = -numext::conj(src.y())*invnm; |
ykuroda | 0:13a5d365ba16 | 172 | perp.coeffRef(1) = numext::conj(src.x())*invnm; |
ykuroda | 0:13a5d365ba16 | 173 | perp.coeffRef(2) = 0; |
ykuroda | 0:13a5d365ba16 | 174 | } |
ykuroda | 0:13a5d365ba16 | 175 | /* if both x and y are close to zero, then the vector is close |
ykuroda | 0:13a5d365ba16 | 176 | * to the z-axis, so it's far from colinear to the x-axis for instance. |
ykuroda | 0:13a5d365ba16 | 177 | * So we take the crossed product with (1,0,0) and normalize it. |
ykuroda | 0:13a5d365ba16 | 178 | */ |
ykuroda | 0:13a5d365ba16 | 179 | else |
ykuroda | 0:13a5d365ba16 | 180 | { |
ykuroda | 0:13a5d365ba16 | 181 | RealScalar invnm = RealScalar(1)/src.template tail<2>().norm(); |
ykuroda | 0:13a5d365ba16 | 182 | perp.coeffRef(0) = 0; |
ykuroda | 0:13a5d365ba16 | 183 | perp.coeffRef(1) = -numext::conj(src.z())*invnm; |
ykuroda | 0:13a5d365ba16 | 184 | perp.coeffRef(2) = numext::conj(src.y())*invnm; |
ykuroda | 0:13a5d365ba16 | 185 | } |
ykuroda | 0:13a5d365ba16 | 186 | |
ykuroda | 0:13a5d365ba16 | 187 | return perp; |
ykuroda | 0:13a5d365ba16 | 188 | } |
ykuroda | 0:13a5d365ba16 | 189 | }; |
ykuroda | 0:13a5d365ba16 | 190 | |
ykuroda | 0:13a5d365ba16 | 191 | template<typename Derived> |
ykuroda | 0:13a5d365ba16 | 192 | struct unitOrthogonal_selector<Derived,2> |
ykuroda | 0:13a5d365ba16 | 193 | { |
ykuroda | 0:13a5d365ba16 | 194 | typedef typename plain_matrix_type<Derived>::type VectorType; |
ykuroda | 0:13a5d365ba16 | 195 | static inline VectorType run(const Derived& src) |
ykuroda | 0:13a5d365ba16 | 196 | { return VectorType(-numext::conj(src.y()), numext::conj(src.x())).normalized(); } |
ykuroda | 0:13a5d365ba16 | 197 | }; |
ykuroda | 0:13a5d365ba16 | 198 | |
ykuroda | 0:13a5d365ba16 | 199 | } // end namespace internal |
ykuroda | 0:13a5d365ba16 | 200 | |
ykuroda | 0:13a5d365ba16 | 201 | /** \returns a unit vector which is orthogonal to \c *this |
ykuroda | 0:13a5d365ba16 | 202 | * |
ykuroda | 0:13a5d365ba16 | 203 | * The size of \c *this must be at least 2. If the size is exactly 2, |
ykuroda | 0:13a5d365ba16 | 204 | * then the returned vector is a counter clock wise rotation of \c *this, i.e., (-y,x).normalized(). |
ykuroda | 0:13a5d365ba16 | 205 | * |
ykuroda | 0:13a5d365ba16 | 206 | * \sa cross() |
ykuroda | 0:13a5d365ba16 | 207 | */ |
ykuroda | 0:13a5d365ba16 | 208 | template<typename Derived> |
ykuroda | 0:13a5d365ba16 | 209 | typename MatrixBase<Derived>::PlainObject |
ykuroda | 0:13a5d365ba16 | 210 | MatrixBase<Derived>::unitOrthogonal() const |
ykuroda | 0:13a5d365ba16 | 211 | { |
ykuroda | 0:13a5d365ba16 | 212 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
ykuroda | 0:13a5d365ba16 | 213 | return internal::unitOrthogonal_selector<Derived>::run(derived()); |
ykuroda | 0:13a5d365ba16 | 214 | } |
ykuroda | 0:13a5d365ba16 | 215 | |
ykuroda | 0:13a5d365ba16 | 216 | } // end namespace Eigen |
ykuroda | 0:13a5d365ba16 | 217 | |
ykuroda | 0:13a5d365ba16 | 218 | #endif // EIGEN_ORTHOMETHODS_H |