1
// This file is part of Eigen, a lightweight C++ template library
4
// Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
6
// Eigen is free software; you can redistribute it and/or
7
// modify it under the terms of the GNU Lesser General Public
8
// License as published by the Free Software Foundation; either
9
// version 3 of the License, or (at your option) any later version.
11
// Alternatively, you can redistribute it and/or
12
// modify it under the terms of the GNU General Public License as
13
// published by the Free Software Foundation; either version 2 of
14
// the License, or (at your option) any later version.
16
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
17
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
19
// GNU General Public License for more details.
21
// You should have received a copy of the GNU Lesser General Public
22
// License and a copy of the GNU General Public License along with
23
// Eigen. If not, see <http://www.gnu.org/licenses/>.
25
#ifndef EIGEN_FORCEALIGNEDACCESS_H
26
#define EIGEN_FORCEALIGNEDACCESS_H
28
/** \class ForceAlignedAccess
29
* \ingroup Core_Module
31
* \brief Enforce aligned packet loads and stores regardless of what is requested
33
* \param ExpressionType the type of the object of which we are forcing aligned packet access
35
* This class is the return type of MatrixBase::forceAlignedAccess()
36
* and most of the time this is the only way it is used.
38
* \sa MatrixBase::forceAlignedAccess()
42
template<typename ExpressionType>
43
struct traits<ForceAlignedAccess<ExpressionType> > : public traits<ExpressionType>
47
template<typename ExpressionType> class ForceAlignedAccess
48
: public internal::dense_xpr_base< ForceAlignedAccess<ExpressionType> >::type
52
typedef typename internal::dense_xpr_base<ForceAlignedAccess>::type Base;
53
EIGEN_DENSE_PUBLIC_INTERFACE(ForceAlignedAccess)
55
inline ForceAlignedAccess(const ExpressionType& matrix) : m_expression(matrix) {}
57
inline Index rows() const { return m_expression.rows(); }
58
inline Index cols() const { return m_expression.cols(); }
59
inline Index outerStride() const { return m_expression.outerStride(); }
60
inline Index innerStride() const { return m_expression.innerStride(); }
62
inline const CoeffReturnType coeff(Index row, Index col) const
64
return m_expression.coeff(row, col);
67
inline Scalar& coeffRef(Index row, Index col)
69
return m_expression.const_cast_derived().coeffRef(row, col);
72
inline const CoeffReturnType coeff(Index index) const
74
return m_expression.coeff(index);
77
inline Scalar& coeffRef(Index index)
79
return m_expression.const_cast_derived().coeffRef(index);
82
template<int LoadMode>
83
inline const PacketScalar packet(Index row, Index col) const
85
return m_expression.template packet<Aligned>(row, col);
88
template<int LoadMode>
89
inline void writePacket(Index row, Index col, const PacketScalar& x)
91
m_expression.const_cast_derived().template writePacket<Aligned>(row, col, x);
94
template<int LoadMode>
95
inline const PacketScalar packet(Index index) const
97
return m_expression.template packet<Aligned>(index);
100
template<int LoadMode>
101
inline void writePacket(Index index, const PacketScalar& x)
103
m_expression.const_cast_derived().template writePacket<Aligned>(index, x);
106
operator const ExpressionType&() const { return m_expression; }
109
const ExpressionType& m_expression;
112
ForceAlignedAccess& operator=(const ForceAlignedAccess&);
115
/** \returns an expression of *this with forced aligned access
116
* \sa forceAlignedAccessIf(),class ForceAlignedAccess
118
template<typename Derived>
119
inline const ForceAlignedAccess<Derived>
120
MatrixBase<Derived>::forceAlignedAccess() const
122
return ForceAlignedAccess<Derived>(derived());
125
/** \returns an expression of *this with forced aligned access
126
* \sa forceAlignedAccessIf(), class ForceAlignedAccess
128
template<typename Derived>
129
inline ForceAlignedAccess<Derived>
130
MatrixBase<Derived>::forceAlignedAccess()
132
return ForceAlignedAccess<Derived>(derived());
135
/** \returns an expression of *this with forced aligned access if \a Enable is true.
136
* \sa forceAlignedAccess(), class ForceAlignedAccess
138
template<typename Derived>
139
template<bool Enable>
140
inline typename internal::add_const_on_value_type<typename internal::conditional<Enable,ForceAlignedAccess<Derived>,Derived&>::type>::type
141
MatrixBase<Derived>::forceAlignedAccessIf() const
146
/** \returns an expression of *this with forced aligned access if \a Enable is true.
147
* \sa forceAlignedAccess(), class ForceAlignedAccess
149
template<typename Derived>
150
template<bool Enable>
151
inline typename internal::conditional<Enable,ForceAlignedAccess<Derived>,Derived&>::type
152
MatrixBase<Derived>::forceAlignedAccessIf()
157
#endif // EIGEN_FORCEALIGNEDACCESS_H