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 namespace Eigen { 14 15 /** \class ForceAlignedAccess 16 * \ingroup Core_Module 17 * 18 * \brief Enforce aligned packet loads and stores regardless of what is requested 19 * 20 * \param ExpressionType the type of the object of which we are forcing aligned packet access 21 * 22 * This class is the return type of MatrixBase::forceAlignedAccess() 23 * and most of the time this is the only way it is used. 24 * 25 * \sa MatrixBase::forceAlignedAccess() 26 */ 27 28 namespace internal { 29 template<typename ExpressionType> 30 struct traits<ForceAlignedAccess<ExpressionType> > : public traits<ExpressionType> 31 {}; 32 } 33 34 template<typename ExpressionType> class ForceAlignedAccess 35 : public internal::dense_xpr_base< ForceAlignedAccess<ExpressionType> >::type 36 { 37 public: 38 39 typedef typename internal::dense_xpr_base<ForceAlignedAccess>::type Base; 40 EIGEN_DENSE_PUBLIC_INTERFACE(ForceAlignedAccess) 41 42 inline ForceAlignedAccess(const ExpressionType& matrix) : m_expression(matrix) {} 43 44 inline Index rows() const { return m_expression.rows(); } 45 inline Index cols() const { return m_expression.cols(); } 46 inline Index outerStride() const { return m_expression.outerStride(); } 47 inline Index innerStride() const { return m_expression.innerStride(); } 48 49 inline const CoeffReturnType coeff(Index row, Index col) const 50 { 51 return m_expression.coeff(row, col); 52 } 53 54 inline Scalar& coeffRef(Index row, Index col) 55 { 56 return m_expression.const_cast_derived().coeffRef(row, col); 57 } 58 59 inline const CoeffReturnType coeff(Index index) const 60 { 61 return m_expression.coeff(index); 62 } 63 64 inline Scalar& coeffRef(Index index) 65 { 66 return m_expression.const_cast_derived().coeffRef(index); 67 } 68 69 template<int LoadMode> 70 inline const PacketScalar packet(Index row, Index col) const 71 { 72 return m_expression.template packet<Aligned>(row, col); 73 } 74 75 template<int LoadMode> 76 inline void writePacket(Index row, Index col, const PacketScalar& x) 77 { 78 m_expression.const_cast_derived().template writePacket<Aligned>(row, col, x); 79 } 80 81 template<int LoadMode> 82 inline const PacketScalar packet(Index index) const 83 { 84 return m_expression.template packet<Aligned>(index); 85 } 86 87 template<int LoadMode> 88 inline void writePacket(Index index, const PacketScalar& x) 89 { 90 m_expression.const_cast_derived().template writePacket<Aligned>(index, x); 91 } 92 93 operator const ExpressionType&() const { return m_expression; } 94 95 protected: 96 const ExpressionType& m_expression; 97 98 private: 99 ForceAlignedAccess& operator=(const ForceAlignedAccess&); 100 }; 101 102 /** \returns an expression of *this with forced aligned access 103 * \sa forceAlignedAccessIf(),class ForceAlignedAccess 104 */ 105 template<typename Derived> 106 inline const ForceAlignedAccess<Derived> 107 MatrixBase<Derived>::forceAlignedAccess() const 108 { 109 return ForceAlignedAccess<Derived>(derived()); 110 } 111 112 /** \returns an expression of *this with forced aligned access 113 * \sa forceAlignedAccessIf(), class ForceAlignedAccess 114 */ 115 template<typename Derived> 116 inline ForceAlignedAccess<Derived> 117 MatrixBase<Derived>::forceAlignedAccess() 118 { 119 return ForceAlignedAccess<Derived>(derived()); 120 } 121 122 /** \returns an expression of *this with forced aligned access if \a Enable is true. 123 * \sa forceAlignedAccess(), class ForceAlignedAccess 124 */ 125 template<typename Derived> 126 template<bool Enable> 127 inline typename internal::add_const_on_value_type<typename internal::conditional<Enable,ForceAlignedAccess<Derived>,Derived&>::type>::type 128 MatrixBase<Derived>::forceAlignedAccessIf() const 129 { 130 return derived(); 131 } 132 133 /** \returns an expression of *this with forced aligned access if \a Enable is true. 134 * \sa forceAlignedAccess(), class ForceAlignedAccess 135 */ 136 template<typename Derived> 137 template<bool Enable> 138 inline typename internal::conditional<Enable,ForceAlignedAccess<Derived>,Derived&>::type 139 MatrixBase<Derived>::forceAlignedAccessIf() 140 { 141 return derived(); 142 } 143 144 } // end namespace Eigen 145 146 #endif // EIGEN_FORCEALIGNEDACCESS_H 147