$darkmode
Eigen  5.0.1-dev
CoreEvaluators.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2011 Benoit Jacob <jacob.benoit.1@gmail.com>
5 // Copyright (C) 2011-2014 Gael Guennebaud <gael.guennebaud@inria.fr>
6 // Copyright (C) 2011-2012 Jitse Niesen <jitse@maths.leeds.ac.uk>
7 //
8 // This Source Code Form is subject to the terms of the Mozilla
9 // Public License v. 2.0. If a copy of the MPL was not distributed
10 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
11 
12 #ifndef EIGEN_COREEVALUATORS_H
13 #define EIGEN_COREEVALUATORS_H
14 
15 // IWYU pragma: private
16 #include "./InternalHeaderCheck.h"
17 
18 namespace Eigen {
19 
20 namespace internal {
21 
22 // This class returns the evaluator kind from the expression storage kind.
23 // Default assumes index based accessors
24 template <typename StorageKind>
25 struct storage_kind_to_evaluator_kind {
26  typedef IndexBased Kind;
27 };
28 
29 // This class returns the evaluator shape from the expression storage kind.
30 // It can be Dense, Sparse, Triangular, Diagonal, SelfAdjoint, Band, etc.
31 template <typename StorageKind>
32 struct storage_kind_to_shape;
33 
34 template <>
35 struct storage_kind_to_shape<Dense> {
36  typedef DenseShape Shape;
37 };
38 template <>
39 struct storage_kind_to_shape<SolverStorage> {
40  typedef SolverShape Shape;
41 };
42 template <>
43 struct storage_kind_to_shape<PermutationStorage> {
44  typedef PermutationShape Shape;
45 };
46 template <>
47 struct storage_kind_to_shape<TranspositionsStorage> {
48  typedef TranspositionsShape Shape;
49 };
50 
51 // Evaluators have to be specialized with respect to various criteria such as:
52 // - storage/structure/shape
53 // - scalar type
54 // - etc.
55 // Therefore, we need specialization of evaluator providing additional template arguments for each kind of evaluators.
56 // We currently distinguish the following kind of evaluators:
57 // - unary_evaluator for expressions taking only one arguments (CwiseUnaryOp, CwiseUnaryView, Transpose,
58 // MatrixWrapper, ArrayWrapper, Reverse, Replicate)
59 // - binary_evaluator for expression taking two arguments (CwiseBinaryOp)
60 // - ternary_evaluator for expression taking three arguments (CwiseTernaryOp)
61 // - product_evaluator for linear algebra products (Product); special case of binary_evaluator because it requires
62 // additional tags for dispatching.
63 // - mapbase_evaluator for Map, Block, Ref
64 // - block_evaluator for Block (special dispatching to a mapbase_evaluator or unary_evaluator)
65 
66 template <typename T, typename Arg1Kind = typename evaluator_traits<typename T::Arg1>::Kind,
67  typename Arg2Kind = typename evaluator_traits<typename T::Arg2>::Kind,
68  typename Arg3Kind = typename evaluator_traits<typename T::Arg3>::Kind,
69  typename Arg1Scalar = typename traits<typename T::Arg1>::Scalar,
70  typename Arg2Scalar = typename traits<typename T::Arg2>::Scalar,
71  typename Arg3Scalar = typename traits<typename T::Arg3>::Scalar>
72 struct ternary_evaluator;
73 
74 template <typename T, typename LhsKind = typename evaluator_traits<typename T::Lhs>::Kind,
75  typename RhsKind = typename evaluator_traits<typename T::Rhs>::Kind,
76  typename LhsScalar = typename traits<typename T::Lhs>::Scalar,
77  typename RhsScalar = typename traits<typename T::Rhs>::Scalar>
78 struct binary_evaluator;
79 
80 template <typename T, typename Kind = typename evaluator_traits<typename T::NestedExpression>::Kind,
81  typename Scalar = typename T::Scalar>
82 struct unary_evaluator;
83 
84 // evaluator_traits<T> contains traits for evaluator<T>
85 
86 template <typename T>
87 struct evaluator_traits_base {
88  // by default, get evaluator kind and shape from storage
89  typedef typename storage_kind_to_evaluator_kind<typename traits<T>::StorageKind>::Kind Kind;
90  typedef typename storage_kind_to_shape<typename traits<T>::StorageKind>::Shape Shape;
91 };
92 
93 // Default evaluator traits
94 template <typename T>
95 struct evaluator_traits : public evaluator_traits_base<T> {};
96 
97 template <typename T, typename Shape = typename evaluator_traits<T>::Shape>
98 struct evaluator_assume_aliasing {
99  static const bool value = false;
100 };
101 
102 // By default, we assume a unary expression:
103 template <typename T>
104 struct evaluator : public unary_evaluator<T> {
105  typedef unary_evaluator<T> Base;
106  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit evaluator(const T& xpr) : Base(xpr) {}
107 };
108 
109 // TODO: Think about const-correctness
110 template <typename T>
111 struct evaluator<const T> : evaluator<T> {
112  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit evaluator(const T& xpr) : evaluator<T>(xpr) {}
113 };
114 
115 // ---------- base class for all evaluators ----------
116 
117 template <typename ExpressionType>
118 struct evaluator_base {
119  // TODO that's not very nice to have to propagate all these traits. They are currently only needed to handle
120  // outer,inner indices.
121  typedef traits<ExpressionType> ExpressionTraits;
122 
123  enum { Alignment = 0 };
124  // noncopyable:
125  // Don't make this class inherit noncopyable as this kills EBO (Empty Base Optimization)
126  // and make complex evaluator much larger than then should do.
127  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr evaluator_base() = default;
128 
129  private:
130  EIGEN_DEVICE_FUNC evaluator_base(const evaluator_base&);
131  EIGEN_DEVICE_FUNC const evaluator_base& operator=(const evaluator_base&);
132 };
133 
134 // -------------------- Matrix and Array --------------------
135 //
136 // evaluator<PlainObjectBase> is a common base class for the
137 // Matrix and Array evaluators.
138 // Here we directly specialize evaluator. This is not really a unary expression, and it is, by definition, dense,
139 // so no need for more sophisticated dispatching.
140 
141 // this helper permits to completely eliminate m_outerStride if it is known at compiletime.
142 template <typename Scalar, int OuterStride>
143 class plainobjectbase_evaluator_data {
144  public:
145  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr plainobjectbase_evaluator_data(const Scalar* ptr, Index outerStride)
146  : data(ptr) {
147 #ifndef EIGEN_INTERNAL_DEBUGGING
148  EIGEN_UNUSED_VARIABLE(outerStride);
149 #endif
150  eigen_internal_assert(outerStride == OuterStride);
151  }
152  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr Index outerStride() const noexcept { return OuterStride; }
153  const Scalar* data;
154 };
155 
156 template <typename Scalar>
157 class plainobjectbase_evaluator_data<Scalar, Dynamic> {
158  public:
159  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr plainobjectbase_evaluator_data(const Scalar* ptr, Index outerStride)
160  : data(ptr), m_outerStride(outerStride) {}
161  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr Index outerStride() const { return m_outerStride; }
162  const Scalar* data;
163 
164  protected:
165  Index m_outerStride;
166 };
167 
168 template <typename Derived>
169 struct evaluator<PlainObjectBase<Derived>> : evaluator_base<Derived> {
170  typedef PlainObjectBase<Derived> PlainObjectType;
171  typedef typename PlainObjectType::Scalar Scalar;
172  typedef typename PlainObjectType::CoeffReturnType CoeffReturnType;
173 
174  enum {
175  IsRowMajor = PlainObjectType::IsRowMajor,
176  IsVectorAtCompileTime = PlainObjectType::IsVectorAtCompileTime,
177  RowsAtCompileTime = PlainObjectType::RowsAtCompileTime,
178  ColsAtCompileTime = PlainObjectType::ColsAtCompileTime,
179 
180  CoeffReadCost = NumTraits<Scalar>::ReadCost,
181  Flags = traits<Derived>::EvaluatorFlags,
182  Alignment = traits<Derived>::Alignment
183  };
184  enum {
185  // We do not need to know the outer stride for vectors
186  OuterStrideAtCompileTime = IsVectorAtCompileTime ? 0
187  : int(IsRowMajor) ? ColsAtCompileTime
188  : RowsAtCompileTime
189  };
190 
191  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr evaluator() : m_d(0, OuterStrideAtCompileTime) {
192  EIGEN_INTERNAL_CHECK_COST_VALUE(CoeffReadCost);
193  }
194 
195  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr explicit evaluator(const PlainObjectType& m)
196  : m_d(m.data(), IsVectorAtCompileTime ? 0 : m.outerStride()) {
197  EIGEN_INTERNAL_CHECK_COST_VALUE(CoeffReadCost);
198  }
199 
200  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr CoeffReturnType coeff(Index row, Index col) const {
201  return coeff(getIndex(row, col));
202  }
203 
204  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr CoeffReturnType coeff(Index index) const { return m_d.data[index]; }
205 
206  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr Scalar& coeffRef(Index row, Index col) {
207  return coeffRef(getIndex(row, col));
208  }
209 
210  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr Scalar& coeffRef(Index index) {
211  return const_cast<Scalar*>(m_d.data)[index];
212  }
213 
214  template <int LoadMode, typename PacketType>
215  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index row, Index col) const {
216  return packet<LoadMode, PacketType>(getIndex(row, col));
217  }
218 
219  template <int LoadMode, typename PacketType>
220  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index index) const {
221  return ploadt<PacketType, LoadMode>(m_d.data + index);
222  }
223 
224  template <int StoreMode, typename PacketType>
225  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacket(Index row, Index col, const PacketType& x) {
226  writePacket<StoreMode, PacketType>(getIndex(row, col), x);
227  }
228 
229  template <int StoreMode, typename PacketType>
230  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacket(Index index, const PacketType& x) {
231  pstoret<Scalar, PacketType, StoreMode>(const_cast<Scalar*>(m_d.data) + index, x);
232  }
233 
234  template <int LoadMode, typename PacketType>
235  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index row, Index col, Index begin, Index count) const {
236  return packetSegment<LoadMode, PacketType>(getIndex(row, col), begin, count);
237  }
238 
239  template <int LoadMode, typename PacketType>
240  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index index, Index begin, Index count) const {
241  return ploadtSegment<PacketType, LoadMode>(m_d.data + index, begin, count);
242  }
243 
244  template <int StoreMode, typename PacketType>
245  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacketSegment(Index row, Index col, const PacketType& x, Index begin,
246  Index count) {
247  writePacketSegment<StoreMode, PacketType>(getIndex(row, col), x, begin, count);
248  }
249 
250  template <int StoreMode, typename PacketType>
251  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacketSegment(Index index, const PacketType& x, Index begin,
252  Index count) {
253  pstoretSegment<Scalar, PacketType, StoreMode>(const_cast<Scalar*>(m_d.data) + index, x, begin, count);
254  }
255 
256  protected:
257  plainobjectbase_evaluator_data<Scalar, OuterStrideAtCompileTime> m_d;
258 
259  private:
260  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index constexpr getIndex(Index row, Index col) const {
261  return IsRowMajor ? row * m_d.outerStride() + col : row + col * m_d.outerStride();
262  }
263 };
264 
265 template <typename Scalar, int Rows, int Cols, int Options, int MaxRows, int MaxCols>
266 struct evaluator<Matrix<Scalar, Rows, Cols, Options, MaxRows, MaxCols>>
267  : evaluator<PlainObjectBase<Matrix<Scalar, Rows, Cols, Options, MaxRows, MaxCols>>> {
268  typedef Matrix<Scalar, Rows, Cols, Options, MaxRows, MaxCols> XprType;
269 
270  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr evaluator() = default;
271 
272  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr explicit evaluator(const XprType& m)
273  : evaluator<PlainObjectBase<XprType>>(m) {}
274 };
275 
276 template <typename Scalar, int Rows, int Cols, int Options, int MaxRows, int MaxCols>
277 struct evaluator<Array<Scalar, Rows, Cols, Options, MaxRows, MaxCols>>
278  : evaluator<PlainObjectBase<Array<Scalar, Rows, Cols, Options, MaxRows, MaxCols>>> {
279  typedef Array<Scalar, Rows, Cols, Options, MaxRows, MaxCols> XprType;
280 
281  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr evaluator() = default;
282 
283  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr explicit evaluator(const XprType& m)
284  : evaluator<PlainObjectBase<XprType>>(m) {}
285 };
286 
287 // -------------------- Transpose --------------------
288 
289 template <typename ArgType>
290 struct unary_evaluator<Transpose<ArgType>, IndexBased> : evaluator_base<Transpose<ArgType>> {
291  typedef Transpose<ArgType> XprType;
292 
293  enum {
294  CoeffReadCost = evaluator<ArgType>::CoeffReadCost,
295  Flags = evaluator<ArgType>::Flags ^ RowMajorBit,
296  Alignment = evaluator<ArgType>::Alignment
297  };
298 
299  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit unary_evaluator(const XprType& t) : m_argImpl(t.nestedExpression()) {}
300 
301  typedef typename XprType::Scalar Scalar;
302  typedef typename XprType::CoeffReturnType CoeffReturnType;
303 
304  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index row, Index col) const {
305  return m_argImpl.coeff(col, row);
306  }
307 
308  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index index) const { return m_argImpl.coeff(index); }
309 
310  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index row, Index col) { return m_argImpl.coeffRef(col, row); }
311 
312  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE typename XprType::Scalar& coeffRef(Index index) {
313  return m_argImpl.coeffRef(index);
314  }
315 
316  template <int LoadMode, typename PacketType>
317  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index row, Index col) const {
318  return m_argImpl.template packet<LoadMode, PacketType>(col, row);
319  }
320 
321  template <int LoadMode, typename PacketType>
322  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index index) const {
323  return m_argImpl.template packet<LoadMode, PacketType>(index);
324  }
325 
326  template <int StoreMode, typename PacketType>
327  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacket(Index row, Index col, const PacketType& x) {
328  m_argImpl.template writePacket<StoreMode, PacketType>(col, row, x);
329  }
330 
331  template <int StoreMode, typename PacketType>
332  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacket(Index index, const PacketType& x) {
333  m_argImpl.template writePacket<StoreMode, PacketType>(index, x);
334  }
335 
336  template <int LoadMode, typename PacketType>
337  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index row, Index col, Index begin, Index count) const {
338  return m_argImpl.template packetSegment<LoadMode, PacketType>(col, row, begin, count);
339  }
340 
341  template <int LoadMode, typename PacketType>
342  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index index, Index begin, Index count) const {
343  return m_argImpl.template packetSegment<LoadMode, PacketType>(index, begin, count);
344  }
345 
346  template <int StoreMode, typename PacketType>
347  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacketSegment(Index row, Index col, const PacketType& x, Index begin,
348  Index count) {
349  m_argImpl.template writePacketSegment<StoreMode, PacketType>(col, row, x, begin, count);
350  }
351 
352  template <int StoreMode, typename PacketType>
353  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacketSegment(Index index, const PacketType& x, Index begin,
354  Index count) {
355  m_argImpl.template writePacketSegment<StoreMode, PacketType>(index, x, begin, count);
356  }
357 
358  protected:
359  evaluator<ArgType> m_argImpl;
360 };
361 
362 // -------------------- CwiseNullaryOp --------------------
363 // Like Matrix and Array, this is not really a unary expression, so we directly specialize evaluator.
364 // Likewise, there is not need to more sophisticated dispatching here.
365 
366 template <typename Scalar, typename NullaryOp, bool has_nullary = has_nullary_operator<NullaryOp>::value,
367  bool has_unary = has_unary_operator<NullaryOp>::value,
368  bool has_binary = has_binary_operator<NullaryOp>::value>
369 struct nullary_wrapper {
370  template <typename IndexType>
371  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar operator()(const NullaryOp& op, IndexType i, IndexType j) const {
372  return op(i, j);
373  }
374  template <typename IndexType>
375  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar operator()(const NullaryOp& op, IndexType i) const {
376  return op(i);
377  }
378 
379  template <typename T, typename IndexType>
380  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE T packetOp(const NullaryOp& op, IndexType i, IndexType j) const {
381  return op.template packetOp<T>(i, j);
382  }
383  template <typename T, typename IndexType>
384  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE T packetOp(const NullaryOp& op, IndexType i) const {
385  return op.template packetOp<T>(i);
386  }
387 };
388 
389 template <typename Scalar, typename NullaryOp>
390 struct nullary_wrapper<Scalar, NullaryOp, true, false, false> {
391  template <typename IndexType>
392  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar operator()(const NullaryOp& op, IndexType = 0, IndexType = 0) const {
393  return op();
394  }
395  template <typename T, typename IndexType>
396  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE T packetOp(const NullaryOp& op, IndexType = 0, IndexType = 0) const {
397  return op.template packetOp<T>();
398  }
399 };
400 
401 template <typename Scalar, typename NullaryOp>
402 struct nullary_wrapper<Scalar, NullaryOp, false, false, true> {
403  template <typename IndexType>
404  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar operator()(const NullaryOp& op, IndexType i, IndexType j = 0) const {
405  return op(i, j);
406  }
407  template <typename T, typename IndexType>
408  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE T packetOp(const NullaryOp& op, IndexType i, IndexType j = 0) const {
409  return op.template packetOp<T>(i, j);
410  }
411 };
412 
413 // We need the following specialization for vector-only functors assigned to a runtime vector,
414 // for instance, using linspace and assigning a RowVectorXd to a MatrixXd or even a row of a MatrixXd.
415 // In this case, i==0 and j is used for the actual iteration.
416 template <typename Scalar, typename NullaryOp>
417 struct nullary_wrapper<Scalar, NullaryOp, false, true, false> {
418  template <typename IndexType>
419  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar operator()(const NullaryOp& op, IndexType i, IndexType j) const {
420  eigen_assert(i == 0 || j == 0);
421  return op(i + j);
422  }
423  template <typename T, typename IndexType>
424  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE T packetOp(const NullaryOp& op, IndexType i, IndexType j) const {
425  eigen_assert(i == 0 || j == 0);
426  return op.template packetOp<T>(i + j);
427  }
428 
429  template <typename IndexType>
430  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar operator()(const NullaryOp& op, IndexType i) const {
431  return op(i);
432  }
433  template <typename T, typename IndexType>
434  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE T packetOp(const NullaryOp& op, IndexType i) const {
435  return op.template packetOp<T>(i);
436  }
437 };
438 
439 template <typename Scalar, typename NullaryOp>
440 struct nullary_wrapper<Scalar, NullaryOp, false, false, false> {};
441 
442 #if 0 && EIGEN_COMP_MSVC > 0
443 // Disable this ugly workaround. This is now handled in traits<Ref>::match,
444 // but this piece of code might still become handly if some other weird compilation
445 // errors pop up again.
446 
447 // MSVC exhibits a weird compilation error when
448 // compiling:
449 // Eigen::MatrixXf A = MatrixXf::Random(3,3);
450 // Ref<const MatrixXf> R = 2.f*A;
451 // and that has_*ary_operator<scalar_constant_op<float>> have not been instantiated yet.
452 // The "problem" is that evaluator<2.f*A> is instantiated by traits<Ref>::match<2.f*A>
453 // and at that time has_*ary_operator<T> returns true regardless of T.
454 // Then nullary_wrapper is badly instantiated as nullary_wrapper<.,.,true,true,true>.
455 // The trick is thus to defer the proper instantiation of nullary_wrapper when coeff(),
456 // and packet() are really instantiated as implemented below:
457 
458 // This is a simple wrapper around Index to enforce the re-instantiation of
459 // has_*ary_operator when needed.
460 template<typename T> struct nullary_wrapper_workaround_msvc {
461  nullary_wrapper_workaround_msvc(const T&);
462  operator T()const;
463 };
464 
465 template<typename Scalar,typename NullaryOp>
466 struct nullary_wrapper<Scalar,NullaryOp,true,true,true>
467 {
468  template <typename IndexType>
469  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar operator()(const NullaryOp& op, IndexType i, IndexType j) const {
470  return nullary_wrapper<Scalar,NullaryOp,
471  has_nullary_operator<NullaryOp,nullary_wrapper_workaround_msvc<IndexType> >::value,
472  has_unary_operator<NullaryOp,nullary_wrapper_workaround_msvc<IndexType> >::value,
473  has_binary_operator<NullaryOp,nullary_wrapper_workaround_msvc<IndexType> >::value>().operator()(op,i,j);
474  }
475  template <typename IndexType>
476  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar operator()(const NullaryOp& op, IndexType i) const {
477  return nullary_wrapper<Scalar,NullaryOp,
478  has_nullary_operator<NullaryOp,nullary_wrapper_workaround_msvc<IndexType> >::value,
479  has_unary_operator<NullaryOp,nullary_wrapper_workaround_msvc<IndexType> >::value,
480  has_binary_operator<NullaryOp,nullary_wrapper_workaround_msvc<IndexType> >::value>().operator()(op,i);
481  }
482 
483  template <typename T, typename IndexType>
484  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE T packetOp(const NullaryOp& op, IndexType i, IndexType j) const {
485  return nullary_wrapper<Scalar,NullaryOp,
486  has_nullary_operator<NullaryOp,nullary_wrapper_workaround_msvc<IndexType> >::value,
487  has_unary_operator<NullaryOp,nullary_wrapper_workaround_msvc<IndexType> >::value,
488  has_binary_operator<NullaryOp,nullary_wrapper_workaround_msvc<IndexType> >::value>().template packetOp<T>(op,i,j);
489  }
490  template <typename T, typename IndexType>
491  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE T packetOp(const NullaryOp& op, IndexType i) const {
492  return nullary_wrapper<Scalar,NullaryOp,
493  has_nullary_operator<NullaryOp,nullary_wrapper_workaround_msvc<IndexType> >::value,
494  has_unary_operator<NullaryOp,nullary_wrapper_workaround_msvc<IndexType> >::value,
495  has_binary_operator<NullaryOp,nullary_wrapper_workaround_msvc<IndexType> >::value>().template packetOp<T>(op,i);
496  }
497 };
498 #endif // MSVC workaround
499 
500 template <typename NullaryOp, typename PlainObjectType>
501 struct evaluator<CwiseNullaryOp<NullaryOp, PlainObjectType>>
502  : evaluator_base<CwiseNullaryOp<NullaryOp, PlainObjectType>> {
503  typedef CwiseNullaryOp<NullaryOp, PlainObjectType> XprType;
504  typedef remove_all_t<PlainObjectType> PlainObjectTypeCleaned;
505 
506  enum {
507  CoeffReadCost = functor_traits<NullaryOp>::Cost,
508 
509  Flags = (evaluator<PlainObjectTypeCleaned>::Flags &
510  (HereditaryBits | (functor_has_linear_access<NullaryOp>::ret ? LinearAccessBit : 0) |
511  (functor_traits<NullaryOp>::PacketAccess ? PacketAccessBit : 0))) |
512  (functor_traits<NullaryOp>::IsRepeatable ? 0 : EvalBeforeNestingBit),
513  Alignment = AlignedMax
514  };
515 
516  EIGEN_DEVICE_FUNC explicit evaluator(const XprType& n) : m_functor(n.functor()), m_wrapper() {
517  EIGEN_INTERNAL_CHECK_COST_VALUE(CoeffReadCost);
518  }
519 
520  typedef typename XprType::CoeffReturnType CoeffReturnType;
521 
522  template <typename IndexType>
523  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(IndexType row, IndexType col) const {
524  return m_wrapper(m_functor, row, col);
525  }
526 
527  template <typename IndexType>
528  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(IndexType index) const {
529  return m_wrapper(m_functor, index);
530  }
531 
532  template <int LoadMode, typename PacketType, typename IndexType>
533  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(IndexType row, IndexType col) const {
534  return m_wrapper.template packetOp<PacketType>(m_functor, row, col);
535  }
536 
537  template <int LoadMode, typename PacketType, typename IndexType>
538  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(IndexType index) const {
539  return m_wrapper.template packetOp<PacketType>(m_functor, index);
540  }
541 
542  template <int LoadMode, typename PacketType, typename IndexType>
543  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(IndexType row, IndexType col, Index /*begin*/,
544  Index /*count*/) const {
545  return packet<LoadMode, PacketType, IndexType>(row, col);
546  }
547 
548  template <int LoadMode, typename PacketType, typename IndexType>
549  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(IndexType index, Index /*begin*/,
550  Index /*count*/) const {
551  return packet<LoadMode, PacketType, IndexType>(index);
552  }
553 
554  protected:
555  const NullaryOp m_functor;
556  const nullary_wrapper<CoeffReturnType, NullaryOp> m_wrapper;
557 };
558 
559 // -------------------- CwiseUnaryOp --------------------
560 
561 template <typename UnaryOp, typename ArgType>
562 struct unary_evaluator<CwiseUnaryOp<UnaryOp, ArgType>, IndexBased> : evaluator_base<CwiseUnaryOp<UnaryOp, ArgType>> {
563  typedef CwiseUnaryOp<UnaryOp, ArgType> XprType;
564 
565  enum {
566  CoeffReadCost = int(evaluator<ArgType>::CoeffReadCost) + int(functor_traits<UnaryOp>::Cost),
567 
568  Flags = evaluator<ArgType>::Flags &
569  (HereditaryBits | LinearAccessBit | (functor_traits<UnaryOp>::PacketAccess ? PacketAccessBit : 0)),
570  Alignment = evaluator<ArgType>::Alignment
571  };
572 
573  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit unary_evaluator(const XprType& op) : m_d(op) {
574  EIGEN_INTERNAL_CHECK_COST_VALUE(functor_traits<UnaryOp>::Cost);
575  EIGEN_INTERNAL_CHECK_COST_VALUE(CoeffReadCost);
576  }
577 
578  typedef typename XprType::CoeffReturnType CoeffReturnType;
579 
580  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index row, Index col) const {
581  return m_d.func()(m_d.argImpl.coeff(row, col));
582  }
583 
584  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index index) const {
585  return m_d.func()(m_d.argImpl.coeff(index));
586  }
587 
588  template <int LoadMode, typename PacketType>
589  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index row, Index col) const {
590  return m_d.func().packetOp(m_d.argImpl.template packet<LoadMode, PacketType>(row, col));
591  }
592 
593  template <int LoadMode, typename PacketType>
594  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index index) const {
595  return m_d.func().packetOp(m_d.argImpl.template packet<LoadMode, PacketType>(index));
596  }
597 
598  template <int LoadMode, typename PacketType>
599  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index row, Index col, Index begin, Index count) const {
600  return m_d.func().packetOp(m_d.argImpl.template packetSegment<LoadMode, PacketType>(row, col, begin, count));
601  }
602 
603  template <int LoadMode, typename PacketType>
604  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index index, Index begin, Index count) const {
605  return m_d.func().packetOp(m_d.argImpl.template packetSegment<LoadMode, PacketType>(index, begin, count));
606  }
607 
608  protected:
609  // this helper permits to completely eliminate the functor if it is empty
610  struct Data {
611  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Data(const XprType& xpr)
612  : op(xpr.functor()), argImpl(xpr.nestedExpression()) {}
613  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const UnaryOp& func() const { return op; }
614  UnaryOp op;
615  evaluator<ArgType> argImpl;
616  };
617 
618  Data m_d;
619 };
620 
621 // ----------------------- Casting ---------------------
622 
623 template <typename SrcType, typename DstType, typename ArgType>
624 struct unary_evaluator<CwiseUnaryOp<core_cast_op<SrcType, DstType>, ArgType>, IndexBased> {
625  using CastOp = core_cast_op<SrcType, DstType>;
626  using XprType = CwiseUnaryOp<CastOp, ArgType>;
627 
628  // Use the largest packet type by default
629  using SrcPacketType = typename packet_traits<SrcType>::type;
630  static constexpr int SrcPacketSize = unpacket_traits<SrcPacketType>::size;
631  static constexpr int SrcPacketBytes = SrcPacketSize * sizeof(SrcType);
632 
633  enum {
634  CoeffReadCost = int(evaluator<ArgType>::CoeffReadCost) + int(functor_traits<CastOp>::Cost),
635  PacketAccess = functor_traits<CastOp>::PacketAccess,
636  ActualPacketAccessBit = PacketAccess ? PacketAccessBit : 0,
637  Flags = evaluator<ArgType>::Flags & (HereditaryBits | LinearAccessBit | ActualPacketAccessBit),
638  IsRowMajor = (evaluator<ArgType>::Flags & RowMajorBit),
639  Alignment = evaluator<ArgType>::Alignment
640  };
641 
642  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit unary_evaluator(const XprType& xpr)
643  : m_argImpl(xpr.nestedExpression()), m_rows(xpr.rows()), m_cols(xpr.cols()) {
644  EIGEN_INTERNAL_CHECK_COST_VALUE(functor_traits<CastOp>::Cost);
645  EIGEN_INTERNAL_CHECK_COST_VALUE(CoeffReadCost);
646  }
647 
648  template <typename DstPacketType>
649  using AltSrcScalarOp = std::enable_if_t<(unpacket_traits<DstPacketType>::size < SrcPacketSize &&
650  !find_packet_by_size<SrcType, unpacket_traits<DstPacketType>::size>::value),
651  bool>;
652  template <typename DstPacketType>
653  using SrcPacketArgs1 =
654  std::enable_if_t<(find_packet_by_size<SrcType, unpacket_traits<DstPacketType>::size>::value), bool>;
655  template <typename DstPacketType>
656  using SrcPacketArgs2 = std::enable_if_t<(unpacket_traits<DstPacketType>::size) == (2 * SrcPacketSize), bool>;
657  template <typename DstPacketType>
658  using SrcPacketArgs4 = std::enable_if_t<(unpacket_traits<DstPacketType>::size) == (4 * SrcPacketSize), bool>;
659  template <typename DstPacketType>
660  using SrcPacketArgs8 = std::enable_if_t<(unpacket_traits<DstPacketType>::size) == (8 * SrcPacketSize), bool>;
661 
662  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE bool check_array_bounds(Index row, Index col, Index begin, Index count) const {
663  return IsRowMajor ? (col + count + begin <= cols()) : (row + count + begin <= rows());
664  }
665  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE bool check_array_bounds(Index index, Index begin, Index count) const {
666  return index + count + begin <= size();
667  }
668 
669  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE SrcType srcCoeff(Index row, Index col, Index offset) const {
670  Index actualRow = IsRowMajor ? row : row + offset;
671  Index actualCol = IsRowMajor ? col + offset : col;
672  return m_argImpl.coeff(actualRow, actualCol);
673  }
674  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE SrcType srcCoeff(Index index, Index offset) const {
675  Index actualIndex = index + offset;
676  return m_argImpl.coeff(actualIndex);
677  }
678 
679  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstType coeff(Index row, Index col) const {
680  return cast<SrcType, DstType>(srcCoeff(row, col, 0));
681  }
682  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstType coeff(Index index) const {
683  return cast<SrcType, DstType>(srcCoeff(index, 0));
684  }
685 
686  template <int LoadMode, typename PacketType = SrcPacketType>
687  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType srcPacket(Index row, Index col, Index offset) const {
688  constexpr int PacketSize = unpacket_traits<PacketType>::size;
689  Index packetOffset = offset * PacketSize;
690  Index actualRow = IsRowMajor ? row : row + packetOffset;
691  Index actualCol = IsRowMajor ? col + packetOffset : col;
692  eigen_assert(check_array_bounds(actualRow, actualCol, 0, PacketSize) && "Array index out of bounds");
693  return m_argImpl.template packet<LoadMode, PacketType>(actualRow, actualCol);
694  }
695  template <int LoadMode, typename PacketType = SrcPacketType>
696  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType srcPacket(Index index, Index offset) const {
697  constexpr int PacketSize = unpacket_traits<PacketType>::size;
698  Index packetOffset = offset * PacketSize;
699  Index actualIndex = index + packetOffset;
700  eigen_assert(check_array_bounds(actualIndex, 0, PacketSize) && "Array index out of bounds");
701  return m_argImpl.template packet<LoadMode, PacketType>(actualIndex);
702  }
703  template <int LoadMode, typename PacketType = SrcPacketType>
704  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType srcPacketSegment(Index row, Index col, Index begin, Index count,
705  Index offset) const {
706  constexpr int PacketSize = unpacket_traits<PacketType>::size;
707  Index packetOffset = offset * PacketSize;
708  Index actualRow = IsRowMajor ? row : row + packetOffset;
709  Index actualCol = IsRowMajor ? col + packetOffset : col;
710  eigen_assert(check_array_bounds(actualRow, actualCol, begin, count) && "Array index out of bounds");
711  return m_argImpl.template packetSegment<LoadMode, PacketType>(actualRow, actualCol, begin, count);
712  }
713  template <int LoadMode, typename PacketType = SrcPacketType>
714  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType srcPacketSegment(Index index, Index begin, Index count,
715  Index offset) const {
716  constexpr int PacketSize = unpacket_traits<PacketType>::size;
717  Index packetOffset = offset * PacketSize;
718  Index actualIndex = index + packetOffset;
719  eigen_assert(check_array_bounds(actualIndex, begin, count) && "Array index out of bounds");
720  return m_argImpl.template packetSegment<LoadMode, PacketType>(actualIndex, begin, count);
721  }
722 
723  template <int NumPackets, int LoadMode, typename PacketType = SrcPacketType>
724  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketBlock<PacketType, NumPackets> srcPacketSegmentHelper(Index row, Index col,
725  Index begin,
726  Index count) const {
727  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
728  PacketBlock<PacketType, NumPackets> packets;
729  for (Index i = 0; i < NumPackets; i++) packets.packet[i] = pzero(PacketType());
730  Index offset = begin / SrcPacketSize;
731  Index actualBegin = begin % SrcPacketSize;
732  for (; offset < NumPackets; offset++) {
733  Index actualCount = numext::mini(SrcPacketSize - actualBegin, count);
734  packets.packet[offset] = srcPacketSegment<SrcLoadMode>(row, col, actualBegin, actualCount, offset);
735  if (count == actualCount) break;
736  actualBegin = 0;
737  count -= actualCount;
738  }
739  return packets;
740  }
741  template <int NumPackets, int LoadMode, typename PacketType = SrcPacketType>
742  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketBlock<PacketType, NumPackets> srcPacketSegmentHelper(Index index,
743  Index begin,
744  Index count) const {
745  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
746  PacketBlock<PacketType, NumPackets> packets;
747  for (Index i = 0; i < NumPackets; i++) packets.packet[i] = pzero(PacketType());
748  Index offset = begin / SrcPacketSize;
749  Index actualBegin = begin % SrcPacketSize;
750  for (; offset < NumPackets; offset++) {
751  Index actualCount = numext::mini(SrcPacketSize - actualBegin, count);
752  packets.packet[offset] = srcPacketSegment<SrcLoadMode>(index, actualBegin, actualCount, offset);
753  if (count == actualCount) break;
754  actualBegin = 0;
755  count -= actualCount;
756  }
757  return packets;
758  }
759 
760  // There is no source packet type with equal or fewer elements than DstPacketType.
761  // This is problematic as the evaluation loop may attempt to access data outside the bounds of the array.
762  // For example, consider the cast utilizing pcast<Packet4f,Packet2d> with an array of size 4: {0.0f,1.0f,2.0f,3.0f}.
763  // The first iteration of the evaluation loop will load 16 bytes: {0.0f,1.0f,2.0f,3.0f} and cast to {0.0,1.0}, which
764  // is acceptable. The second iteration will load 16 bytes: {2.0f,3.0f,?,?}, which is outside the bounds of the array.
765  template <int LoadMode, typename DstPacketType, AltSrcScalarOp<DstPacketType> = true>
766  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packet(Index row, Index col) const {
767  constexpr int DstPacketSize = unpacket_traits<DstPacketType>::size;
768  constexpr int SrcBytesIncrement = DstPacketSize * sizeof(SrcType);
769  constexpr int SrcLoadMode = plain_enum_min(SrcBytesIncrement, LoadMode);
770  return pcast<SrcPacketType, DstPacketType>(srcPacketSegment<SrcLoadMode>(row, col, 0, DstPacketSize, 0));
771  }
772  // Use the source packet type with the same size as DstPacketType, if it exists
773  template <int LoadMode, typename DstPacketType, SrcPacketArgs1<DstPacketType> = true>
774  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packet(Index row, Index col) const {
775  constexpr int DstPacketSize = unpacket_traits<DstPacketType>::size;
776  using SizedSrcPacketType = typename find_packet_by_size<SrcType, DstPacketSize>::type;
777  constexpr int SrcBytesIncrement = DstPacketSize * sizeof(SrcType);
778  constexpr int SrcLoadMode = plain_enum_min(SrcBytesIncrement, LoadMode);
779  return pcast<SizedSrcPacketType, DstPacketType>(srcPacket<SrcLoadMode, SizedSrcPacketType>(row, col, 0));
780  }
781  // unpacket_traits<DstPacketType>::size == 2 * SrcPacketSize
782  template <int LoadMode, typename DstPacketType, SrcPacketArgs2<DstPacketType> = true>
783  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packet(Index row, Index col) const {
784  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
785  return pcast<SrcPacketType, DstPacketType>(srcPacket<SrcLoadMode>(row, col, 0),
786  srcPacket<SrcLoadMode>(row, col, 1));
787  }
788  // unpacket_traits<DstPacketType>::size == 4 * SrcPacketSize
789  template <int LoadMode, typename DstPacketType, SrcPacketArgs4<DstPacketType> = true>
790  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packet(Index row, Index col) const {
791  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
792  return pcast<SrcPacketType, DstPacketType>(srcPacket<SrcLoadMode>(row, col, 0), srcPacket<SrcLoadMode>(row, col, 1),
793  srcPacket<SrcLoadMode>(row, col, 2),
794  srcPacket<SrcLoadMode>(row, col, 3));
795  }
796  // unpacket_traits<DstPacketType>::size == 8 * SrcPacketSize
797  template <int LoadMode, typename DstPacketType, SrcPacketArgs8<DstPacketType> = true>
798  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packet(Index row, Index col) const {
799  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
800  return pcast<SrcPacketType, DstPacketType>(
801  srcPacket<SrcLoadMode>(row, col, 0), srcPacket<SrcLoadMode>(row, col, 1), srcPacket<SrcLoadMode>(row, col, 2),
802  srcPacket<SrcLoadMode>(row, col, 3), srcPacket<SrcLoadMode>(row, col, 4), srcPacket<SrcLoadMode>(row, col, 5),
803  srcPacket<SrcLoadMode>(row, col, 6), srcPacket<SrcLoadMode>(row, col, 7));
804  }
805 
806  // packetSegment variants
807  template <int LoadMode, typename DstPacketType, AltSrcScalarOp<DstPacketType> = true>
808  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packetSegment(Index row, Index col, Index begin,
809  Index count) const {
810  constexpr int DstPacketSize = unpacket_traits<DstPacketType>::size;
811  constexpr int SrcBytesIncrement = DstPacketSize * sizeof(SrcType);
812  constexpr int SrcLoadMode = plain_enum_min(SrcBytesIncrement, LoadMode);
813  return pcast<SrcPacketType, DstPacketType>(srcPacketSegment<SrcLoadMode>(row, col, begin, count, 0));
814  }
815  // Use the source packet type with the same size as DstPacketType, if it exists
816  template <int LoadMode, typename DstPacketType, SrcPacketArgs1<DstPacketType> = true>
817  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packetSegment(Index row, Index col, Index begin,
818  Index count) const {
819  constexpr int DstPacketSize = unpacket_traits<DstPacketType>::size;
820  using SizedSrcPacketType = typename find_packet_by_size<SrcType, DstPacketSize>::type;
821  constexpr int SrcBytesIncrement = DstPacketSize * sizeof(SrcType);
822  constexpr int SrcLoadMode = plain_enum_min(SrcBytesIncrement, LoadMode);
823  return pcast<SizedSrcPacketType, DstPacketType>(
824  srcPacketSegment<SrcLoadMode, SizedSrcPacketType>(row, col, begin, count, 0));
825  }
826  // unpacket_traits<DstPacketType>::size == 2 * SrcPacketSize
827  template <int LoadMode, typename DstPacketType, SrcPacketArgs2<DstPacketType> = true>
828  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packetSegment(Index row, Index col, Index begin,
829  Index count) const {
830  constexpr int NumPackets = 2;
831  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
832  PacketBlock<SrcPacketType, NumPackets> packets =
833  srcPacketSegmentHelper<NumPackets, SrcLoadMode>(row, col, begin, count);
834  return pcast<SrcPacketType, DstPacketType>(packets.packet[0], packets.packet[1]);
835  }
836  // unpacket_traits<DstPacketType>::size == 4 * SrcPacketSize
837  template <int LoadMode, typename DstPacketType, SrcPacketArgs4<DstPacketType> = true>
838  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packetSegment(Index row, Index col, Index begin,
839  Index count) const {
840  constexpr int NumPackets = 4;
841  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
842  PacketBlock<SrcPacketType, NumPackets> packets =
843  srcPacketSegmentHelper<NumPackets, SrcLoadMode>(row, col, begin, count);
844  return pcast<SrcPacketType, DstPacketType>(packets.packet[0], packets.packet[1], packets.packet[2],
845  packets.packet[3]);
846  }
847  // unpacket_traits<DstPacketType>::size == 8 * SrcPacketSize
848  template <int LoadMode, typename DstPacketType, SrcPacketArgs8<DstPacketType> = true>
849  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packetSegment(Index row, Index col, Index begin,
850  Index count) const {
851  constexpr int NumPackets = 8;
852  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
853  PacketBlock<SrcPacketType, NumPackets> packets =
854  srcPacketSegmentHelper<NumPackets, SrcLoadMode>(row, col, begin, count);
855  return pcast<SrcPacketType, DstPacketType>(packets.packet[0], packets.packet[1], packets.packet[2],
856  packets.packet[3], packets.packet[4], packets.packet[5],
857  packets.packet[6], packets.packet[7]);
858  }
859 
860  // Analogous routines for linear access.
861  template <int LoadMode, typename DstPacketType, AltSrcScalarOp<DstPacketType> = true>
862  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packet(Index index) const {
863  constexpr int DstPacketSize = unpacket_traits<DstPacketType>::size;
864  constexpr int SrcBytesIncrement = DstPacketSize * sizeof(SrcType);
865  constexpr int SrcLoadMode = plain_enum_min(SrcBytesIncrement, LoadMode);
866  return pcast<SrcPacketType, DstPacketType>(srcPacketSegment<SrcLoadMode>(index, 0, DstPacketSize, 0));
867  }
868  template <int LoadMode, typename DstPacketType, SrcPacketArgs1<DstPacketType> = true>
869  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packet(Index index) const {
870  constexpr int DstPacketSize = unpacket_traits<DstPacketType>::size;
871  using SizedSrcPacketType = typename find_packet_by_size<SrcType, DstPacketSize>::type;
872  constexpr int SrcBytesIncrement = DstPacketSize * sizeof(SrcType);
873  constexpr int SrcLoadMode = plain_enum_min(SrcBytesIncrement, LoadMode);
874  return pcast<SizedSrcPacketType, DstPacketType>(srcPacket<SrcLoadMode, SizedSrcPacketType>(index, 0));
875  }
876  template <int LoadMode, typename DstPacketType, SrcPacketArgs2<DstPacketType> = true>
877  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packet(Index index) const {
878  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
879  return pcast<SrcPacketType, DstPacketType>(srcPacket<SrcLoadMode>(index, 0), srcPacket<SrcLoadMode>(index, 1));
880  }
881  template <int LoadMode, typename DstPacketType, SrcPacketArgs4<DstPacketType> = true>
882  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packet(Index index) const {
883  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
884  return pcast<SrcPacketType, DstPacketType>(srcPacket<SrcLoadMode>(index, 0), srcPacket<SrcLoadMode>(index, 1),
885  srcPacket<SrcLoadMode>(index, 2), srcPacket<SrcLoadMode>(index, 3));
886  }
887  template <int LoadMode, typename DstPacketType, SrcPacketArgs8<DstPacketType> = true>
888  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packet(Index index) const {
889  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
890  return pcast<SrcPacketType, DstPacketType>(srcPacket<SrcLoadMode>(index, 0), srcPacket<SrcLoadMode>(index, 1),
891  srcPacket<SrcLoadMode>(index, 2), srcPacket<SrcLoadMode>(index, 3),
892  srcPacket<SrcLoadMode>(index, 4), srcPacket<SrcLoadMode>(index, 5),
893  srcPacket<SrcLoadMode>(index, 6), srcPacket<SrcLoadMode>(index, 7));
894  }
895 
896  // packetSegment variants
897  template <int LoadMode, typename DstPacketType, AltSrcScalarOp<DstPacketType> = true>
898  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packetSegment(Index index, Index begin, Index count) const {
899  constexpr int DstPacketSize = unpacket_traits<DstPacketType>::size;
900  constexpr int SrcBytesIncrement = DstPacketSize * sizeof(SrcType);
901  constexpr int SrcLoadMode = plain_enum_min(SrcBytesIncrement, LoadMode);
902  return pcast<SrcPacketType, DstPacketType>(srcPacketSegment<SrcLoadMode>(index, begin, count, 0));
903  }
904  // Use the source packet type with the same size as DstPacketType, if it exists
905  template <int LoadMode, typename DstPacketType, SrcPacketArgs1<DstPacketType> = true>
906  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packetSegment(Index index, Index begin, Index count) const {
907  constexpr int DstPacketSize = unpacket_traits<DstPacketType>::size;
908  using SizedSrcPacketType = typename find_packet_by_size<SrcType, DstPacketSize>::type;
909  constexpr int SrcBytesIncrement = DstPacketSize * sizeof(SrcType);
910  constexpr int SrcLoadMode = plain_enum_min(SrcBytesIncrement, LoadMode);
911  return pcast<SizedSrcPacketType, DstPacketType>(
912  srcPacketSegment<SrcLoadMode, SizedSrcPacketType>(index, begin, count, 0));
913  }
914  // unpacket_traits<DstPacketType>::size == 2 * SrcPacketSize
915  template <int LoadMode, typename DstPacketType, SrcPacketArgs2<DstPacketType> = true>
916  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packetSegment(Index index, Index begin, Index count) const {
917  constexpr int NumPackets = 2;
918  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
919  PacketBlock<SrcPacketType, NumPackets> packets =
920  srcPacketSegmentHelper<NumPackets, SrcLoadMode>(index, begin, count);
921  return pcast<SrcPacketType, DstPacketType>(packets.packet[0], packets.packet[1]);
922  }
923  // unpacket_traits<DstPacketType>::size == 4 * SrcPacketSize
924  template <int LoadMode, typename DstPacketType, SrcPacketArgs4<DstPacketType> = true>
925  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packetSegment(Index index, Index begin, Index count) const {
926  constexpr int NumPackets = 4;
927  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
928  PacketBlock<SrcPacketType, NumPackets> packets =
929  srcPacketSegmentHelper<NumPackets, SrcLoadMode>(index, begin, count);
930  return pcast<SrcPacketType, DstPacketType>(packets.packet[0], packets.packet[1], packets.packet[2],
931  packets.packet[3]);
932  }
933  // unpacket_traits<DstPacketType>::size == 8 * SrcPacketSize
934  template <int LoadMode, typename DstPacketType, SrcPacketArgs8<DstPacketType> = true>
935  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DstPacketType packetSegment(Index index, Index begin, Index count) const {
936  constexpr int NumPackets = 8;
937  constexpr int SrcLoadMode = plain_enum_min(SrcPacketBytes, LoadMode);
938  PacketBlock<SrcPacketType, NumPackets> packets =
939  srcPacketSegmentHelper<NumPackets, SrcLoadMode>(index, begin, count);
940  return pcast<SrcPacketType, DstPacketType>(packets.packet[0], packets.packet[1], packets.packet[2],
941  packets.packet[3], packets.packet[4], packets.packet[5],
942  packets.packet[6], packets.packet[7]);
943  }
944 
945  constexpr EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index rows() const { return m_rows; }
946  constexpr EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index cols() const { return m_cols; }
947  constexpr EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index size() const { return m_rows * m_cols; }
948 
949  protected:
950  const evaluator<ArgType> m_argImpl;
951  const variable_if_dynamic<Index, XprType::RowsAtCompileTime> m_rows;
952  const variable_if_dynamic<Index, XprType::ColsAtCompileTime> m_cols;
953 };
954 
955 // -------------------- CwiseTernaryOp --------------------
956 
957 // this is a ternary expression
958 template <typename TernaryOp, typename Arg1, typename Arg2, typename Arg3>
959 struct evaluator<CwiseTernaryOp<TernaryOp, Arg1, Arg2, Arg3>>
960  : public ternary_evaluator<CwiseTernaryOp<TernaryOp, Arg1, Arg2, Arg3>> {
961  typedef CwiseTernaryOp<TernaryOp, Arg1, Arg2, Arg3> XprType;
962  typedef ternary_evaluator<CwiseTernaryOp<TernaryOp, Arg1, Arg2, Arg3>> Base;
963 
964  EIGEN_DEVICE_FUNC explicit evaluator(const XprType& xpr) : Base(xpr) {}
965 };
966 
967 template <typename TernaryOp, typename Arg1, typename Arg2, typename Arg3>
968 struct ternary_evaluator<CwiseTernaryOp<TernaryOp, Arg1, Arg2, Arg3>, IndexBased, IndexBased>
969  : evaluator_base<CwiseTernaryOp<TernaryOp, Arg1, Arg2, Arg3>> {
970  typedef CwiseTernaryOp<TernaryOp, Arg1, Arg2, Arg3> XprType;
971 
972  enum {
973  CoeffReadCost = int(evaluator<Arg1>::CoeffReadCost) + int(evaluator<Arg2>::CoeffReadCost) +
974  int(evaluator<Arg3>::CoeffReadCost) + int(functor_traits<TernaryOp>::Cost),
975 
976  Arg1Flags = evaluator<Arg1>::Flags,
977  Arg2Flags = evaluator<Arg2>::Flags,
978  Arg3Flags = evaluator<Arg3>::Flags,
979  SameType = is_same<typename Arg1::Scalar, typename Arg2::Scalar>::value &&
980  is_same<typename Arg1::Scalar, typename Arg3::Scalar>::value,
981  StorageOrdersAgree = (int(Arg1Flags) & RowMajorBit) == (int(Arg2Flags) & RowMajorBit) &&
982  (int(Arg1Flags) & RowMajorBit) == (int(Arg3Flags) & RowMajorBit),
983  Flags0 = (int(Arg1Flags) | int(Arg2Flags) | int(Arg3Flags)) &
984  (HereditaryBits |
985  (int(Arg1Flags) & int(Arg2Flags) & int(Arg3Flags) &
986  ((StorageOrdersAgree ? LinearAccessBit : 0) |
987  (functor_traits<TernaryOp>::PacketAccess && StorageOrdersAgree && SameType ? PacketAccessBit : 0)))),
988  Flags = (Flags0 & ~RowMajorBit) | (Arg1Flags & RowMajorBit),
989  Alignment = plain_enum_min(plain_enum_min(evaluator<Arg1>::Alignment, evaluator<Arg2>::Alignment),
990  evaluator<Arg3>::Alignment)
991  };
992 
993  EIGEN_DEVICE_FUNC explicit ternary_evaluator(const XprType& xpr) : m_d(xpr) {
994  EIGEN_INTERNAL_CHECK_COST_VALUE(functor_traits<TernaryOp>::Cost);
995  EIGEN_INTERNAL_CHECK_COST_VALUE(CoeffReadCost);
996  }
997 
998  typedef typename XprType::CoeffReturnType CoeffReturnType;
999 
1000  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index row, Index col) const {
1001  return m_d.func()(m_d.arg1Impl.coeff(row, col), m_d.arg2Impl.coeff(row, col), m_d.arg3Impl.coeff(row, col));
1002  }
1003 
1004  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index index) const {
1005  return m_d.func()(m_d.arg1Impl.coeff(index), m_d.arg2Impl.coeff(index), m_d.arg3Impl.coeff(index));
1006  }
1007 
1008  template <int LoadMode, typename PacketType>
1009  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index row, Index col) const {
1010  return m_d.func().packetOp(m_d.arg1Impl.template packet<LoadMode, PacketType>(row, col),
1011  m_d.arg2Impl.template packet<LoadMode, PacketType>(row, col),
1012  m_d.arg3Impl.template packet<LoadMode, PacketType>(row, col));
1013  }
1014 
1015  template <int LoadMode, typename PacketType>
1016  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index index) const {
1017  return m_d.func().packetOp(m_d.arg1Impl.template packet<LoadMode, PacketType>(index),
1018  m_d.arg2Impl.template packet<LoadMode, PacketType>(index),
1019  m_d.arg3Impl.template packet<LoadMode, PacketType>(index));
1020  }
1021 
1022  template <int LoadMode, typename PacketType>
1023  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index row, Index col, Index begin, Index count) const {
1024  return m_d.func().packetOp(m_d.arg1Impl.template packetSegment<LoadMode, PacketType>(row, col, begin, count),
1025  m_d.arg2Impl.template packetSegment<LoadMode, PacketType>(row, col, begin, count),
1026  m_d.arg3Impl.template packetSegment<LoadMode, PacketType>(row, col, begin, count));
1027  }
1028 
1029  template <int LoadMode, typename PacketType>
1030  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index index, Index begin, Index count) const {
1031  return m_d.func().packetOp(m_d.arg1Impl.template packetSegment<LoadMode, PacketType>(index, begin, count),
1032  m_d.arg2Impl.template packetSegment<LoadMode, PacketType>(index, begin, count),
1033  m_d.arg3Impl.template packetSegment<LoadMode, PacketType>(index, begin, count));
1034  }
1035 
1036  protected:
1037  // this helper permits to completely eliminate the functor if it is empty
1038  struct Data {
1039  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Data(const XprType& xpr)
1040  : op(xpr.functor()), arg1Impl(xpr.arg1()), arg2Impl(xpr.arg2()), arg3Impl(xpr.arg3()) {}
1041  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const TernaryOp& func() const { return op; }
1042  TernaryOp op;
1043  evaluator<Arg1> arg1Impl;
1044  evaluator<Arg2> arg2Impl;
1045  evaluator<Arg3> arg3Impl;
1046  };
1047 
1048  Data m_d;
1049 };
1050 
1051 template <typename Arg1, typename Arg2, typename Scalar, typename CmpLhsType, typename CmpRhsType, ComparisonName cmp>
1052 struct scalar_boolean_select_spec {
1053  using DummyTernaryOp = scalar_boolean_select_op<Scalar, Scalar, bool>;
1054  using DummyArg3 = CwiseBinaryOp<scalar_cmp_op<Scalar, Scalar, cmp, false>, CmpLhsType, CmpRhsType>;
1055  using DummyXprType = CwiseTernaryOp<DummyTernaryOp, Arg1, Arg2, DummyArg3>;
1056 
1057  // only use the typed comparison if it is vectorized
1058  static constexpr bool UseTyped = functor_traits<scalar_cmp_op<Scalar, Scalar, cmp, true>>::PacketAccess;
1059  using CondScalar = std::conditional_t<UseTyped, Scalar, bool>;
1060 
1061  using TernaryOp = scalar_boolean_select_op<Scalar, Scalar, CondScalar>;
1062  using Arg3 = CwiseBinaryOp<scalar_cmp_op<Scalar, Scalar, cmp, UseTyped>, CmpLhsType, CmpRhsType>;
1063  using XprType = CwiseTernaryOp<TernaryOp, Arg1, Arg2, Arg3>;
1064 
1065  using Base = ternary_evaluator<XprType>;
1066 };
1067 
1068 // specialization for expressions like (a < b).select(c, d) to enable full vectorization
1069 template <typename Arg1, typename Arg2, typename Scalar, typename CmpLhsType, typename CmpRhsType, ComparisonName cmp>
1070 struct evaluator<CwiseTernaryOp<scalar_boolean_select_op<Scalar, Scalar, bool>, Arg1, Arg2,
1071  CwiseBinaryOp<scalar_cmp_op<Scalar, Scalar, cmp, false>, CmpLhsType, CmpRhsType>>>
1072  : public scalar_boolean_select_spec<Arg1, Arg2, Scalar, CmpLhsType, CmpRhsType, cmp>::Base {
1073  using Helper = scalar_boolean_select_spec<Arg1, Arg2, Scalar, CmpLhsType, CmpRhsType, cmp>;
1074  using Base = typename Helper::Base;
1075  using DummyXprType = typename Helper::DummyXprType;
1076  using Arg3 = typename Helper::Arg3;
1077  using XprType = typename Helper::XprType;
1078 
1079  EIGEN_DEVICE_FUNC explicit evaluator(const DummyXprType& xpr)
1080  : Base(XprType(xpr.arg1(), xpr.arg2(), Arg3(xpr.arg3().lhs(), xpr.arg3().rhs()))) {}
1081 };
1082 
1083 // -------------------- CwiseBinaryOp --------------------
1084 
1085 // this is a binary expression
1086 template <typename BinaryOp, typename Lhs, typename Rhs>
1087 struct evaluator<CwiseBinaryOp<BinaryOp, Lhs, Rhs>> : public binary_evaluator<CwiseBinaryOp<BinaryOp, Lhs, Rhs>> {
1088  typedef CwiseBinaryOp<BinaryOp, Lhs, Rhs> XprType;
1089  typedef binary_evaluator<CwiseBinaryOp<BinaryOp, Lhs, Rhs>> Base;
1090 
1091  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit evaluator(const XprType& xpr) : Base(xpr) {}
1092 };
1093 
1094 template <typename BinaryOp, typename Lhs, typename Rhs>
1095 struct binary_evaluator<CwiseBinaryOp<BinaryOp, Lhs, Rhs>, IndexBased, IndexBased>
1096  : evaluator_base<CwiseBinaryOp<BinaryOp, Lhs, Rhs>> {
1097  typedef CwiseBinaryOp<BinaryOp, Lhs, Rhs> XprType;
1098 
1099  enum {
1100  CoeffReadCost =
1101  int(evaluator<Lhs>::CoeffReadCost) + int(evaluator<Rhs>::CoeffReadCost) + int(functor_traits<BinaryOp>::Cost),
1102 
1103  LhsFlags = evaluator<Lhs>::Flags,
1104  RhsFlags = evaluator<Rhs>::Flags,
1105  SameType = is_same<typename Lhs::Scalar, typename Rhs::Scalar>::value,
1106  StorageOrdersAgree = (int(LhsFlags) & RowMajorBit) == (int(RhsFlags) & RowMajorBit),
1107  Flags0 = (int(LhsFlags) | int(RhsFlags)) &
1108  (HereditaryBits |
1109  (int(LhsFlags) & int(RhsFlags) &
1110  ((StorageOrdersAgree ? LinearAccessBit : 0) |
1111  (functor_traits<BinaryOp>::PacketAccess && StorageOrdersAgree && SameType ? PacketAccessBit : 0)))),
1112  Flags = (Flags0 & ~RowMajorBit) | (LhsFlags & RowMajorBit),
1113  Alignment = plain_enum_min(evaluator<Lhs>::Alignment, evaluator<Rhs>::Alignment)
1114  };
1115 
1116  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit binary_evaluator(const XprType& xpr) : m_d(xpr) {
1117  EIGEN_INTERNAL_CHECK_COST_VALUE(functor_traits<BinaryOp>::Cost);
1118  EIGEN_INTERNAL_CHECK_COST_VALUE(CoeffReadCost);
1119  }
1120 
1121  typedef typename XprType::CoeffReturnType CoeffReturnType;
1122 
1123  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index row, Index col) const {
1124  return m_d.func()(m_d.lhsImpl.coeff(row, col), m_d.rhsImpl.coeff(row, col));
1125  }
1126 
1127  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index index) const {
1128  return m_d.func()(m_d.lhsImpl.coeff(index), m_d.rhsImpl.coeff(index));
1129  }
1130 
1131  template <int LoadMode, typename PacketType>
1132  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index row, Index col) const {
1133  return m_d.func().packetOp(m_d.lhsImpl.template packet<LoadMode, PacketType>(row, col),
1134  m_d.rhsImpl.template packet<LoadMode, PacketType>(row, col));
1135  }
1136 
1137  template <int LoadMode, typename PacketType>
1138  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index index) const {
1139  return m_d.func().packetOp(m_d.lhsImpl.template packet<LoadMode, PacketType>(index),
1140  m_d.rhsImpl.template packet<LoadMode, PacketType>(index));
1141  }
1142 
1143  template <int LoadMode, typename PacketType>
1144  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index row, Index col, Index begin, Index count) const {
1145  return m_d.func().packetOp(m_d.lhsImpl.template packetSegment<LoadMode, PacketType>(row, col, begin, count),
1146  m_d.rhsImpl.template packetSegment<LoadMode, PacketType>(row, col, begin, count));
1147  }
1148 
1149  template <int LoadMode, typename PacketType>
1150  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index index, Index begin, Index count) const {
1151  return m_d.func().packetOp(m_d.lhsImpl.template packetSegment<LoadMode, PacketType>(index, begin, count),
1152  m_d.rhsImpl.template packetSegment<LoadMode, PacketType>(index, begin, count));
1153  }
1154 
1155  protected:
1156  // this helper permits to completely eliminate the functor if it is empty
1157  struct Data {
1158  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Data(const XprType& xpr)
1159  : op(xpr.functor()), lhsImpl(xpr.lhs()), rhsImpl(xpr.rhs()) {}
1160  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const BinaryOp& func() const { return op; }
1161  BinaryOp op;
1162  evaluator<Lhs> lhsImpl;
1163  evaluator<Rhs> rhsImpl;
1164  };
1165 
1166  Data m_d;
1167 };
1168 
1169 // -------------------- CwiseUnaryView --------------------
1170 
1171 template <typename UnaryOp, typename ArgType, typename StrideType>
1172 struct unary_evaluator<CwiseUnaryView<UnaryOp, ArgType, StrideType>, IndexBased>
1173  : evaluator_base<CwiseUnaryView<UnaryOp, ArgType, StrideType>> {
1174  typedef CwiseUnaryView<UnaryOp, ArgType, StrideType> XprType;
1175 
1176  enum {
1177  CoeffReadCost = int(evaluator<ArgType>::CoeffReadCost) + int(functor_traits<UnaryOp>::Cost),
1178 
1179  Flags = (evaluator<ArgType>::Flags & (HereditaryBits | LinearAccessBit | DirectAccessBit)),
1180 
1181  Alignment = 0 // FIXME it is not very clear why alignment is necessarily lost...
1182  };
1183 
1184  EIGEN_DEVICE_FUNC explicit unary_evaluator(const XprType& op) : m_d(op) {
1185  EIGEN_INTERNAL_CHECK_COST_VALUE(functor_traits<UnaryOp>::Cost);
1186  EIGEN_INTERNAL_CHECK_COST_VALUE(CoeffReadCost);
1187  }
1188 
1189  typedef typename XprType::Scalar Scalar;
1190  typedef typename XprType::CoeffReturnType CoeffReturnType;
1191 
1192  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index row, Index col) const {
1193  return m_d.func()(m_d.argImpl.coeff(row, col));
1194  }
1195 
1196  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index index) const {
1197  return m_d.func()(m_d.argImpl.coeff(index));
1198  }
1199 
1200  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index row, Index col) {
1201  return m_d.func()(m_d.argImpl.coeffRef(row, col));
1202  }
1203 
1204  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index index) {
1205  return m_d.func()(m_d.argImpl.coeffRef(index));
1206  }
1207 
1208  protected:
1209  // this helper permits to completely eliminate the functor if it is empty
1210  struct Data {
1211  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Data(const XprType& xpr)
1212  : op(xpr.functor()), argImpl(xpr.nestedExpression()) {}
1213  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const UnaryOp& func() const { return op; }
1214  UnaryOp op;
1215  evaluator<ArgType> argImpl;
1216  };
1217 
1218  Data m_d;
1219 };
1220 
1221 // -------------------- Map --------------------
1222 
1223 // FIXME perhaps the PlainObjectType could be provided by Derived::PlainObject ?
1224 // but that might complicate template specialization
1225 template <typename Derived, typename PlainObjectType>
1226 struct mapbase_evaluator;
1227 
1228 template <typename Derived, typename PlainObjectType>
1229 struct mapbase_evaluator : evaluator_base<Derived> {
1230  typedef Derived XprType;
1231  typedef typename XprType::PointerType PointerType;
1232  typedef typename XprType::Scalar Scalar;
1233  typedef typename XprType::CoeffReturnType CoeffReturnType;
1234 
1235  enum {
1236  IsRowMajor = XprType::RowsAtCompileTime,
1237  ColsAtCompileTime = XprType::ColsAtCompileTime,
1238  CoeffReadCost = NumTraits<Scalar>::ReadCost
1239  };
1240 
1241  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit mapbase_evaluator(const XprType& map)
1242  : m_data(const_cast<PointerType>(map.data())),
1243  m_innerStride(map.innerStride()),
1244  m_outerStride(map.outerStride()) {
1245  EIGEN_STATIC_ASSERT(check_implication((evaluator<Derived>::Flags & PacketAccessBit) != 0,
1246  inner_stride_at_compile_time<Derived>::ret == 1),
1247  PACKET_ACCESS_REQUIRES_TO_HAVE_INNER_STRIDE_FIXED_TO_1);
1248  EIGEN_INTERNAL_CHECK_COST_VALUE(CoeffReadCost);
1249  }
1250 
1251  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index row, Index col) const {
1252  return m_data[col * colStride() + row * rowStride()];
1253  }
1254 
1255  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index index) const {
1256  return m_data[index * m_innerStride.value()];
1257  }
1258 
1259  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index row, Index col) {
1260  return m_data[col * colStride() + row * rowStride()];
1261  }
1262 
1263  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index index) { return m_data[index * m_innerStride.value()]; }
1264 
1265  template <int LoadMode, typename PacketType>
1266  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index row, Index col) const {
1267  PointerType ptr = m_data + row * rowStride() + col * colStride();
1268  return ploadt<PacketType, LoadMode>(ptr);
1269  }
1270 
1271  template <int LoadMode, typename PacketType>
1272  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index index) const {
1273  return ploadt<PacketType, LoadMode>(m_data + index * m_innerStride.value());
1274  }
1275 
1276  template <int StoreMode, typename PacketType>
1277  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacket(Index row, Index col, const PacketType& x) {
1278  PointerType ptr = m_data + row * rowStride() + col * colStride();
1279  pstoret<Scalar, PacketType, StoreMode>(ptr, x);
1280  }
1281 
1282  template <int StoreMode, typename PacketType>
1283  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacket(Index index, const PacketType& x) {
1284  pstoret<Scalar, PacketType, StoreMode>(m_data + index * m_innerStride.value(), x);
1285  }
1286 
1287  template <int LoadMode, typename PacketType>
1288  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index row, Index col, Index begin, Index count) const {
1289  PointerType ptr = m_data + row * rowStride() + col * colStride();
1290  return ploadtSegment<PacketType, LoadMode>(ptr, begin, count);
1291  }
1292 
1293  template <int LoadMode, typename PacketType>
1294  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index index, Index begin, Index count) const {
1295  return ploadtSegment<PacketType, LoadMode>(m_data + index * m_innerStride.value(), begin, count);
1296  }
1297 
1298  template <int StoreMode, typename PacketType>
1299  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacketSegment(Index row, Index col, const PacketType& x, Index begin,
1300  Index count) {
1301  PointerType ptr = m_data + row * rowStride() + col * colStride();
1302  pstoretSegment<Scalar, PacketType, StoreMode>(ptr, x, begin, count);
1303  }
1304 
1305  template <int StoreMode, typename PacketType>
1306  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacketSegment(Index index, const PacketType& x, Index begin,
1307  Index count) {
1308  pstoretSegment<Scalar, PacketType, StoreMode>(m_data + index * m_innerStride.value(), x, begin, count);
1309  }
1310 
1311  protected:
1312  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr Index rowStride() const noexcept {
1313  return XprType::IsRowMajor ? m_outerStride.value() : m_innerStride.value();
1314  }
1315  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr Index colStride() const noexcept {
1316  return XprType::IsRowMajor ? m_innerStride.value() : m_outerStride.value();
1317  }
1318 
1319  PointerType m_data;
1320  const variable_if_dynamic<Index, XprType::InnerStrideAtCompileTime> m_innerStride;
1321  const variable_if_dynamic<Index, XprType::OuterStrideAtCompileTime> m_outerStride;
1322 };
1323 
1324 template <typename PlainObjectType, int MapOptions, typename StrideType>
1325 struct evaluator<Map<PlainObjectType, MapOptions, StrideType>>
1326  : public mapbase_evaluator<Map<PlainObjectType, MapOptions, StrideType>, PlainObjectType> {
1327  typedef Map<PlainObjectType, MapOptions, StrideType> XprType;
1328  typedef typename XprType::Scalar Scalar;
1329  // TODO: should check for smaller packet types once we can handle multi-sized packet types
1330  typedef typename packet_traits<Scalar>::type PacketScalar;
1331 
1332  enum {
1333  InnerStrideAtCompileTime = StrideType::InnerStrideAtCompileTime == 0
1334  ? int(PlainObjectType::InnerStrideAtCompileTime)
1335  : int(StrideType::InnerStrideAtCompileTime),
1336  OuterStrideAtCompileTime = StrideType::OuterStrideAtCompileTime == 0
1337  ? int(PlainObjectType::OuterStrideAtCompileTime)
1338  : int(StrideType::OuterStrideAtCompileTime),
1339  HasNoInnerStride = InnerStrideAtCompileTime == 1,
1340  HasNoOuterStride = StrideType::OuterStrideAtCompileTime == 0,
1341  HasNoStride = HasNoInnerStride && HasNoOuterStride,
1342  IsDynamicSize = PlainObjectType::SizeAtCompileTime == Dynamic,
1343 
1344  PacketAccessMask = bool(HasNoInnerStride) ? ~int(0) : ~int(PacketAccessBit),
1345  LinearAccessMask =
1346  bool(HasNoStride) || bool(PlainObjectType::IsVectorAtCompileTime) ? ~int(0) : ~int(LinearAccessBit),
1347  Flags = int(evaluator<PlainObjectType>::Flags) & (LinearAccessMask & PacketAccessMask),
1348 
1349  Alignment = int(MapOptions) & int(AlignedMask)
1350  };
1351 
1352  EIGEN_DEVICE_FUNC explicit evaluator(const XprType& map) : mapbase_evaluator<XprType, PlainObjectType>(map) {}
1353 };
1354 
1355 // -------------------- Ref --------------------
1356 
1357 template <typename PlainObjectType, int RefOptions, typename StrideType>
1358 struct evaluator<Ref<PlainObjectType, RefOptions, StrideType>>
1359  : public mapbase_evaluator<Ref<PlainObjectType, RefOptions, StrideType>, PlainObjectType> {
1360  typedef Ref<PlainObjectType, RefOptions, StrideType> XprType;
1361 
1362  enum {
1363  Flags = evaluator<Map<PlainObjectType, RefOptions, StrideType>>::Flags,
1364  Alignment = evaluator<Map<PlainObjectType, RefOptions, StrideType>>::Alignment
1365  };
1366 
1367  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit evaluator(const XprType& ref)
1368  : mapbase_evaluator<XprType, PlainObjectType>(ref) {}
1369 };
1370 
1371 // -------------------- Block --------------------
1372 
1373 template <typename ArgType, int BlockRows, int BlockCols, bool InnerPanel,
1374  bool HasDirectAccess = has_direct_access<ArgType>::ret>
1375 struct block_evaluator;
1376 
1377 template <typename ArgType, int BlockRows, int BlockCols, bool InnerPanel>
1378 struct evaluator<Block<ArgType, BlockRows, BlockCols, InnerPanel>>
1379  : block_evaluator<ArgType, BlockRows, BlockCols, InnerPanel> {
1380  typedef Block<ArgType, BlockRows, BlockCols, InnerPanel> XprType;
1381  typedef typename XprType::Scalar Scalar;
1382  // TODO: should check for smaller packet types once we can handle multi-sized packet types
1383  typedef typename packet_traits<Scalar>::type PacketScalar;
1384 
1385  enum {
1386  CoeffReadCost = evaluator<ArgType>::CoeffReadCost,
1387 
1388  RowsAtCompileTime = traits<XprType>::RowsAtCompileTime,
1389  ColsAtCompileTime = traits<XprType>::ColsAtCompileTime,
1390  MaxRowsAtCompileTime = traits<XprType>::MaxRowsAtCompileTime,
1391  MaxColsAtCompileTime = traits<XprType>::MaxColsAtCompileTime,
1392 
1393  ArgTypeIsRowMajor = (int(evaluator<ArgType>::Flags) & RowMajorBit) != 0,
1394  IsRowMajor = (MaxRowsAtCompileTime == 1 && MaxColsAtCompileTime != 1) ? 1
1395  : (MaxColsAtCompileTime == 1 && MaxRowsAtCompileTime != 1) ? 0
1396  : ArgTypeIsRowMajor,
1397  HasSameStorageOrderAsArgType = (IsRowMajor == ArgTypeIsRowMajor),
1398  InnerSize = IsRowMajor ? int(ColsAtCompileTime) : int(RowsAtCompileTime),
1399  InnerStrideAtCompileTime = HasSameStorageOrderAsArgType ? int(inner_stride_at_compile_time<ArgType>::ret)
1400  : int(outer_stride_at_compile_time<ArgType>::ret),
1401  OuterStrideAtCompileTime = HasSameStorageOrderAsArgType ? int(outer_stride_at_compile_time<ArgType>::ret)
1402  : int(inner_stride_at_compile_time<ArgType>::ret),
1403  MaskPacketAccessBit = (InnerStrideAtCompileTime == 1 || HasSameStorageOrderAsArgType) ? PacketAccessBit : 0,
1404 
1405  FlagsLinearAccessBit = (RowsAtCompileTime == 1 || ColsAtCompileTime == 1 ||
1406  (InnerPanel && (evaluator<ArgType>::Flags & LinearAccessBit)))
1407  ? LinearAccessBit
1408  : 0,
1409  FlagsRowMajorBit = XprType::Flags & RowMajorBit,
1410  Flags0 = evaluator<ArgType>::Flags & ((HereditaryBits & ~RowMajorBit) | DirectAccessBit | MaskPacketAccessBit),
1411  Flags = Flags0 | FlagsLinearAccessBit | FlagsRowMajorBit,
1412 
1413  PacketAlignment = unpacket_traits<PacketScalar>::alignment,
1414  Alignment0 = (InnerPanel && (OuterStrideAtCompileTime != Dynamic) && (OuterStrideAtCompileTime != 0) &&
1415  (((OuterStrideAtCompileTime * int(sizeof(Scalar))) % int(PacketAlignment)) == 0))
1416  ? int(PacketAlignment)
1417  : 0,
1418  Alignment = plain_enum_min(evaluator<ArgType>::Alignment, Alignment0)
1419  };
1420  typedef block_evaluator<ArgType, BlockRows, BlockCols, InnerPanel> block_evaluator_type;
1421  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit evaluator(const XprType& block) : block_evaluator_type(block) {
1422  EIGEN_INTERNAL_CHECK_COST_VALUE(CoeffReadCost);
1423  }
1424 };
1425 
1426 // no direct-access => dispatch to a unary evaluator
1427 template <typename ArgType, int BlockRows, int BlockCols, bool InnerPanel>
1428 struct block_evaluator<ArgType, BlockRows, BlockCols, InnerPanel, /*HasDirectAccess*/ false>
1429  : unary_evaluator<Block<ArgType, BlockRows, BlockCols, InnerPanel>> {
1430  typedef Block<ArgType, BlockRows, BlockCols, InnerPanel> XprType;
1431 
1432  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit block_evaluator(const XprType& block)
1433  : unary_evaluator<XprType>(block) {}
1434 };
1435 
1436 template <typename ArgType, int BlockRows, int BlockCols, bool InnerPanel>
1437 struct unary_evaluator<Block<ArgType, BlockRows, BlockCols, InnerPanel>, IndexBased>
1438  : evaluator_base<Block<ArgType, BlockRows, BlockCols, InnerPanel>> {
1439  typedef Block<ArgType, BlockRows, BlockCols, InnerPanel> XprType;
1440 
1441  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit unary_evaluator(const XprType& block)
1442  : m_argImpl(block.nestedExpression()),
1443  m_startRow(block.startRow()),
1444  m_startCol(block.startCol()),
1445  m_linear_offset(ForwardLinearAccess
1446  ? (ArgType::IsRowMajor
1447  ? block.startRow() * block.nestedExpression().cols() + block.startCol()
1448  : block.startCol() * block.nestedExpression().rows() + block.startRow())
1449  : 0) {}
1450 
1451  typedef typename XprType::Scalar Scalar;
1452  typedef typename XprType::CoeffReturnType CoeffReturnType;
1453 
1454  enum {
1455  RowsAtCompileTime = XprType::RowsAtCompileTime,
1456  ForwardLinearAccess = (InnerPanel || int(XprType::IsRowMajor) == int(ArgType::IsRowMajor)) &&
1457  bool(evaluator<ArgType>::Flags & LinearAccessBit)
1458  };
1459 
1460  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index row, Index col) const {
1461  return m_argImpl.coeff(m_startRow.value() + row, m_startCol.value() + col);
1462  }
1463 
1464  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index index) const {
1465  return linear_coeff_impl(index, bool_constant<ForwardLinearAccess>());
1466  }
1467 
1468  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index row, Index col) {
1469  return m_argImpl.coeffRef(m_startRow.value() + row, m_startCol.value() + col);
1470  }
1471 
1472  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index index) {
1473  return linear_coeffRef_impl(index, bool_constant<ForwardLinearAccess>());
1474  }
1475 
1476  template <int LoadMode, typename PacketType>
1477  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index row, Index col) const {
1478  return m_argImpl.template packet<LoadMode, PacketType>(m_startRow.value() + row, m_startCol.value() + col);
1479  }
1480 
1481  template <int LoadMode, typename PacketType>
1482  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index index) const {
1483  if (ForwardLinearAccess)
1484  return m_argImpl.template packet<LoadMode, PacketType>(m_linear_offset.value() + index);
1485  else
1486  return packet<LoadMode, PacketType>(RowsAtCompileTime == 1 ? 0 : index, RowsAtCompileTime == 1 ? index : 0);
1487  }
1488 
1489  template <int StoreMode, typename PacketType>
1490  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacket(Index row, Index col, const PacketType& x) {
1491  return m_argImpl.template writePacket<StoreMode, PacketType>(m_startRow.value() + row, m_startCol.value() + col, x);
1492  }
1493 
1494  template <int StoreMode, typename PacketType>
1495  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacket(Index index, const PacketType& x) {
1496  if (ForwardLinearAccess)
1497  return m_argImpl.template writePacket<StoreMode, PacketType>(m_linear_offset.value() + index, x);
1498  else
1499  return writePacket<StoreMode, PacketType>(RowsAtCompileTime == 1 ? 0 : index, RowsAtCompileTime == 1 ? index : 0,
1500  x);
1501  }
1502 
1503  template <int LoadMode, typename PacketType>
1504  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index row, Index col, Index begin, Index count) const {
1505  return m_argImpl.template packetSegment<LoadMode, PacketType>(m_startRow.value() + row, m_startCol.value() + col,
1506  begin, count);
1507  }
1508 
1509  template <int LoadMode, typename PacketType>
1510  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index index, Index begin, Index count) const {
1511  if (ForwardLinearAccess)
1512  return m_argImpl.template packetSegment<LoadMode, PacketType>(m_linear_offset.value() + index, begin, count);
1513  else
1514  return packetSegment<LoadMode, PacketType>(RowsAtCompileTime == 1 ? 0 : index, RowsAtCompileTime == 1 ? index : 0,
1515  begin, count);
1516  }
1517 
1518  template <int StoreMode, typename PacketType>
1519  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacketSegment(Index row, Index col, const PacketType& x, Index begin,
1520  Index count) {
1521  return m_argImpl.template writePacketSegment<StoreMode, PacketType>(m_startRow.value() + row,
1522  m_startCol.value() + col, x, begin, count);
1523  }
1524 
1525  template <int StoreMode, typename PacketType>
1526  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacketSegment(Index index, const PacketType& x, Index begin,
1527  Index count) {
1528  if (ForwardLinearAccess)
1529  return m_argImpl.template writePacketSegment<StoreMode, PacketType>(m_linear_offset.value() + index, x, begin,
1530  count);
1531  else
1532  return writePacketSegment<StoreMode, PacketType>(RowsAtCompileTime == 1 ? 0 : index,
1533  RowsAtCompileTime == 1 ? index : 0, x, begin, count);
1534  }
1535 
1536  protected:
1537  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType
1538  linear_coeff_impl(Index index, internal::true_type /* ForwardLinearAccess */) const {
1539  return m_argImpl.coeff(m_linear_offset.value() + index);
1540  }
1541  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType
1542  linear_coeff_impl(Index index, internal::false_type /* not ForwardLinearAccess */) const {
1543  return coeff(RowsAtCompileTime == 1 ? 0 : index, RowsAtCompileTime == 1 ? index : 0);
1544  }
1545 
1546  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& linear_coeffRef_impl(Index index,
1547  internal::true_type /* ForwardLinearAccess */) {
1548  return m_argImpl.coeffRef(m_linear_offset.value() + index);
1549  }
1550  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& linear_coeffRef_impl(
1551  Index index, internal::false_type /* not ForwardLinearAccess */) {
1552  return coeffRef(RowsAtCompileTime == 1 ? 0 : index, RowsAtCompileTime == 1 ? index : 0);
1553  }
1554 
1555  evaluator<ArgType> m_argImpl;
1556  const variable_if_dynamic<Index, (ArgType::RowsAtCompileTime == 1 && BlockRows == 1) ? 0 : Dynamic> m_startRow;
1557  const variable_if_dynamic<Index, (ArgType::ColsAtCompileTime == 1 && BlockCols == 1) ? 0 : Dynamic> m_startCol;
1558  const variable_if_dynamic<Index, ForwardLinearAccess ? Dynamic : 0> m_linear_offset;
1559 };
1560 
1561 // TODO: This evaluator does not actually use the child evaluator;
1562 // all action is via the data() as returned by the Block expression.
1563 
1564 template <typename ArgType, int BlockRows, int BlockCols, bool InnerPanel>
1565 struct block_evaluator<ArgType, BlockRows, BlockCols, InnerPanel, /* HasDirectAccess */ true>
1566  : mapbase_evaluator<Block<ArgType, BlockRows, BlockCols, InnerPanel>,
1567  typename Block<ArgType, BlockRows, BlockCols, InnerPanel>::PlainObject> {
1568  typedef Block<ArgType, BlockRows, BlockCols, InnerPanel> XprType;
1569  typedef typename XprType::Scalar Scalar;
1570 
1571  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit block_evaluator(const XprType& block)
1572  : mapbase_evaluator<XprType, typename XprType::PlainObject>(block) {
1573  eigen_internal_assert((internal::is_constant_evaluated() ||
1574  (std::uintptr_t(block.data()) % plain_enum_max(1, evaluator<XprType>::Alignment)) == 0) &&
1575  "data is not aligned");
1576  }
1577 };
1578 
1579 // -------------------- Replicate --------------------
1580 
1581 template <typename ArgType, int RowFactor, int ColFactor>
1582 struct unary_evaluator<Replicate<ArgType, RowFactor, ColFactor>>
1583  : evaluator_base<Replicate<ArgType, RowFactor, ColFactor>> {
1584  typedef Replicate<ArgType, RowFactor, ColFactor> XprType;
1585  typedef typename XprType::CoeffReturnType CoeffReturnType;
1586  enum { Factor = (RowFactor == Dynamic || ColFactor == Dynamic) ? Dynamic : RowFactor * ColFactor };
1587  typedef typename nested_eval<ArgType, Factor>::type ArgTypeNested;
1588  typedef remove_all_t<ArgTypeNested> ArgTypeNestedCleaned;
1589 
1590  enum {
1591  CoeffReadCost = evaluator<ArgTypeNestedCleaned>::CoeffReadCost,
1592  LinearAccessMask = XprType::IsVectorAtCompileTime ? LinearAccessBit : 0,
1593  Flags = (evaluator<ArgTypeNestedCleaned>::Flags & (HereditaryBits | LinearAccessMask) & ~RowMajorBit) |
1594  (traits<XprType>::Flags & RowMajorBit),
1595 
1596  Alignment = evaluator<ArgTypeNestedCleaned>::Alignment
1597  };
1598 
1599  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit unary_evaluator(const XprType& replicate)
1600  : m_arg(replicate.nestedExpression()),
1601  m_argImpl(m_arg),
1602  m_rows(replicate.nestedExpression().rows()),
1603  m_cols(replicate.nestedExpression().cols()) {}
1604 
1605  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index row, Index col) const {
1606  // try to avoid using modulo; this is a pure optimization strategy
1607  const Index actual_row = traits<XprType>::RowsAtCompileTime == 1 ? 0 : RowFactor == 1 ? row : row % m_rows.value();
1608  const Index actual_col = traits<XprType>::ColsAtCompileTime == 1 ? 0 : ColFactor == 1 ? col : col % m_cols.value();
1609 
1610  return m_argImpl.coeff(actual_row, actual_col);
1611  }
1612 
1613  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index index) const {
1614  // try to avoid using modulo; this is a pure optimization strategy
1615  const Index actual_index = traits<XprType>::RowsAtCompileTime == 1
1616  ? (ColFactor == 1 ? index : index % m_cols.value())
1617  : (RowFactor == 1 ? index : index % m_rows.value());
1618 
1619  return m_argImpl.coeff(actual_index);
1620  }
1621 
1622  template <int LoadMode, typename PacketType>
1623  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index row, Index col) const {
1624  const Index actual_row = traits<XprType>::RowsAtCompileTime == 1 ? 0 : RowFactor == 1 ? row : row % m_rows.value();
1625  const Index actual_col = traits<XprType>::ColsAtCompileTime == 1 ? 0 : ColFactor == 1 ? col : col % m_cols.value();
1626 
1627  return m_argImpl.template packet<LoadMode, PacketType>(actual_row, actual_col);
1628  }
1629 
1630  template <int LoadMode, typename PacketType>
1631  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index index) const {
1632  const Index actual_index = traits<XprType>::RowsAtCompileTime == 1
1633  ? (ColFactor == 1 ? index : index % m_cols.value())
1634  : (RowFactor == 1 ? index : index % m_rows.value());
1635 
1636  return m_argImpl.template packet<LoadMode, PacketType>(actual_index);
1637  }
1638 
1639  template <int LoadMode, typename PacketType>
1640  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index row, Index col, Index begin, Index count) const {
1641  const Index actual_row = traits<XprType>::RowsAtCompileTime == 1 ? 0 : RowFactor == 1 ? row : row % m_rows.value();
1642  const Index actual_col = traits<XprType>::ColsAtCompileTime == 1 ? 0 : ColFactor == 1 ? col : col % m_cols.value();
1643 
1644  return m_argImpl.template packetSegment<LoadMode, PacketType>(actual_row, actual_col, begin, count);
1645  }
1646 
1647  template <int LoadMode, typename PacketType>
1648  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index index, Index begin, Index count) const {
1649  const Index actual_index = traits<XprType>::RowsAtCompileTime == 1
1650  ? (ColFactor == 1 ? index : index % m_cols.value())
1651  : (RowFactor == 1 ? index : index % m_rows.value());
1652 
1653  return m_argImpl.template packetSegment<LoadMode, PacketType>(actual_index, begin, count);
1654  }
1655 
1656  protected:
1657  const ArgTypeNested m_arg;
1658  evaluator<ArgTypeNestedCleaned> m_argImpl;
1659  const variable_if_dynamic<Index, ArgType::RowsAtCompileTime> m_rows;
1660  const variable_if_dynamic<Index, ArgType::ColsAtCompileTime> m_cols;
1661 };
1662 
1663 // -------------------- MatrixWrapper and ArrayWrapper --------------------
1664 //
1665 // evaluator_wrapper_base<T> is a common base class for the
1666 // MatrixWrapper and ArrayWrapper evaluators.
1667 
1668 template <typename XprType>
1669 struct evaluator_wrapper_base : evaluator_base<XprType> {
1670  typedef remove_all_t<typename XprType::NestedExpressionType> ArgType;
1671  enum {
1672  CoeffReadCost = evaluator<ArgType>::CoeffReadCost,
1673  Flags = evaluator<ArgType>::Flags,
1674  Alignment = evaluator<ArgType>::Alignment
1675  };
1676 
1677  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit evaluator_wrapper_base(const ArgType& arg) : m_argImpl(arg) {}
1678 
1679  typedef typename ArgType::Scalar Scalar;
1680  typedef typename ArgType::CoeffReturnType CoeffReturnType;
1681 
1682  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index row, Index col) const {
1683  return m_argImpl.coeff(row, col);
1684  }
1685 
1686  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index index) const { return m_argImpl.coeff(index); }
1687 
1688  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index row, Index col) { return m_argImpl.coeffRef(row, col); }
1689 
1690  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index index) { return m_argImpl.coeffRef(index); }
1691 
1692  template <int LoadMode, typename PacketType>
1693  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index row, Index col) const {
1694  return m_argImpl.template packet<LoadMode, PacketType>(row, col);
1695  }
1696 
1697  template <int LoadMode, typename PacketType>
1698  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index index) const {
1699  return m_argImpl.template packet<LoadMode, PacketType>(index);
1700  }
1701 
1702  template <int StoreMode, typename PacketType>
1703  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacket(Index row, Index col, const PacketType& x) {
1704  m_argImpl.template writePacket<StoreMode>(row, col, x);
1705  }
1706 
1707  template <int StoreMode, typename PacketType>
1708  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacket(Index index, const PacketType& x) {
1709  m_argImpl.template writePacket<StoreMode>(index, x);
1710  }
1711 
1712  template <int LoadMode, typename PacketType>
1713  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index row, Index col, Index begin, Index count) const {
1714  return m_argImpl.template packetSegment<LoadMode, PacketType>(row, col, begin, count);
1715  }
1716 
1717  template <int LoadMode, typename PacketType>
1718  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index index, Index begin, Index count) const {
1719  return m_argImpl.template packetSegment<LoadMode, PacketType>(index, begin, count);
1720  }
1721 
1722  template <int StoreMode, typename PacketType>
1723  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacketSegment(Index row, Index col, const PacketType& x, Index begin,
1724  Index count) {
1725  m_argImpl.template writePacketSegment<StoreMode>(row, col, x, begin, count);
1726  }
1727 
1728  template <int StoreMode, typename PacketType>
1729  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacketSegment(Index index, const PacketType& x, Index begin,
1730  Index count) {
1731  m_argImpl.template writePacketSegment<StoreMode>(index, x, begin, count);
1732  }
1733 
1734  protected:
1735  evaluator<ArgType> m_argImpl;
1736 };
1737 
1738 template <typename TArgType>
1739 struct unary_evaluator<MatrixWrapper<TArgType>> : evaluator_wrapper_base<MatrixWrapper<TArgType>> {
1740  typedef MatrixWrapper<TArgType> XprType;
1741 
1742  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit unary_evaluator(const XprType& wrapper)
1743  : evaluator_wrapper_base<MatrixWrapper<TArgType>>(wrapper.nestedExpression()) {}
1744 };
1745 
1746 template <typename TArgType>
1747 struct unary_evaluator<ArrayWrapper<TArgType>> : evaluator_wrapper_base<ArrayWrapper<TArgType>> {
1748  typedef ArrayWrapper<TArgType> XprType;
1749 
1750  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit unary_evaluator(const XprType& wrapper)
1751  : evaluator_wrapper_base<ArrayWrapper<TArgType>>(wrapper.nestedExpression()) {}
1752 };
1753 
1754 // -------------------- Reverse --------------------
1755 
1756 // defined in Reverse.h:
1757 template <typename PacketType, bool ReversePacket>
1758 struct reverse_packet_cond;
1759 
1760 template <typename ArgType, int Direction>
1761 struct unary_evaluator<Reverse<ArgType, Direction>> : evaluator_base<Reverse<ArgType, Direction>> {
1762  typedef Reverse<ArgType, Direction> XprType;
1763  typedef typename XprType::Scalar Scalar;
1764  typedef typename XprType::CoeffReturnType CoeffReturnType;
1765 
1766  enum {
1767  IsRowMajor = XprType::IsRowMajor,
1768  IsColMajor = !IsRowMajor,
1769  ReverseRow = (Direction == Vertical) || (Direction == BothDirections),
1770  ReverseCol = (Direction == Horizontal) || (Direction == BothDirections),
1771  ReversePacket = (Direction == BothDirections) || ((Direction == Vertical) && IsColMajor) ||
1772  ((Direction == Horizontal) && IsRowMajor),
1773 
1774  CoeffReadCost = evaluator<ArgType>::CoeffReadCost,
1775 
1776  // let's enable LinearAccess only with vectorization because of the product overhead
1777  // FIXME enable DirectAccess with negative strides?
1778  Flags0 = evaluator<ArgType>::Flags,
1779  LinearAccess =
1780  ((Direction == BothDirections) && (int(Flags0) & PacketAccessBit)) ||
1781  ((ReverseRow && XprType::ColsAtCompileTime == 1) || (ReverseCol && XprType::RowsAtCompileTime == 1))
1782  ? LinearAccessBit
1783  : 0,
1784 
1785  Flags = int(Flags0) & (HereditaryBits | PacketAccessBit | LinearAccess),
1786 
1787  Alignment = 0 // FIXME in some rare cases, Alignment could be preserved, like a Vector4f.
1788  };
1789 
1790  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit unary_evaluator(const XprType& reverse)
1791  : m_argImpl(reverse.nestedExpression()),
1792  m_rows(ReverseRow ? reverse.nestedExpression().rows() : 1),
1793  m_cols(ReverseCol ? reverse.nestedExpression().cols() : 1) {}
1794 
1795  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index row, Index col) const {
1796  return m_argImpl.coeff(ReverseRow ? m_rows.value() - row - 1 : row, ReverseCol ? m_cols.value() - col - 1 : col);
1797  }
1798 
1799  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index index) const {
1800  return m_argImpl.coeff(m_rows.value() * m_cols.value() - index - 1);
1801  }
1802 
1803  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index row, Index col) {
1804  return m_argImpl.coeffRef(ReverseRow ? m_rows.value() - row - 1 : row, ReverseCol ? m_cols.value() - col - 1 : col);
1805  }
1806 
1807  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index index) {
1808  return m_argImpl.coeffRef(m_rows.value() * m_cols.value() - index - 1);
1809  }
1810 
1811  template <int LoadMode, typename PacketType>
1812  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index row, Index col) const {
1813  static constexpr int PacketSize = unpacket_traits<PacketType>::size;
1814  static constexpr int OffsetRow = ReverseRow && IsColMajor ? PacketSize : 1;
1815  static constexpr int OffsetCol = ReverseCol && IsRowMajor ? PacketSize : 1;
1816  using reverse_packet = reverse_packet_cond<PacketType, ReversePacket>;
1817 
1818  Index actualRow = ReverseRow ? m_rows.value() - row - OffsetRow : row;
1819  Index actualCol = ReverseCol ? m_cols.value() - col - OffsetCol : col;
1820 
1821  return reverse_packet::run(m_argImpl.template packet<LoadMode, PacketType>(actualRow, actualCol));
1822  }
1823 
1824  template <int LoadMode, typename PacketType>
1825  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packet(Index index) const {
1826  static constexpr int PacketSize = unpacket_traits<PacketType>::size;
1827 
1828  Index actualIndex = m_rows.value() * m_cols.value() - index - PacketSize;
1829 
1830  return preverse(m_argImpl.template packet<LoadMode, PacketType>(actualIndex));
1831  }
1832 
1833  template <int LoadMode, typename PacketType>
1834  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacket(Index row, Index col, const PacketType& x) {
1835  static constexpr int PacketSize = unpacket_traits<PacketType>::size;
1836  static constexpr int OffsetRow = ReverseRow && IsColMajor ? PacketSize : 1;
1837  static constexpr int OffsetCol = ReverseCol && IsRowMajor ? PacketSize : 1;
1838  using reverse_packet = reverse_packet_cond<PacketType, ReversePacket>;
1839 
1840  Index actualRow = ReverseRow ? m_rows.value() - row - OffsetRow : row;
1841  Index actualCol = ReverseCol ? m_cols.value() - col - OffsetCol : col;
1842 
1843  m_argImpl.template writePacket<LoadMode>(actualRow, actualCol, reverse_packet::run(x));
1844  }
1845 
1846  template <int LoadMode, typename PacketType>
1847  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacket(Index index, const PacketType& x) {
1848  static constexpr int PacketSize = unpacket_traits<PacketType>::size;
1849 
1850  Index actualIndex = m_rows.value() * m_cols.value() - index - PacketSize;
1851 
1852  m_argImpl.template writePacket<LoadMode>(actualIndex, preverse(x));
1853  }
1854 
1855  template <int LoadMode, typename PacketType>
1856  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index row, Index col, Index begin, Index count) const {
1857  static constexpr int PacketSize = unpacket_traits<PacketType>::size;
1858  static constexpr int OffsetRow = ReverseRow && IsColMajor ? PacketSize : 1;
1859  static constexpr int OffsetCol = ReverseCol && IsRowMajor ? PacketSize : 1;
1860  using reverse_packet = reverse_packet_cond<PacketType, ReversePacket>;
1861 
1862  Index actualRow = ReverseRow ? m_rows.value() - row - OffsetRow : row;
1863  Index actualCol = ReverseCol ? m_cols.value() - col - OffsetCol : col;
1864  Index actualBegin = ReversePacket ? (PacketSize - count - begin) : begin;
1865 
1866  return reverse_packet::run(
1867  m_argImpl.template packetSegment<LoadMode, PacketType>(actualRow, actualCol, actualBegin, count));
1868  }
1869 
1870  template <int LoadMode, typename PacketType>
1871  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE PacketType packetSegment(Index index, Index begin, Index count) const {
1872  static constexpr int PacketSize = unpacket_traits<PacketType>::size;
1873 
1874  Index actualIndex = m_rows.value() * m_cols.value() - index - PacketSize;
1875  Index actualBegin = PacketSize - count - begin;
1876 
1877  return preverse(m_argImpl.template packetSegment<LoadMode, PacketType>(actualIndex, actualBegin, count));
1878  }
1879 
1880  template <int LoadMode, typename PacketType>
1881  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacketSegment(Index row, Index col, const PacketType& x, Index begin,
1882  Index count) {
1883  static constexpr int PacketSize = unpacket_traits<PacketType>::size;
1884  static constexpr int OffsetRow = ReverseRow && IsColMajor ? PacketSize : 1;
1885  static constexpr int OffsetCol = ReverseCol && IsRowMajor ? PacketSize : 1;
1886  using reverse_packet = reverse_packet_cond<PacketType, ReversePacket>;
1887 
1888  Index actualRow = ReverseRow ? m_rows.value() - row - OffsetRow : row;
1889  Index actualCol = ReverseCol ? m_cols.value() - col - OffsetCol : col;
1890  Index actualBegin = ReversePacket ? (PacketSize - count - begin) : begin;
1891 
1892  m_argImpl.template writePacketSegment<LoadMode>(actualRow, actualCol, reverse_packet::run(x), actualBegin, count);
1893  }
1894 
1895  template <int LoadMode, typename PacketType>
1896  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void writePacketSegment(Index index, const PacketType& x, Index begin,
1897  Index count) {
1898  static constexpr int PacketSize = unpacket_traits<PacketType>::size;
1899 
1900  Index actualIndex = m_rows.value() * m_cols.value() - index - PacketSize;
1901  Index actualBegin = PacketSize - count - begin;
1902 
1903  m_argImpl.template writePacketSegment<LoadMode>(actualIndex, preverse(x), actualBegin, count);
1904  }
1905 
1906  protected:
1907  evaluator<ArgType> m_argImpl;
1908 
1909  // If we do not reverse rows, then we do not need to know the number of rows; same for columns
1910  // Nonetheless, in this case it is important to set to 1 such that the coeff(index) method works fine for vectors.
1911  const variable_if_dynamic<Index, ReverseRow ? ArgType::RowsAtCompileTime : 1> m_rows;
1912  const variable_if_dynamic<Index, ReverseCol ? ArgType::ColsAtCompileTime : 1> m_cols;
1913 };
1914 
1915 // -------------------- Diagonal --------------------
1916 
1917 template <typename ArgType, int DiagIndex>
1918 struct evaluator<Diagonal<ArgType, DiagIndex>> : evaluator_base<Diagonal<ArgType, DiagIndex>> {
1919  typedef Diagonal<ArgType, DiagIndex> XprType;
1920 
1921  enum {
1922  CoeffReadCost = evaluator<ArgType>::CoeffReadCost,
1923 
1924  Flags =
1925  (unsigned int)(evaluator<ArgType>::Flags & (HereditaryBits | DirectAccessBit) & ~RowMajorBit) | LinearAccessBit,
1926 
1927  Alignment = 0
1928  };
1929 
1930  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE explicit evaluator(const XprType& diagonal)
1931  : m_argImpl(diagonal.nestedExpression()), m_index(diagonal.index()) {}
1932 
1933  typedef typename XprType::Scalar Scalar;
1934  typedef typename XprType::CoeffReturnType CoeffReturnType;
1935 
1936  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index row, Index) const {
1937  return m_argImpl.coeff(row + rowOffset(), row + colOffset());
1938  }
1939 
1940  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE CoeffReturnType coeff(Index index) const {
1941  return m_argImpl.coeff(index + rowOffset(), index + colOffset());
1942  }
1943 
1944  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index row, Index) {
1945  return m_argImpl.coeffRef(row + rowOffset(), row + colOffset());
1946  }
1947 
1948  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar& coeffRef(Index index) {
1949  return m_argImpl.coeffRef(index + rowOffset(), index + colOffset());
1950  }
1951 
1952  protected:
1953  evaluator<ArgType> m_argImpl;
1954  const variable_if_dynamicindex<Index, XprType::DiagIndex> m_index;
1955 
1956  private:
1957  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr Index rowOffset() const {
1958  return m_index.value() > 0 ? 0 : -m_index.value();
1959  }
1960  EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE constexpr Index colOffset() const {
1961  return m_index.value() > 0 ? m_index.value() : 0;
1962  }
1963 };
1964 
1965 //----------------------------------------------------------------------
1966 // deprecated code
1967 //----------------------------------------------------------------------
1968 
1969 // -------------------- EvalToTemp --------------------
1970 
1971 // expression class for evaluating nested expression to a temporary
1972 
1973 template <typename ArgType>
1974 class EvalToTemp;
1975 
1976 template <typename ArgType>
1977 struct traits<EvalToTemp<ArgType>> : public traits<ArgType> {};
1978 
1979 template <typename ArgType>
1980 class EvalToTemp : public dense_xpr_base<EvalToTemp<ArgType>>::type {
1981  public:
1982  typedef typename dense_xpr_base<EvalToTemp>::type Base;
1983  EIGEN_GENERIC_PUBLIC_INTERFACE(EvalToTemp)
1984 
1985  explicit EvalToTemp(const ArgType& arg) : m_arg(arg) {}
1986 
1987  const ArgType& arg() const { return m_arg; }
1988 
1989  constexpr Index rows() const noexcept { return m_arg.rows(); }
1990 
1991  constexpr Index cols() const noexcept { return m_arg.cols(); }
1992 
1993  private:
1994  const ArgType& m_arg;
1995 };
1996 
1997 template <typename ArgType>
1998 struct evaluator<EvalToTemp<ArgType>> : public evaluator<typename ArgType::PlainObject> {
1999  typedef EvalToTemp<ArgType> XprType;
2000  typedef typename ArgType::PlainObject PlainObject;
2001  typedef evaluator<PlainObject> Base;
2002 
2003  EIGEN_DEVICE_FUNC explicit evaluator(const XprType& xpr) : m_result(xpr.arg()) {
2004  internal::construct_at<Base>(this, m_result);
2005  }
2006 
2007  // This constructor is used when nesting an EvalTo evaluator in another evaluator
2008  EIGEN_DEVICE_FUNC evaluator(const ArgType& arg) : m_result(arg) { internal::construct_at<Base>(this, m_result); }
2009 
2010  protected:
2011  PlainObject m_result;
2012 };
2013 
2014 } // namespace internal
2015 
2016 } // end namespace Eigen
2017 
2018 #endif // EIGEN_COREEVALUATORS_H
Definition: Constants.h:266
const unsigned int DirectAccessBit
Definition: Constants.h:159
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_arg_op< typename Derived::Scalar >, const Derived > arg(const Eigen::ArrayBase< Derived > &x)
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
const unsigned int RowMajorBit
Definition: Constants.h:70
const unsigned int PacketAccessBit
Definition: Constants.h:97
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Definition: Constants.h:269
Definition: Constants.h:272
const int Dynamic
Definition: Constants.h:25
const unsigned int EvalBeforeNestingBit
Definition: Constants.h:74
const unsigned int ActualPacketAccessBit
Definition: Constants.h:108
const unsigned int LinearAccessBit
Definition: Constants.h:133