11#ifndef EIGEN_ORTHOMETHODS_H
12#define EIGEN_ORTHOMETHODS_H
15#include "./InternalHeaderCheck.h"
22template <
typename Derived,
typename OtherDerived,
int Size>
24 typedef typename ScalarBinaryOpTraits<typename internal::traits<Derived>::Scalar,
25 typename internal::traits<OtherDerived>::Scalar>::ReturnType Scalar;
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)
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)));
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;
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));
88template <
typename Derived>
89template <
typename OtherDerived>
90EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
typename internal::cross_impl<Derived, OtherDerived>::return_type
92 return internal::cross_impl<Derived, OtherDerived>::run(*
this, other);
97template <
int Arch,
typename VectorLhs,
typename VectorRhs,
typename Scalar =
typename VectorLhs::Scalar,
99 bool((
int(evaluator<VectorLhs>::Flags) &
int(evaluator<VectorRhs>::Flags)) &
PacketAccessBit)>
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);
121template <
typename Derived>
122template <
typename OtherDerived>
124 const MatrixBase<OtherDerived>& other)
const {
125 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived, 4)
126 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived, 4)
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());
133 return internal::cross3_impl<Architecture::Target, internal::remove_all_t<DerivedNested>,
134 internal::remove_all_t<OtherDerivedNested>>::run(lhs, rhs);
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)
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)
155 typename internal::nested_eval<ExpressionType, 2>::type mat(_expression());
156 typename internal::nested_eval<OtherDerived, 2>::type vec(other.derived());
158 CrossReturnType res(_expression().rows(), _expression().cols());
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();
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();
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());
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;
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) {
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;
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;
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();
249template <
typename Derived>
251 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
252 return internal::unitOrthogonal_selector<Derived>::run(derived());
@ 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