Eigen  3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
OrthoMethods.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_ORTHOMETHODS_H
12 #define EIGEN_ORTHOMETHODS_H
13 
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
29 template<typename Derived>
30 template<typename OtherDerived>
31 #ifndef EIGEN_PARSED_BY_DOXYGEN
32 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
33 typename MatrixBase<Derived>::template cross_product_return_type<OtherDerived>::type
34 #else
35 typename MatrixBase<Derived>::PlainObject
36 #endif
38 {
39  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,3)
40  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3)
41 
42  // Note that there is no need for an expression here since the compiler
43  // optimize such a small temporary very well (even within a complex expression)
44  typename internal::nested_eval<Derived,2>::type lhs(derived());
45  typename internal::nested_eval<OtherDerived,2>::type rhs(other.derived());
46  return typename cross_product_return_type<OtherDerived>::type(
47  numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),
48  numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),
49  numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0))
50  );
51 }
52 
53 namespace internal {
54 
55 template< int Arch,typename VectorLhs,typename VectorRhs,
56  typename Scalar = typename VectorLhs::Scalar,
57  bool Vectorizable = bool((VectorLhs::Flags&VectorRhs::Flags)&PacketAccessBit)>
58 struct cross3_impl {
59  EIGEN_DEVICE_FUNC static inline typename internal::plain_matrix_type<VectorLhs>::type
60  run(const VectorLhs& lhs, const VectorRhs& rhs)
61  {
62  return typename internal::plain_matrix_type<VectorLhs>::type(
63  numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),
64  numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),
65  numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0)),
66  0
67  );
68  }
69 };
70 
71 }
72 
82 template<typename Derived>
83 template<typename OtherDerived>
84 EIGEN_DEVICE_FUNC inline typename MatrixBase<Derived>::PlainObject
86 {
87  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,4)
88  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,4)
89 
90  typedef typename internal::nested_eval<Derived,2>::type DerivedNested;
91  typedef typename internal::nested_eval<OtherDerived,2>::type OtherDerivedNested;
92  DerivedNested lhs(derived());
93  OtherDerivedNested rhs(other.derived());
94 
95  return internal::cross3_impl<Architecture::Target,
96  internal::remove_all_t<DerivedNested>,
97  internal::remove_all_t<OtherDerivedNested>>::run(lhs,rhs);
98 }
99 
109 template<typename ExpressionType, int Direction>
110 template<typename OtherDerived>
111 EIGEN_DEVICE_FUNC
112 const typename VectorwiseOp<ExpressionType,Direction>::CrossReturnType
114 {
115  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3)
116  EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename OtherDerived::Scalar>::value),
117  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
118 
119  typename internal::nested_eval<ExpressionType,2>::type mat(_expression());
120  typename internal::nested_eval<OtherDerived,2>::type vec(other.derived());
121 
122  CrossReturnType res(_expression().rows(),_expression().cols());
123  if(Direction==Vertical)
124  {
125  eigen_assert(CrossReturnType::RowsAtCompileTime==3 && "the matrix must have exactly 3 rows");
126  res.row(0) = (mat.row(1) * vec.coeff(2) - mat.row(2) * vec.coeff(1)).conjugate();
127  res.row(1) = (mat.row(2) * vec.coeff(0) - mat.row(0) * vec.coeff(2)).conjugate();
128  res.row(2) = (mat.row(0) * vec.coeff(1) - mat.row(1) * vec.coeff(0)).conjugate();
129  }
130  else
131  {
132  eigen_assert(CrossReturnType::ColsAtCompileTime==3 && "the matrix must have exactly 3 columns");
133  res.col(0) = (mat.col(1) * vec.coeff(2) - mat.col(2) * vec.coeff(1)).conjugate();
134  res.col(1) = (mat.col(2) * vec.coeff(0) - mat.col(0) * vec.coeff(2)).conjugate();
135  res.col(2) = (mat.col(0) * vec.coeff(1) - mat.col(1) * vec.coeff(0)).conjugate();
136  }
137  return res;
138 }
139 
140 namespace internal {
141 
142 template<typename Derived, int Size = Derived::SizeAtCompileTime>
143 struct unitOrthogonal_selector
144 {
145  typedef typename plain_matrix_type<Derived>::type VectorType;
146  typedef typename traits<Derived>::Scalar Scalar;
147  typedef typename NumTraits<Scalar>::Real RealScalar;
148  typedef Matrix<Scalar,2,1> Vector2;
149  EIGEN_DEVICE_FUNC
150  static inline VectorType run(const Derived& src)
151  {
152  VectorType perp = VectorType::Zero(src.size());
153  Index maxi = 0;
154  Index sndi = 0;
155  src.cwiseAbs().maxCoeff(&maxi);
156  if (maxi==0)
157  sndi = 1;
158  RealScalar invnm = RealScalar(1)/(Vector2() << src.coeff(sndi),src.coeff(maxi)).finished().norm();
159  perp.coeffRef(maxi) = -numext::conj(src.coeff(sndi)) * invnm;
160  perp.coeffRef(sndi) = numext::conj(src.coeff(maxi)) * invnm;
161 
162  return perp;
163  }
164 };
165 
166 template<typename Derived>
167 struct unitOrthogonal_selector<Derived,3>
168 {
169  typedef typename plain_matrix_type<Derived>::type VectorType;
170  typedef typename traits<Derived>::Scalar Scalar;
171  typedef typename NumTraits<Scalar>::Real RealScalar;
172  EIGEN_DEVICE_FUNC
173  static inline VectorType run(const Derived& src)
174  {
175  VectorType perp;
176  /* Let us compute the crossed product of *this with a vector
177  * that is not too close to being colinear to *this.
178  */
179 
180  /* unless the x and y coords are both close to zero, we can
181  * simply take ( -y, x, 0 ) and normalize it.
182  */
183  if((!isMuchSmallerThan(src.x(), src.z()))
184  || (!isMuchSmallerThan(src.y(), src.z())))
185  {
186  RealScalar invnm = RealScalar(1)/src.template head<2>().norm();
187  perp.coeffRef(0) = -numext::conj(src.y())*invnm;
188  perp.coeffRef(1) = numext::conj(src.x())*invnm;
189  perp.coeffRef(2) = 0;
190  }
191  /* if both x and y are close to zero, then the vector is close
192  * to the z-axis, so it's far from colinear to the x-axis for instance.
193  * So we take the crossed product with (1,0,0) and normalize it.
194  */
195  else
196  {
197  RealScalar invnm = RealScalar(1)/src.template tail<2>().norm();
198  perp.coeffRef(0) = 0;
199  perp.coeffRef(1) = -numext::conj(src.z())*invnm;
200  perp.coeffRef(2) = numext::conj(src.y())*invnm;
201  }
202 
203  return perp;
204  }
205 };
206 
207 template<typename Derived>
208 struct unitOrthogonal_selector<Derived,2>
209 {
210  typedef typename plain_matrix_type<Derived>::type VectorType;
211  EIGEN_DEVICE_FUNC
212  static inline VectorType run(const Derived& src)
213  { return VectorType(-numext::conj(src.y()), numext::conj(src.x())).normalized(); }
214 };
215 
216 } // end namespace internal
217 
227 template<typename Derived>
228 EIGEN_DEVICE_FUNC typename MatrixBase<Derived>::PlainObject
230 {
231  EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
232  return internal::unitOrthogonal_selector<Derived>::run(derived());
233 }
234 
235 } // end namespace Eigen
236 
237 #endif // EIGEN_ORTHOMETHODS_H
Derived & derived()
Definition: EigenBase.h:48
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182
PlainObject cross(const MatrixBase< OtherDerived > &other) const
Definition: OrthoMethods.h:37
const CrossReturnType cross(const MatrixBase< OtherDerived > &other) const
Definition: OrthoMethods.h:113
PlainObject unitOrthogonal(void) const
Definition: OrthoMethods.h:229
PlainObject cross3(const MatrixBase< OtherDerived > &other) const
Definition: OrthoMethods.h:85
@ Vertical
Definition: Constants.h:266
const unsigned int PacketAccessBit
Definition: Constants.h:96
Matrix< Type, 2, 1 > Vector2
[c++11] 2×1 vector of type Type.
Definition: Matrix.h:533
Namespace containing all symbols from the Eigen library.
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231