OrthoMethods.h
Go to the documentation of this file.
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 
23 template<typename Derived>
24 template<typename OtherDerived>
25 inline typename MatrixBase<Derived>::template cross_product_return_type<OtherDerived>::type
27 {
30 
31  // Note that there is no need for an expression here since the compiler
32  // optimize such a small temporary very well (even within a complex expression)
33  typename internal::nested<Derived,2>::type lhs(derived());
34  typename internal::nested<OtherDerived,2>::type rhs(other.derived());
36  numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),
37  numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),
38  numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0))
39  );
40 }
41 
42 namespace internal {
43 
44 template< int Arch,typename VectorLhs,typename VectorRhs,
45  typename Scalar = typename VectorLhs::Scalar,
46  bool Vectorizable = bool((VectorLhs::Flags&VectorRhs::Flags)&PacketAccessBit)>
47 struct cross3_impl {
49  run(const VectorLhs& lhs, const VectorRhs& rhs)
50  {
52  numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),
53  numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),
54  numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0)),
55  0
56  );
57  }
58 };
59 
60 }
61 
71 template<typename Derived>
72 template<typename OtherDerived>
75 {
78 
79  typedef typename internal::nested<Derived,2>::type DerivedNested;
80  typedef typename internal::nested<OtherDerived,2>::type OtherDerivedNested;
81  DerivedNested lhs(derived());
82  OtherDerivedNested rhs(other.derived());
83 
87 }
88 
98 template<typename ExpressionType, int Direction>
99 template<typename OtherDerived>
102 {
105  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
106 
107  CrossReturnType res(_expression().rows(),_expression().cols());
108  if(Direction==Vertical)
109  {
110  eigen_assert(CrossReturnType::RowsAtCompileTime==3 && "the matrix must have exactly 3 rows");
111  res.row(0) = (_expression().row(1) * other.coeff(2) - _expression().row(2) * other.coeff(1)).conjugate();
112  res.row(1) = (_expression().row(2) * other.coeff(0) - _expression().row(0) * other.coeff(2)).conjugate();
113  res.row(2) = (_expression().row(0) * other.coeff(1) - _expression().row(1) * other.coeff(0)).conjugate();
114  }
115  else
116  {
117  eigen_assert(CrossReturnType::ColsAtCompileTime==3 && "the matrix must have exactly 3 columns");
118  res.col(0) = (_expression().col(1) * other.coeff(2) - _expression().col(2) * other.coeff(1)).conjugate();
119  res.col(1) = (_expression().col(2) * other.coeff(0) - _expression().col(0) * other.coeff(2)).conjugate();
120  res.col(2) = (_expression().col(0) * other.coeff(1) - _expression().col(1) * other.coeff(0)).conjugate();
121  }
122  return res;
123 }
124 
125 namespace internal {
126 
127 template<typename Derived, int Size = Derived::SizeAtCompileTime>
129 {
131  typedef typename traits<Derived>::Scalar Scalar;
133  typedef typename Derived::Index Index;
135  static inline VectorType run(const Derived& src)
136  {
137  VectorType perp = VectorType::Zero(src.size());
138  Index maxi = 0;
139  Index sndi = 0;
140  src.cwiseAbs().maxCoeff(&maxi);
141  if (maxi==0)
142  sndi = 1;
143  RealScalar invnm = RealScalar(1)/(Vector2() << src.coeff(sndi),src.coeff(maxi)).finished().norm();
144  perp.coeffRef(maxi) = -numext::conj(src.coeff(sndi)) * invnm;
145  perp.coeffRef(sndi) = numext::conj(src.coeff(maxi)) * invnm;
146 
147  return perp;
148  }
149 };
150 
151 template<typename Derived>
152 struct unitOrthogonal_selector<Derived,3>
153 {
155  typedef typename traits<Derived>::Scalar Scalar;
157  static inline VectorType run(const Derived& src)
158  {
159  VectorType perp;
160  /* Let us compute the crossed product of *this with a vector
161  * that is not too close to being colinear to *this.
162  */
163 
164  /* unless the x and y coords are both close to zero, we can
165  * simply take ( -y, x, 0 ) and normalize it.
166  */
167  if((!isMuchSmallerThan(src.x(), src.z()))
168  || (!isMuchSmallerThan(src.y(), src.z())))
169  {
170  RealScalar invnm = RealScalar(1)/src.template head<2>().norm();
171  perp.coeffRef(0) = -numext::conj(src.y())*invnm;
172  perp.coeffRef(1) = numext::conj(src.x())*invnm;
173  perp.coeffRef(2) = 0;
174  }
175  /* if both x and y are close to zero, then the vector is close
176  * to the z-axis, so it's far from colinear to the x-axis for instance.
177  * So we take the crossed product with (1,0,0) and normalize it.
178  */
179  else
180  {
181  RealScalar invnm = RealScalar(1)/src.template tail<2>().norm();
182  perp.coeffRef(0) = 0;
183  perp.coeffRef(1) = -numext::conj(src.z())*invnm;
184  perp.coeffRef(2) = numext::conj(src.y())*invnm;
185  }
186 
187  return perp;
188  }
189 };
190 
191 template<typename Derived>
192 struct unitOrthogonal_selector<Derived,2>
193 {
195  static inline VectorType run(const Derived& src)
196  { return VectorType(-numext::conj(src.y()), numext::conj(src.x())).normalized(); }
197 };
198 
199 } // end namespace internal
200 
208 template<typename Derived>
211 {
214 }
215 
216 } // end namespace Eigen
217 
218 #endif // EIGEN_ORTHOMETHODS_H
ColXpr col(Index i)
Definition: DenseBase.h:709
PlainObject unitOrthogonal(void) const
Definition: OrthoMethods.h:210
const CrossReturnType cross(const MatrixBase< OtherDerived > &other) const
Definition: OrthoMethods.h:101
NumTraits< Scalar >::Real RealScalar
Definition: OrthoMethods.h:132
static internal::plain_matrix_type< VectorLhs >::type run(const VectorLhs &lhs, const VectorRhs &rhs)
Definition: OrthoMethods.h:49
traits< Derived >::Scalar Scalar
Definition: OrthoMethods.h:131
plain_matrix_type< Derived >::type VectorType
Definition: OrthoMethods.h:154
Definition: LDLT.h:16
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
RowXpr row(Index i)
Definition: DenseBase.h:726
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Definition: StaticAssert.h:111
bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, typename NumTraits< Scalar >::Real precision=NumTraits< Scalar >::dummy_precision())
internal::traits< Derived >::Scalar Scalar
Definition: MatrixBase.h:56
const unsigned int PacketAccessBit
Definition: Constants.h:81
PlainObject cross3(const MatrixBase< OtherDerived > &other) const
Definition: OrthoMethods.h:74
ConjugateReturnType conjugate() const
plain_matrix_type< Derived >::type VectorType
Definition: OrthoMethods.h:130
static VectorType run(const Derived &src)
Definition: OrthoMethods.h:135
cross_product_return_type< OtherDerived >::type cross(const MatrixBase< OtherDerived > &other) const
static VectorType run(const Derived &src)
Definition: OrthoMethods.h:157
plain_matrix_type< Derived >::type VectorType
Definition: OrthoMethods.h:194
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:127
#define eigen_assert(x)
#define EIGEN_STATIC_ASSERT_VECTOR_ONLY(TYPE)
Definition: StaticAssert.h:126
ExpressionType::PlainObject CrossReturnType
Definition: VectorwiseOp.h:557
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
static VectorType run(const Derived &src)
Definition: OrthoMethods.h:195
#define EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(TYPE, SIZE)
Definition: StaticAssert.h:141


tuw_aruco
Author(s): Lukas Pfeifhofer
autogenerated on Mon Jun 10 2019 15:40:54