Eigen  5.0.1-dev+60122df6
 
Loading...
Searching...
No Matches
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// IWYU pragma: private
15#include "./InternalHeaderCheck.h"
16
17namespace Eigen {
18
19namespace internal {
20
21// Vector3 version (default)
22template <typename Derived, typename OtherDerived, int Size>
23struct cross_impl {
24 typedef typename ScalarBinaryOpTraits<typename internal::traits<Derived>::Scalar,
25 typename internal::traits<OtherDerived>::Scalar>::ReturnType Scalar;
26 typedef Matrix<Scalar, MatrixBase<Derived>::RowsAtCompileTime, MatrixBase<Derived>::ColsAtCompileTime> return_type;
27
28 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE return_type run(const MatrixBase<Derived>& first,
29 const MatrixBase<OtherDerived>& second) {
30 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived, 3)
31 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, 3)
32
33 // Note that there is no need for an expression here since the compiler
34 // optimize such a small temporary very well (even within a complex expression)
35 typename internal::nested_eval<Derived, 2>::type lhs(first.derived());
36 typename internal::nested_eval<OtherDerived, 2>::type rhs(second.derived());
37 return return_type(numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),
38 numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),
39 numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0)));
40 }
41};
42
43// Vector2 version
44template <typename Derived, typename OtherDerived>
45struct cross_impl<Derived, OtherDerived, 2> {
46 typedef typename ScalarBinaryOpTraits<typename internal::traits<Derived>::Scalar,
47 typename internal::traits<OtherDerived>::Scalar>::ReturnType Scalar;
48 typedef Scalar return_type;
49
50 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE return_type run(const MatrixBase<Derived>& first,
51 const MatrixBase<OtherDerived>& second) {
52 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived, 2);
53 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, 2);
54 typename internal::nested_eval<Derived, 2>::type lhs(first.derived());
55 typename internal::nested_eval<OtherDerived, 2>::type rhs(second.derived());
56 return numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0));
57 }
58};
59
60} // end namespace internal
61
88template <typename Derived>
89template <typename OtherDerived>
90EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE typename internal::cross_impl<Derived, OtherDerived>::return_type
91MatrixBase<Derived>::cross(const MatrixBase<OtherDerived>& other) const {
92 return internal::cross_impl<Derived, OtherDerived>::run(*this, other);
93}
94
95namespace internal {
96
97template <int Arch, typename VectorLhs, typename VectorRhs, typename Scalar = typename VectorLhs::Scalar,
98 bool Vectorizable =
99 bool((int(evaluator<VectorLhs>::Flags) & int(evaluator<VectorRhs>::Flags)) & PacketAccessBit)>
100struct cross3_impl {
101 EIGEN_DEVICE_FUNC static inline typename internal::plain_matrix_type<VectorLhs>::type run(const VectorLhs& lhs,
102 const VectorRhs& rhs) {
103 return typename internal::plain_matrix_type<VectorLhs>::type(
104 numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),
105 numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),
106 numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0)), 0);
107 }
108};
109
110} // namespace internal
111
121template <typename Derived>
122template <typename OtherDerived>
123EIGEN_DEVICE_FUNC inline typename MatrixBase<Derived>::PlainObject MatrixBase<Derived>::cross3(
124 const MatrixBase<OtherDerived>& other) const {
125 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived, 4)
126 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, 4)
127
128 typedef typename internal::nested_eval<Derived, 2>::type DerivedNested;
129 typedef typename internal::nested_eval<OtherDerived, 2>::type OtherDerivedNested;
130 DerivedNested lhs(derived());
131 OtherDerivedNested rhs(other.derived());
132
133 return internal::cross3_impl<Architecture::Target, internal::remove_all_t<DerivedNested>,
134 internal::remove_all_t<OtherDerivedNested>>::run(lhs, rhs);
135}
136
146template <typename ExpressionType, int Direction>
147template <typename OtherDerived>
148EIGEN_DEVICE_FUNC const typename VectorwiseOp<ExpressionType, Direction>::CrossReturnType
150 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, 3)
151 EIGEN_STATIC_ASSERT(
152 (internal::is_same<Scalar, typename OtherDerived::Scalar>::value),
153 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
154
155 typename internal::nested_eval<ExpressionType, 2>::type mat(_expression());
156 typename internal::nested_eval<OtherDerived, 2>::type vec(other.derived());
157
158 CrossReturnType res(_expression().rows(), _expression().cols());
159 if (Direction == Vertical) {
160 eigen_assert(CrossReturnType::RowsAtCompileTime == 3 && "the matrix must have exactly 3 rows");
161 res.row(0) = (mat.row(1) * vec.coeff(2) - mat.row(2) * vec.coeff(1)).conjugate();
162 res.row(1) = (mat.row(2) * vec.coeff(0) - mat.row(0) * vec.coeff(2)).conjugate();
163 res.row(2) = (mat.row(0) * vec.coeff(1) - mat.row(1) * vec.coeff(0)).conjugate();
164 } else {
165 eigen_assert(CrossReturnType::ColsAtCompileTime == 3 && "the matrix must have exactly 3 columns");
166 res.col(0) = (mat.col(1) * vec.coeff(2) - mat.col(2) * vec.coeff(1)).conjugate();
167 res.col(1) = (mat.col(2) * vec.coeff(0) - mat.col(0) * vec.coeff(2)).conjugate();
168 res.col(2) = (mat.col(0) * vec.coeff(1) - mat.col(1) * vec.coeff(0)).conjugate();
169 }
170 return res;
171}
172
173namespace internal {
174
175template <typename Derived, int Size = Derived::SizeAtCompileTime>
176struct unitOrthogonal_selector {
177 typedef typename plain_matrix_type<Derived>::type VectorType;
178 typedef typename traits<Derived>::Scalar Scalar;
179 typedef typename NumTraits<Scalar>::Real RealScalar;
181 EIGEN_DEVICE_FUNC static inline VectorType run(const Derived& src) {
182 VectorType perp = VectorType::Zero(src.size());
183 Index maxi = 0;
184 Index sndi = 0;
185 src.cwiseAbs().maxCoeff(&maxi);
186 if (maxi == 0) sndi = 1;
187 RealScalar invnm = RealScalar(1) / (Vector2() << src.coeff(sndi), src.coeff(maxi)).finished().norm();
188 perp.coeffRef(maxi) = -numext::conj(src.coeff(sndi)) * invnm;
189 perp.coeffRef(sndi) = numext::conj(src.coeff(maxi)) * invnm;
190
191 return perp;
192 }
193};
194
195template <typename Derived>
196struct unitOrthogonal_selector<Derived, 3> {
197 typedef typename plain_matrix_type<Derived>::type VectorType;
198 typedef typename traits<Derived>::Scalar Scalar;
199 typedef typename NumTraits<Scalar>::Real RealScalar;
200 EIGEN_DEVICE_FUNC static inline VectorType run(const Derived& src) {
201 VectorType perp;
202 /* Let us compute the crossed product of *this with a vector
203 * that is not too close to being colinear to *this.
204 */
205
206 /* unless the x and y coords are both close to zero, we can
207 * simply take ( -y, x, 0 ) and normalize it.
208 */
209 if ((!isMuchSmallerThan(src.x(), src.z())) || (!isMuchSmallerThan(src.y(), src.z()))) {
210 RealScalar invnm = RealScalar(1) / src.template head<2>().norm();
211 perp.coeffRef(0) = -numext::conj(src.y()) * invnm;
212 perp.coeffRef(1) = numext::conj(src.x()) * invnm;
213 perp.coeffRef(2) = 0;
214 }
215 /* if both x and y are close to zero, then the vector is close
216 * to the z-axis, so it's far from colinear to the x-axis for instance.
217 * So we take the crossed product with (1,0,0) and normalize it.
218 */
219 else {
220 RealScalar invnm = RealScalar(1) / src.template tail<2>().norm();
221 perp.coeffRef(0) = 0;
222 perp.coeffRef(1) = -numext::conj(src.z()) * invnm;
223 perp.coeffRef(2) = numext::conj(src.y()) * invnm;
224 }
225
226 return perp;
227 }
228};
229
230template <typename Derived>
231struct unitOrthogonal_selector<Derived, 2> {
232 typedef typename plain_matrix_type<Derived>::type VectorType;
233 EIGEN_DEVICE_FUNC static inline VectorType run(const Derived& src) {
234 return VectorType(-numext::conj(src.y()), numext::conj(src.x())).normalized();
235 }
236};
237
238} // end namespace internal
239
249template <typename Derived>
250EIGEN_DEVICE_FUNC typename MatrixBase<Derived>::PlainObject MatrixBase<Derived>::unitOrthogonal() const {
251 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
252 return internal::unitOrthogonal_selector<Derived>::run(derived());
253}
254
255} // end namespace Eigen
256
257#endif // EIGEN_ORTHOMETHODS_H
@ ColsAtCompileTime
Definition DenseBase.h:102
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:186
const CrossReturnType cross(const MatrixBase< OtherDerived > &other) const
Definition OrthoMethods.h:149
PlainObject unitOrthogonal(void) const
Definition OrthoMethods.h:250
PlainObject cross3(const MatrixBase< OtherDerived > &other) const
Definition OrthoMethods.h:123
internal::cross_impl< Derived, OtherDerived >::return_type cross(const MatrixBase< OtherDerived > &other) const
Definition OrthoMethods.h:91
@ Vertical
Definition Constants.h:266
const unsigned int PacketAccessBit
Definition Constants.h:97
Matrix< Type, 2, 1 > Vector2
[c++11] 2×1 vector of type Type.
Definition Matrix.h:511
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82