Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
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 <[email protected]>
5// Copyright (C) 2006-2008 Benoit Jacob <[email protected]>
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
86template <typename Derived>
87template <typename OtherDerived>
88EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
89#ifndef EIGEN_PARSED_BY_DOXYGEN
90 typename internal::cross_impl<Derived, OtherDerived>::return_type
91#else
92 inline std::conditional_t<SizeAtCompileTime == 2, Scalar, PlainObject>
93#endif
95 return internal::cross_impl<Derived, OtherDerived>::run(*this, other);
96}
97
98namespace internal {
99
100template <int Arch, typename VectorLhs, typename VectorRhs, typename Scalar = typename VectorLhs::Scalar,
101 bool Vectorizable = bool((evaluator<VectorLhs>::Flags & evaluator<VectorRhs>::Flags) & PacketAccessBit)>
102struct cross3_impl {
103 EIGEN_DEVICE_FUNC static inline typename internal::plain_matrix_type<VectorLhs>::type run(const VectorLhs& lhs,
104 const VectorRhs& rhs) {
105 return typename internal::plain_matrix_type<VectorLhs>::type(
106 numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),
107 numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),
108 numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0)), 0);
109 }
110};
111
112} // namespace internal
113
123template <typename Derived>
124template <typename OtherDerived>
125EIGEN_DEVICE_FUNC inline typename MatrixBase<Derived>::PlainObject MatrixBase<Derived>::cross3(
126 const MatrixBase<OtherDerived>& other) const {
127 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived, 4)
128 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, 4)
129
130 typedef typename internal::nested_eval<Derived, 2>::type DerivedNested;
131 typedef typename internal::nested_eval<OtherDerived, 2>::type OtherDerivedNested;
132 DerivedNested lhs(derived());
133 OtherDerivedNested rhs(other.derived());
134
135 return internal::cross3_impl<Architecture::Target, internal::remove_all_t<DerivedNested>,
136 internal::remove_all_t<OtherDerivedNested>>::run(lhs, rhs);
137}
138
148template <typename ExpressionType, int Direction>
149template <typename OtherDerived>
150EIGEN_DEVICE_FUNC const typename VectorwiseOp<ExpressionType, Direction>::CrossReturnType
152 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, 3)
153 EIGEN_STATIC_ASSERT(
154 (internal::is_same<Scalar, typename OtherDerived::Scalar>::value),
155 YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
156
157 typename internal::nested_eval<ExpressionType, 2>::type mat(_expression());
158 typename internal::nested_eval<OtherDerived, 2>::type vec(other.derived());
159
160 CrossReturnType res(_expression().rows(), _expression().cols());
161 if (Direction == Vertical) {
162 eigen_assert(CrossReturnType::RowsAtCompileTime == 3 && "the matrix must have exactly 3 rows");
163 res.row(0) = (mat.row(1) * vec.coeff(2) - mat.row(2) * vec.coeff(1)).conjugate();
164 res.row(1) = (mat.row(2) * vec.coeff(0) - mat.row(0) * vec.coeff(2)).conjugate();
165 res.row(2) = (mat.row(0) * vec.coeff(1) - mat.row(1) * vec.coeff(0)).conjugate();
166 } else {
167 eigen_assert(CrossReturnType::ColsAtCompileTime == 3 && "the matrix must have exactly 3 columns");
168 res.col(0) = (mat.col(1) * vec.coeff(2) - mat.col(2) * vec.coeff(1)).conjugate();
169 res.col(1) = (mat.col(2) * vec.coeff(0) - mat.col(0) * vec.coeff(2)).conjugate();
170 res.col(2) = (mat.col(0) * vec.coeff(1) - mat.col(1) * vec.coeff(0)).conjugate();
171 }
172 return res;
173}
174
175namespace internal {
176
177template <typename Derived, int Size = Derived::SizeAtCompileTime>
178struct unitOrthogonal_selector {
179 typedef typename plain_matrix_type<Derived>::type VectorType;
180 typedef typename traits<Derived>::Scalar Scalar;
181 typedef typename NumTraits<Scalar>::Real RealScalar;
183 EIGEN_DEVICE_FUNC static inline VectorType run(const Derived& src) {
184 VectorType perp = VectorType::Zero(src.size());
185 Index maxi = 0;
186 Index sndi = 0;
187 src.cwiseAbs().maxCoeff(&maxi);
188 if (maxi == 0) sndi = 1;
189 RealScalar invnm = RealScalar(1) / (Vector2() << src.coeff(sndi), src.coeff(maxi)).finished().norm();
190 perp.coeffRef(maxi) = -numext::conj(src.coeff(sndi)) * invnm;
191 perp.coeffRef(sndi) = numext::conj(src.coeff(maxi)) * invnm;
192
193 return perp;
194 }
195};
196
197template <typename Derived>
198struct unitOrthogonal_selector<Derived, 3> {
199 typedef typename plain_matrix_type<Derived>::type VectorType;
200 typedef typename traits<Derived>::Scalar Scalar;
201 typedef typename NumTraits<Scalar>::Real RealScalar;
202 EIGEN_DEVICE_FUNC static inline VectorType run(const Derived& src) {
203 VectorType perp;
204 /* Let us compute the crossed product of *this with a vector
205 * that is not too close to being colinear to *this.
206 */
207
208 /* unless the x and y coords are both close to zero, we can
209 * simply take ( -y, x, 0 ) and normalize it.
210 */
211 if ((!isMuchSmallerThan(src.x(), src.z())) || (!isMuchSmallerThan(src.y(), src.z()))) {
212 RealScalar invnm = RealScalar(1) / src.template head<2>().norm();
213 perp.coeffRef(0) = -numext::conj(src.y()) * invnm;
214 perp.coeffRef(1) = numext::conj(src.x()) * invnm;
215 perp.coeffRef(2) = 0;
216 }
217 /* if both x and y are close to zero, then the vector is close
218 * to the z-axis, so it's far from colinear to the x-axis for instance.
219 * So we take the crossed product with (1,0,0) and normalize it.
220 */
221 else {
222 RealScalar invnm = RealScalar(1) / src.template tail<2>().norm();
223 perp.coeffRef(0) = 0;
224 perp.coeffRef(1) = -numext::conj(src.z()) * invnm;
225 perp.coeffRef(2) = numext::conj(src.y()) * invnm;
226 }
227
228 return perp;
229 }
230};
231
232template <typename Derived>
233struct unitOrthogonal_selector<Derived, 2> {
234 typedef typename plain_matrix_type<Derived>::type VectorType;
235 EIGEN_DEVICE_FUNC static inline VectorType run(const Derived& src) {
236 return VectorType(-numext::conj(src.y()), numext::conj(src.x())).normalized();
237 }
238};
239
240} // end namespace internal
241
251template <typename Derived>
252EIGEN_DEVICE_FUNC typename MatrixBase<Derived>::PlainObject MatrixBase<Derived>::unitOrthogonal() const {
253 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
254 return internal::unitOrthogonal_selector<Derived>::run(derived());
255}
256
257} // end namespace Eigen
258
259#endif // EIGEN_ORTHOMETHODS_H
@ ColsAtCompileTime
Definition DenseBase.h:102
Derived & derived()
Definition EigenBase.h:49
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
constexpr const Scalar & coeff(Index rowId, Index colId) const
Definition PlainObjectBase.h:198
Pseudo expression providing broadcasting and partial reduction operations.
Definition VectorwiseOp.h:192
PlainObject unitOrthogonal(void) const
Definition OrthoMethods.h:252
@ Vertical
Definition Constants.h:266
Namespace containing all symbols from the Eigen library.
Definition Core:137
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:83