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