$darkmode
Eigen  5.0.1-dev
ForceAlignedAccess.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_FORCEALIGNEDACCESS_H
11 #define EIGEN_FORCEALIGNEDACCESS_H
12 
13 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
31 namespace internal {
32 template <typename ExpressionType>
33 struct traits<ForceAlignedAccess<ExpressionType>> : public traits<ExpressionType> {};
34 } // namespace internal
35 
36 template <typename ExpressionType>
37 class ForceAlignedAccess : public internal::dense_xpr_base<ForceAlignedAccess<ExpressionType>>::type {
38  public:
39  typedef typename internal::dense_xpr_base<ForceAlignedAccess>::type Base;
40  EIGEN_DENSE_PUBLIC_INTERFACE(ForceAlignedAccess)
41 
42  EIGEN_DEVICE_FUNC explicit inline ForceAlignedAccess(const ExpressionType& matrix) : m_expression(matrix) {}
43 
44  EIGEN_DEVICE_FUNC constexpr Index rows() const noexcept { return m_expression.rows(); }
45  EIGEN_DEVICE_FUNC constexpr Index cols() const noexcept { return m_expression.cols(); }
46  EIGEN_DEVICE_FUNC constexpr Index outerStride() const noexcept { return m_expression.outerStride(); }
47  EIGEN_DEVICE_FUNC constexpr Index innerStride() const noexcept { return m_expression.innerStride(); }
48 
49  EIGEN_DEVICE_FUNC inline const CoeffReturnType coeff(Index row, Index col) const {
50  return m_expression.coeff(row, col);
51  }
52 
53  EIGEN_DEVICE_FUNC inline Scalar& coeffRef(Index row, Index col) {
54  return m_expression.const_cast_derived().coeffRef(row, col);
55  }
56 
57  EIGEN_DEVICE_FUNC inline const CoeffReturnType coeff(Index index) const { return m_expression.coeff(index); }
58 
59  EIGEN_DEVICE_FUNC inline Scalar& coeffRef(Index index) { return m_expression.const_cast_derived().coeffRef(index); }
60 
61  template <int LoadMode>
62  inline const PacketScalar packet(Index row, Index col) const {
63  return m_expression.template packet<Aligned>(row, col);
64  }
65 
66  template <int LoadMode>
67  inline void writePacket(Index row, Index col, const PacketScalar& x) {
68  m_expression.const_cast_derived().template writePacket<Aligned>(row, col, x);
69  }
70 
71  template <int LoadMode>
72  inline const PacketScalar packet(Index index) const {
73  return m_expression.template packet<Aligned>(index);
74  }
75 
76  template <int LoadMode>
77  inline void writePacket(Index index, const PacketScalar& x) {
78  m_expression.const_cast_derived().template writePacket<Aligned>(index, x);
79  }
80 
81  EIGEN_DEVICE_FUNC operator const ExpressionType&() const { return m_expression; }
82 
83  protected:
84  const ExpressionType& m_expression;
85 
86  private:
87  ForceAlignedAccess& operator=(const ForceAlignedAccess&);
88 };
89 
93 template <typename Derived>
95  return ForceAlignedAccess<Derived>(derived());
96 }
97 
101 template <typename Derived>
102 inline ForceAlignedAccess<Derived> MatrixBase<Derived>::forceAlignedAccess() {
103  return ForceAlignedAccess<Derived>(derived());
104 }
105 
109 template <typename Derived>
110 template <bool Enable>
111 inline add_const_on_value_type_t<std::conditional_t<Enable, ForceAlignedAccess<Derived>, Derived&>>
113  return derived(); // FIXME This should not work but apparently is never used
114 }
115 
119 template <typename Derived>
120 template <bool Enable>
121 inline std::conditional_t<Enable, ForceAlignedAccess<Derived>, Derived&> MatrixBase<Derived>::forceAlignedAccessIf() {
122  return derived(); // FIXME This should not work but apparently is never used
123 }
124 
125 } // end namespace Eigen
126 
127 #endif // EIGEN_FORCEALIGNEDACCESS_H
Enforce aligned packet loads and stores regardless of what is requested.
Definition: ForceAlignedAccess.h:37
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
const Derived & forceAlignedAccess() const
Definition: MatrixBase.h:299
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52