Eigen  3.2.93
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 namespace Eigen {
15 
27 template<typename Derived>
28 template<typename OtherDerived>
29 #ifndef EIGEN_PARSED_BY_DOXYGEN
30 inline typename MatrixBase<Derived>::template cross_product_return_type<OtherDerived>::type
31 #else
32 inline typename MatrixBase<Derived>::PlainObject
33 #endif
35 {
36  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,3)
37  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3)
38 
39  // Note that there is no need for an expression here since the compiler
40  // optimize such a small temporary very well (even within a complex expression)
41  typename internal::nested_eval<Derived,2>::type lhs(derived());
42  typename internal::nested_eval<OtherDerived,2>::type rhs(other.derived());
43  return typename cross_product_return_type<OtherDerived>::type(
44  numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),
45  numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),
46  numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0))
47  );
48 }
49 
50 namespace internal {
51 
52 template< int Arch,typename VectorLhs,typename VectorRhs,
53  typename Scalar = typename VectorLhs::Scalar,
54  bool Vectorizable = bool((VectorLhs::Flags&VectorRhs::Flags)&PacketAccessBit)>
55 struct cross3_impl {
56  static inline typename internal::plain_matrix_type<VectorLhs>::type
57  run(const VectorLhs& lhs, const VectorRhs& rhs)
58  {
59  return typename internal::plain_matrix_type<VectorLhs>::type(
60  numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),
61  numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),
62  numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0)),
63  0
64  );
65  }
66 };
67 
68 }
69 
79 template<typename Derived>
80 template<typename OtherDerived>
81 inline typename MatrixBase<Derived>::PlainObject
83 {
84  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,4)
85  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,4)
86 
87  typedef typename internal::nested_eval<Derived,2>::type DerivedNested;
88  typedef typename internal::nested_eval<OtherDerived,2>::type OtherDerivedNested;
89  DerivedNested lhs(derived());
90  OtherDerivedNested rhs(other.derived());
91 
92  return internal::cross3_impl<Architecture::Target,
93  typename internal::remove_all<DerivedNested>::type,
94  typename internal::remove_all<OtherDerivedNested>::type>::run(lhs,rhs);
95 }
96 
106 template<typename ExpressionType, int Direction>
107 template<typename OtherDerived>
108 const typename VectorwiseOp<ExpressionType,Direction>::CrossReturnType
110 {
111  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3)
112  EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename OtherDerived::Scalar>::value),
113  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
114 
115  typename internal::nested_eval<ExpressionType,2>::type mat(_expression());
116  typename internal::nested_eval<OtherDerived,2>::type vec(other.derived());
117 
118  CrossReturnType res(_expression().rows(),_expression().cols());
119  if(Direction==Vertical)
120  {
121  eigen_assert(CrossReturnType::RowsAtCompileTime==3 && "the matrix must have exactly 3 rows");
122  res.row(0) = (mat.row(1) * vec.coeff(2) - mat.row(2) * vec.coeff(1)).conjugate();
123  res.row(1) = (mat.row(2) * vec.coeff(0) - mat.row(0) * vec.coeff(2)).conjugate();
124  res.row(2) = (mat.row(0) * vec.coeff(1) - mat.row(1) * vec.coeff(0)).conjugate();
125  }
126  else
127  {
128  eigen_assert(CrossReturnType::ColsAtCompileTime==3 && "the matrix must have exactly 3 columns");
129  res.col(0) = (mat.col(1) * vec.coeff(2) - mat.col(2) * vec.coeff(1)).conjugate();
130  res.col(1) = (mat.col(2) * vec.coeff(0) - mat.col(0) * vec.coeff(2)).conjugate();
131  res.col(2) = (mat.col(0) * vec.coeff(1) - mat.col(1) * vec.coeff(0)).conjugate();
132  }
133  return res;
134 }
135 
136 namespace internal {
137 
138 template<typename Derived, int Size = Derived::SizeAtCompileTime>
139 struct unitOrthogonal_selector
140 {
141  typedef typename plain_matrix_type<Derived>::type VectorType;
142  typedef typename traits<Derived>::Scalar Scalar;
143  typedef typename NumTraits<Scalar>::Real RealScalar;
144  typedef Matrix<Scalar,2,1> Vector2;
145  EIGEN_DEVICE_FUNC
146  static inline VectorType run(const Derived& src)
147  {
148  VectorType perp = VectorType::Zero(src.size());
149  Index maxi = 0;
150  Index sndi = 0;
151  src.cwiseAbs().maxCoeff(&maxi);
152  if (maxi==0)
153  sndi = 1;
154  RealScalar invnm = RealScalar(1)/(Vector2() << src.coeff(sndi),src.coeff(maxi)).finished().norm();
155  perp.coeffRef(maxi) = -numext::conj(src.coeff(sndi)) * invnm;
156  perp.coeffRef(sndi) = numext::conj(src.coeff(maxi)) * invnm;
157 
158  return perp;
159  }
160 };
161 
162 template<typename Derived>
163 struct unitOrthogonal_selector<Derived,3>
164 {
165  typedef typename plain_matrix_type<Derived>::type VectorType;
166  typedef typename traits<Derived>::Scalar Scalar;
167  typedef typename NumTraits<Scalar>::Real RealScalar;
168  EIGEN_DEVICE_FUNC
169  static inline VectorType run(const Derived& src)
170  {
171  VectorType perp;
172  /* Let us compute the crossed product of *this with a vector
173  * that is not too close to being colinear to *this.
174  */
175 
176  /* unless the x and y coords are both close to zero, we can
177  * simply take ( -y, x, 0 ) and normalize it.
178  */
179  if((!isMuchSmallerThan(src.x(), src.z()))
180  || (!isMuchSmallerThan(src.y(), src.z())))
181  {
182  RealScalar invnm = RealScalar(1)/src.template head<2>().norm();
183  perp.coeffRef(0) = -numext::conj(src.y())*invnm;
184  perp.coeffRef(1) = numext::conj(src.x())*invnm;
185  perp.coeffRef(2) = 0;
186  }
187  /* if both x and y are close to zero, then the vector is close
188  * to the z-axis, so it's far from colinear to the x-axis for instance.
189  * So we take the crossed product with (1,0,0) and normalize it.
190  */
191  else
192  {
193  RealScalar invnm = RealScalar(1)/src.template tail<2>().norm();
194  perp.coeffRef(0) = 0;
195  perp.coeffRef(1) = -numext::conj(src.z())*invnm;
196  perp.coeffRef(2) = numext::conj(src.y())*invnm;
197  }
198 
199  return perp;
200  }
201 };
202 
203 template<typename Derived>
204 struct unitOrthogonal_selector<Derived,2>
205 {
206  typedef typename plain_matrix_type<Derived>::type VectorType;
207  EIGEN_DEVICE_FUNC
208  static inline VectorType run(const Derived& src)
209  { return VectorType(-numext::conj(src.y()), numext::conj(src.x())).normalized(); }
210 };
211 
212 } // end namespace internal
213 
223 template<typename Derived>
224 typename MatrixBase<Derived>::PlainObject
226 {
227  EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
228  return internal::unitOrthogonal_selector<Derived>::run(derived());
229 }
230 
231 } // end namespace Eigen
232 
233 #endif // EIGEN_ORTHOMETHODS_H
PlainObject unitOrthogonal(void) const
Definition: OrthoMethods.h:225
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:66
Definition: Constants.h:265
Eigen::Index Index
Definition: VectorwiseOp.h:162
Namespace containing all symbols from the Eigen library.
Definition: Core:271
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:167
Derived & derived()
Definition: EigenBase.h:44
const unsigned int PacketAccessBit
Definition: Constants.h:89
PlainObject cross3(const MatrixBase< OtherDerived > &other) const
Definition: OrthoMethods.h:82
const CrossReturnType cross(const MatrixBase< OtherDerived > &other) const
Definition: OrthoMethods.h:109
Definition: Eigen_Colamd.h:50
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:178
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
PlainObject cross(const MatrixBase< OtherDerived > &other) const
Definition: OrthoMethods.h:34