aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Core/arch/Default/ConjHelper.h
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2017-06-15 10:16:30 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2017-06-15 10:16:30 +0200
commitb240080e6443e4fb203ac1cba0ad5bd4fdad56b3 (patch)
treeb4b68258306115e23bad6bd4a73803c326bd86bd /Eigen/src/Core/arch/Default/ConjHelper.h
parentb8e805497e446e7159f231238b4a8fd22fe70749 (diff)
bug #1436: fix compilation of Jacobi rotations with ARM NEON, some specializations of internal::conj_helper were missing.
Diffstat (limited to 'Eigen/src/Core/arch/Default/ConjHelper.h')
-rw-r--r--Eigen/src/Core/arch/Default/ConjHelper.h29
1 files changed, 29 insertions, 0 deletions
diff --git a/Eigen/src/Core/arch/Default/ConjHelper.h b/Eigen/src/Core/arch/Default/ConjHelper.h
new file mode 100644
index 000000000..4cfe34e05
--- /dev/null
+++ b/Eigen/src/Core/arch/Default/ConjHelper.h
@@ -0,0 +1,29 @@
+
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2017 Gael Guennebaud <gael.guennebaud@inria.fr>
+//
+// This Source Code Form is subject to the terms of the Mozilla
+// Public License v. 2.0. If a copy of the MPL was not distributed
+// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
+
+#ifndef EIGEN_ARCH_CONJ_HELPER_H
+#define EIGEN_ARCH_CONJ_HELPER_H
+
+#define EIGEN_MAKE_CONJ_HELPER_CPLX_REAL(PACKET_CPLX, PACKET_REAL) \
+ template<> struct conj_helper<PACKET_REAL, PACKET_CPLX, false,false> { \
+ EIGEN_STRONG_INLINE PACKET_CPLX pmadd(const PACKET_REAL& x, const PACKET_CPLX& y, const PACKET_CPLX& c) const \
+ { return padd(c, pmul(x,y)); } \
+ EIGEN_STRONG_INLINE PACKET_CPLX pmul(const PACKET_REAL& x, const PACKET_CPLX& y) const \
+ { return PACKET_CPLX(Eigen::internal::pmul<PACKET_REAL>(x, y.v)); } \
+ }; \
+ \
+ template<> struct conj_helper<PACKET_CPLX, PACKET_REAL, false,false> { \
+ EIGEN_STRONG_INLINE PACKET_CPLX pmadd(const PACKET_CPLX& x, const PACKET_REAL& y, const PACKET_CPLX& c) const \
+ { return padd(c, pmul(x,y)); } \
+ EIGEN_STRONG_INLINE PACKET_CPLX pmul(const PACKET_CPLX& x, const PACKET_REAL& y) const \
+ { return PACKET_CPLX(Eigen::internal::pmul<PACKET_REAL>(x.v, y)); } \
+ };
+
+#endif // EIGEN_ARCH_CONJ_HELPER_H