aboutsummaryrefslogtreecommitdiff
path: root/Eigen/src/Core/ForceAlignedAccess.h
blob: 807c7a29346f47747f4feb70bfd4e81a7d138e53 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2009-2010 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_FORCEALIGNEDACCESS_H
#define EIGEN_FORCEALIGNEDACCESS_H

namespace Eigen {

/** \class ForceAlignedAccess
  * \ingroup Core_Module
  *
  * \brief Enforce aligned packet loads and stores regardless of what is requested
  *
  * \param ExpressionType the type of the object of which we are forcing aligned packet access
  *
  * This class is the return type of MatrixBase::forceAlignedAccess()
  * and most of the time this is the only way it is used.
  *
  * \sa MatrixBase::forceAlignedAccess()
  */

namespace internal {
template<typename ExpressionType>
struct traits<ForceAlignedAccess<ExpressionType> > : public traits<ExpressionType>
{};
}

template<typename ExpressionType> class ForceAlignedAccess
  : public internal::dense_xpr_base< ForceAlignedAccess<ExpressionType> >::type
{
  public:

    typedef typename internal::dense_xpr_base<ForceAlignedAccess>::type Base;
    EIGEN_DENSE_PUBLIC_INTERFACE(ForceAlignedAccess)

    inline ForceAlignedAccess(const ExpressionType& matrix) : m_expression(matrix) {}

    inline Index rows() const { return m_expression.rows(); }
    inline Index cols() const { return m_expression.cols(); }
    inline Index outerStride() const { return m_expression.outerStride(); }
    inline Index innerStride() const { return m_expression.innerStride(); }

    inline const CoeffReturnType coeff(Index row, Index col) const
    {
      return m_expression.coeff(row, col);
    }

    inline Scalar& coeffRef(Index row, Index col)
    {
      return m_expression.const_cast_derived().coeffRef(row, col);
    }

    inline const CoeffReturnType coeff(Index index) const
    {
      return m_expression.coeff(index);
    }

    inline Scalar& coeffRef(Index index)
    {
      return m_expression.const_cast_derived().coeffRef(index);
    }

    template<int LoadMode>
    inline const PacketScalar packet(Index row, Index col) const
    {
      return m_expression.template packet<Aligned>(row, col);
    }

    template<int LoadMode>
    inline void writePacket(Index row, Index col, const PacketScalar& x)
    {
      m_expression.const_cast_derived().template writePacket<Aligned>(row, col, x);
    }

    template<int LoadMode>
    inline const PacketScalar packet(Index index) const
    {
      return m_expression.template packet<Aligned>(index);
    }

    template<int LoadMode>
    inline void writePacket(Index index, const PacketScalar& x)
    {
      m_expression.const_cast_derived().template writePacket<Aligned>(index, x);
    }

    operator const ExpressionType&() const { return m_expression; }

  protected:
    const ExpressionType& m_expression;

  private:
    ForceAlignedAccess& operator=(const ForceAlignedAccess&);
};

/** \returns an expression of *this with forced aligned access
  * \sa forceAlignedAccessIf(),class ForceAlignedAccess
  */
template<typename Derived>
inline const ForceAlignedAccess<Derived>
MatrixBase<Derived>::forceAlignedAccess() const
{
  return ForceAlignedAccess<Derived>(derived());
}

/** \returns an expression of *this with forced aligned access
  * \sa forceAlignedAccessIf(), class ForceAlignedAccess
  */
template<typename Derived>
inline ForceAlignedAccess<Derived>
MatrixBase<Derived>::forceAlignedAccess()
{
  return ForceAlignedAccess<Derived>(derived());
}

/** \returns an expression of *this with forced aligned access if \a Enable is true.
  * \sa forceAlignedAccess(), class ForceAlignedAccess
  */
template<typename Derived>
template<bool Enable>
inline typename internal::add_const_on_value_type<typename internal::conditional<Enable,ForceAlignedAccess<Derived>,Derived&>::type>::type
MatrixBase<Derived>::forceAlignedAccessIf() const
{
  return derived();
}

/** \returns an expression of *this with forced aligned access if \a Enable is true.
  * \sa forceAlignedAccess(), class ForceAlignedAccess
  */
template<typename Derived>
template<bool Enable>
inline typename internal::conditional<Enable,ForceAlignedAccess<Derived>,Derived&>::type
MatrixBase<Derived>::forceAlignedAccessIf()
{
  return derived();
}

} // end namespace Eigen

#endif // EIGEN_FORCEALIGNEDACCESS_H