Eigne Matrix Class Library
Dependents: Eigen_test Odometry_test AttitudeEstimation_usingTicker MPU9250_Quaternion_Binary_Serial ... more
OrthoMethods.h
00001 // This file is part of Eigen, a lightweight C++ template library 00002 // for linear algebra. 00003 // 00004 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr> 00005 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com> 00006 // 00007 // This Source Code Form is subject to the terms of the Mozilla 00008 // Public License v. 2.0. If a copy of the MPL was not distributed 00009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 00010 00011 #ifndef EIGEN_ORTHOMETHODS_H 00012 #define EIGEN_ORTHOMETHODS_H 00013 00014 namespace Eigen { 00015 00016 /** \geometry_module 00017 * 00018 * \returns the cross product of \c *this and \a other 00019 * 00020 * Here is a very good explanation of cross-product: http://xkcd.com/199/ 00021 * \sa MatrixBase::cross3() 00022 */ 00023 template<typename Derived> 00024 template<typename OtherDerived> 00025 inline typename MatrixBase<Derived>::template cross_product_return_type<OtherDerived>::type 00026 MatrixBase<Derived>::cross (const MatrixBase<OtherDerived>& other) const 00027 { 00028 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,3) 00029 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3) 00030 00031 // Note that there is no need for an expression here since the compiler 00032 // optimize such a small temporary very well (even within a complex expression) 00033 typename internal::nested<Derived,2>::type lhs(derived()); 00034 typename internal::nested<OtherDerived,2>::type rhs(other.derived()); 00035 return typename cross_product_return_type<OtherDerived>::type( 00036 numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)), 00037 numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)), 00038 numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0)) 00039 ); 00040 } 00041 00042 namespace internal { 00043 00044 template< int Arch,typename VectorLhs,typename VectorRhs, 00045 typename Scalar = typename VectorLhs::Scalar, 00046 bool Vectorizable = bool((VectorLhs::Flags&VectorRhs::Flags)&PacketAccessBit)> 00047 struct cross3_impl { 00048 static inline typename internal::plain_matrix_type<VectorLhs>::type 00049 run(const VectorLhs& lhs, const VectorRhs& rhs) 00050 { 00051 return typename internal::plain_matrix_type<VectorLhs>::type( 00052 numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)), 00053 numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)), 00054 numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0)), 00055 0 00056 ); 00057 } 00058 }; 00059 00060 } 00061 00062 /** \geometry_module 00063 * 00064 * \returns the cross product of \c *this and \a other using only the x, y, and z coefficients 00065 * 00066 * The size of \c *this and \a other must be four. This function is especially useful 00067 * when using 4D vectors instead of 3D ones to get advantage of SSE/AltiVec vectorization. 00068 * 00069 * \sa MatrixBase::cross() 00070 */ 00071 template<typename Derived> 00072 template<typename OtherDerived> 00073 inline typename MatrixBase<Derived>::PlainObject 00074 MatrixBase<Derived>::cross3 (const MatrixBase<OtherDerived>& other) const 00075 { 00076 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,4) 00077 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,4) 00078 00079 typedef typename internal::nested<Derived,2>::type DerivedNested; 00080 typedef typename internal::nested<OtherDerived,2>::type OtherDerivedNested; 00081 DerivedNested lhs(derived()); 00082 OtherDerivedNested rhs(other.derived()); 00083 00084 return internal::cross3_impl<Architecture::Target, 00085 typename internal::remove_all<DerivedNested>::type, 00086 typename internal::remove_all<OtherDerivedNested>::type>::run(lhs,rhs); 00087 } 00088 00089 /** \returns a matrix expression of the cross product of each column or row 00090 * of the referenced expression with the \a other vector. 00091 * 00092 * The referenced matrix must have one dimension equal to 3. 00093 * The result matrix has the same dimensions than the referenced one. 00094 * 00095 * \geometry_module 00096 * 00097 * \sa MatrixBase::cross() */ 00098 template<typename ExpressionType, int Direction> 00099 template<typename OtherDerived> 00100 const typename VectorwiseOp<ExpressionType,Direction>::CrossReturnType 00101 VectorwiseOp<ExpressionType,Direction>::cross (const MatrixBase<OtherDerived>& other) const 00102 { 00103 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3) 00104 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename OtherDerived::Scalar>::value), 00105 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY) 00106 00107 CrossReturnType res(_expression().rows(),_expression().cols()); 00108 if(Direction==Vertical) 00109 { 00110 eigen_assert(CrossReturnType::RowsAtCompileTime==3 && "the matrix must have exactly 3 rows"); 00111 res.row(0) = (_expression().row(1) * other.coeff(2) - _expression().row(2) * other.coeff(1)).conjugate(); 00112 res.row(1) = (_expression().row(2) * other.coeff(0) - _expression().row(0) * other.coeff(2)).conjugate(); 00113 res.row(2) = (_expression().row(0) * other.coeff(1) - _expression().row(1) * other.coeff(0)).conjugate(); 00114 } 00115 else 00116 { 00117 eigen_assert(CrossReturnType::ColsAtCompileTime==3 && "the matrix must have exactly 3 columns"); 00118 res.col(0) = (_expression().col(1) * other.coeff(2) - _expression().col(2) * other.coeff(1)).conjugate(); 00119 res.col(1) = (_expression().col(2) * other.coeff(0) - _expression().col(0) * other.coeff(2)).conjugate(); 00120 res.col(2) = (_expression().col(0) * other.coeff(1) - _expression().col(1) * other.coeff(0)).conjugate(); 00121 } 00122 return res; 00123 } 00124 00125 namespace internal { 00126 00127 template<typename Derived, int Size = Derived::SizeAtCompileTime> 00128 struct unitOrthogonal_selector 00129 { 00130 typedef typename plain_matrix_type<Derived>::type VectorType; 00131 typedef typename traits<Derived>::Scalar Scalar; 00132 typedef typename NumTraits<Scalar>::Real RealScalar; 00133 typedef typename Derived::Index Index; 00134 typedef Matrix<Scalar,2,1> Vector2; 00135 static inline VectorType run(const Derived& src) 00136 { 00137 VectorType perp = VectorType::Zero(src.size()); 00138 Index maxi = 0; 00139 Index sndi = 0; 00140 src.cwiseAbs().maxCoeff(&maxi); 00141 if (maxi==0) 00142 sndi = 1; 00143 RealScalar invnm = RealScalar(1)/(Vector2() << src.coeff(sndi),src.coeff(maxi)).finished().norm(); 00144 perp.coeffRef(maxi) = -numext::conj(src.coeff(sndi)) * invnm; 00145 perp.coeffRef(sndi) = numext::conj(src.coeff(maxi)) * invnm; 00146 00147 return perp; 00148 } 00149 }; 00150 00151 template<typename Derived> 00152 struct unitOrthogonal_selector<Derived,3> 00153 { 00154 typedef typename plain_matrix_type<Derived>::type VectorType; 00155 typedef typename traits<Derived>::Scalar Scalar; 00156 typedef typename NumTraits<Scalar>::Real RealScalar; 00157 static inline VectorType run(const Derived& src) 00158 { 00159 VectorType perp; 00160 /* Let us compute the crossed product of *this with a vector 00161 * that is not too close to being colinear to *this. 00162 */ 00163 00164 /* unless the x and y coords are both close to zero, we can 00165 * simply take ( -y, x, 0 ) and normalize it. 00166 */ 00167 if((!isMuchSmallerThan(src.x(), src.z())) 00168 || (!isMuchSmallerThan(src.y(), src.z()))) 00169 { 00170 RealScalar invnm = RealScalar(1)/src.template head<2>().norm(); 00171 perp.coeffRef(0) = -numext::conj(src.y())*invnm; 00172 perp.coeffRef(1) = numext::conj(src.x())*invnm; 00173 perp.coeffRef(2) = 0; 00174 } 00175 /* if both x and y are close to zero, then the vector is close 00176 * to the z-axis, so it's far from colinear to the x-axis for instance. 00177 * So we take the crossed product with (1,0,0) and normalize it. 00178 */ 00179 else 00180 { 00181 RealScalar invnm = RealScalar(1)/src.template tail<2>().norm(); 00182 perp.coeffRef(0) = 0; 00183 perp.coeffRef(1) = -numext::conj(src.z())*invnm; 00184 perp.coeffRef(2) = numext::conj(src.y())*invnm; 00185 } 00186 00187 return perp; 00188 } 00189 }; 00190 00191 template<typename Derived> 00192 struct unitOrthogonal_selector<Derived,2> 00193 { 00194 typedef typename plain_matrix_type<Derived>::type VectorType; 00195 static inline VectorType run(const Derived& src) 00196 { return VectorType(-numext::conj(src.y()), numext::conj(src.x())).normalized(); } 00197 }; 00198 00199 } // end namespace internal 00200 00201 /** \returns a unit vector which is orthogonal to \c *this 00202 * 00203 * The size of \c *this must be at least 2. If the size is exactly 2, 00204 * then the returned vector is a counter clock wise rotation of \c *this, i.e., (-y,x).normalized(). 00205 * 00206 * \sa cross() 00207 */ 00208 template<typename Derived> 00209 typename MatrixBase<Derived>::PlainObject 00210 MatrixBase<Derived>::unitOrthogonal () const 00211 { 00212 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) 00213 return internal::unitOrthogonal_selector<Derived>::run(derived()); 00214 } 00215 00216 } // end namespace Eigen 00217 00218 #endif // EIGEN_ORTHOMETHODS_H
Generated on Tue Jul 12 2022 17:46:58 by 1.7.2