Eigen  3.3.9
 
Loading...
Searching...
No Matches
BlasUtil.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_BLASUTIL_H
11#define EIGEN_BLASUTIL_H
12
13// This file contains many lightweight helper classes used to
14// implement and control fast level 2 and level 3 BLAS-like routines.
15
16namespace Eigen {
17
18namespace internal {
19
20// forward declarations
21template<typename LhsScalar, typename RhsScalar, typename Index, typename DataMapper, int mr, int nr, bool ConjugateLhs=false, bool ConjugateRhs=false>
22struct gebp_kernel;
23
24template<typename Scalar, typename Index, typename DataMapper, int nr, int StorageOrder, bool Conjugate = false, bool PanelMode=false>
25struct gemm_pack_rhs;
26
27template<typename Scalar, typename Index, typename DataMapper, int Pack1, int Pack2, int StorageOrder, bool Conjugate = false, bool PanelMode = false>
28struct gemm_pack_lhs;
29
30template<
31 typename Index,
32 typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
33 typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs,
34 int ResStorageOrder, int ResInnerStride>
35struct general_matrix_matrix_product;
36
37template<typename Index,
38 typename LhsScalar, typename LhsMapper, int LhsStorageOrder, bool ConjugateLhs,
39 typename RhsScalar, typename RhsMapper, bool ConjugateRhs, int Version=Specialized>
40struct general_matrix_vector_product;
41
42
43template<bool Conjugate> struct conj_if;
44
45template<> struct conj_if<true> {
46 template<typename T>
47 inline T operator()(const T& x) const { return numext::conj(x); }
48 template<typename T>
49 inline T pconj(const T& x) const { return internal::pconj(x); }
50};
51
52template<> struct conj_if<false> {
53 template<typename T>
54 inline const T& operator()(const T& x) const { return x; }
55 template<typename T>
56 inline const T& pconj(const T& x) const { return x; }
57};
58
59// Generic implementation for custom complex types.
60template<typename LhsScalar, typename RhsScalar, bool ConjLhs, bool ConjRhs>
61struct conj_helper
62{
63 typedef typename ScalarBinaryOpTraits<LhsScalar,RhsScalar>::ReturnType Scalar;
64
65 EIGEN_STRONG_INLINE Scalar pmadd(const LhsScalar& x, const RhsScalar& y, const Scalar& c) const
66 { return padd(c, pmul(x,y)); }
67
68 EIGEN_STRONG_INLINE Scalar pmul(const LhsScalar& x, const RhsScalar& y) const
69 { return conj_if<ConjLhs>()(x) * conj_if<ConjRhs>()(y); }
70};
71
72template<typename Scalar> struct conj_helper<Scalar,Scalar,false,false>
73{
74 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar pmadd(const Scalar& x, const Scalar& y, const Scalar& c) const { return internal::pmadd(x,y,c); }
75 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar pmul(const Scalar& x, const Scalar& y) const { return internal::pmul(x,y); }
76};
77
78template<typename RealScalar> struct conj_helper<std::complex<RealScalar>, std::complex<RealScalar>, false,true>
79{
80 typedef std::complex<RealScalar> Scalar;
81 EIGEN_STRONG_INLINE Scalar pmadd(const Scalar& x, const Scalar& y, const Scalar& c) const
82 { return c + pmul(x,y); }
83
84 EIGEN_STRONG_INLINE Scalar pmul(const Scalar& x, const Scalar& y) const
85 { return Scalar(numext::real(x)*numext::real(y) + numext::imag(x)*numext::imag(y), numext::imag(x)*numext::real(y) - numext::real(x)*numext::imag(y)); }
86};
87
88template<typename RealScalar> struct conj_helper<std::complex<RealScalar>, std::complex<RealScalar>, true,false>
89{
90 typedef std::complex<RealScalar> Scalar;
91 EIGEN_STRONG_INLINE Scalar pmadd(const Scalar& x, const Scalar& y, const Scalar& c) const
92 { return c + pmul(x,y); }
93
94 EIGEN_STRONG_INLINE Scalar pmul(const Scalar& x, const Scalar& y) const
95 { return Scalar(numext::real(x)*numext::real(y) + numext::imag(x)*numext::imag(y), numext::real(x)*numext::imag(y) - numext::imag(x)*numext::real(y)); }
96};
97
98template<typename RealScalar> struct conj_helper<std::complex<RealScalar>, std::complex<RealScalar>, true,true>
99{
100 typedef std::complex<RealScalar> Scalar;
101 EIGEN_STRONG_INLINE Scalar pmadd(const Scalar& x, const Scalar& y, const Scalar& c) const
102 { return c + pmul(x,y); }
103
104 EIGEN_STRONG_INLINE Scalar pmul(const Scalar& x, const Scalar& y) const
105 { return Scalar(numext::real(x)*numext::real(y) - numext::imag(x)*numext::imag(y), - numext::real(x)*numext::imag(y) - numext::imag(x)*numext::real(y)); }
106};
107
108template<typename RealScalar,bool Conj> struct conj_helper<std::complex<RealScalar>, RealScalar, Conj,false>
109{
110 typedef std::complex<RealScalar> Scalar;
111 EIGEN_STRONG_INLINE Scalar pmadd(const Scalar& x, const RealScalar& y, const Scalar& c) const
112 { return padd(c, pmul(x,y)); }
113 EIGEN_STRONG_INLINE Scalar pmul(const Scalar& x, const RealScalar& y) const
114 { return conj_if<Conj>()(x)*y; }
115};
116
117template<typename RealScalar,bool Conj> struct conj_helper<RealScalar, std::complex<RealScalar>, false,Conj>
118{
119 typedef std::complex<RealScalar> Scalar;
120 EIGEN_STRONG_INLINE Scalar pmadd(const RealScalar& x, const Scalar& y, const Scalar& c) const
121 { return padd(c, pmul(x,y)); }
122 EIGEN_STRONG_INLINE Scalar pmul(const RealScalar& x, const Scalar& y) const
123 { return x*conj_if<Conj>()(y); }
124};
125
126template<typename From,typename To> struct get_factor {
127 EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE To run(const From& x) { return To(x); }
128};
129
130template<typename Scalar> struct get_factor<Scalar,typename NumTraits<Scalar>::Real> {
131 EIGEN_DEVICE_FUNC
132 static EIGEN_STRONG_INLINE typename NumTraits<Scalar>::Real run(const Scalar& x) { return numext::real(x); }
133};
134
135
136template<typename Scalar, typename Index>
137class BlasVectorMapper {
138 public:
139 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE BlasVectorMapper(Scalar *data) : m_data(data) {}
140
141 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE Scalar operator()(Index i) const {
142 return m_data[i];
143 }
144 template <typename Packet, int AlignmentType>
145 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE Packet load(Index i) const {
146 return ploadt<Packet, AlignmentType>(m_data + i);
147 }
148
149 template <typename Packet>
150 EIGEN_DEVICE_FUNC bool aligned(Index i) const {
151 return (UIntPtr(m_data+i)%sizeof(Packet))==0;
152 }
153
154 protected:
155 Scalar* m_data;
156};
157
158template<typename Scalar, typename Index, int AlignmentType, int Incr=1>
159class BlasLinearMapper;
160
161template<typename Scalar, typename Index, int AlignmentType>
162class BlasLinearMapper<Scalar,Index,AlignmentType,1> {
163 public:
164 typedef typename packet_traits<Scalar>::type Packet;
165 typedef typename packet_traits<Scalar>::half HalfPacket;
166
167 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE BlasLinearMapper(Scalar *data, Index incr=1)
168 : m_data(data)
169 {
170 EIGEN_ONLY_USED_FOR_DEBUG(incr);
171 eigen_assert(incr==1);
172 }
173
174 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE void prefetch(int i) const {
175 internal::prefetch(&operator()(i));
176 }
177
178 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE Scalar& operator()(Index i) const {
179 return m_data[i];
180 }
181
182 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE Packet loadPacket(Index i) const {
183 return ploadt<Packet, AlignmentType>(m_data + i);
184 }
185
186 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE HalfPacket loadHalfPacket(Index i) const {
187 return ploadt<HalfPacket, AlignmentType>(m_data + i);
188 }
189
190 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE void storePacket(Index i, const Packet &p) const {
191 pstoret<Scalar, Packet, AlignmentType>(m_data + i, p);
192 }
193
194 protected:
195 Scalar *m_data;
196};
197
198// Lightweight helper class to access matrix coefficients.
199template<typename Scalar, typename Index, int StorageOrder, int AlignmentType = Unaligned, int Incr = 1>
200class blas_data_mapper;
201
202template<typename Scalar, typename Index, int StorageOrder, int AlignmentType>
203class blas_data_mapper<Scalar,Index,StorageOrder,AlignmentType,1>
204{
205public:
206 typedef typename packet_traits<Scalar>::type Packet;
207 typedef typename packet_traits<Scalar>::half HalfPacket;
208
209 typedef BlasLinearMapper<Scalar, Index, AlignmentType> LinearMapper;
210 typedef BlasVectorMapper<Scalar, Index> VectorMapper;
211
212 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE blas_data_mapper(Scalar* data, Index stride, Index incr=1)
213 : m_data(data), m_stride(stride)
214 {
215 EIGEN_ONLY_USED_FOR_DEBUG(incr);
216 eigen_assert(incr==1);
217 }
218
219 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE blas_data_mapper<Scalar, Index, StorageOrder, AlignmentType>
220 getSubMapper(Index i, Index j) const {
221 return blas_data_mapper<Scalar, Index, StorageOrder, AlignmentType>(&operator()(i, j), m_stride);
222 }
223
224 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE LinearMapper getLinearMapper(Index i, Index j) const {
225 return LinearMapper(&operator()(i, j));
226 }
227
228 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE VectorMapper getVectorMapper(Index i, Index j) const {
229 return VectorMapper(&operator()(i, j));
230 }
231
232
233 EIGEN_DEVICE_FUNC
234 EIGEN_ALWAYS_INLINE Scalar& operator()(Index i, Index j) const {
235 return m_data[StorageOrder==RowMajor ? j + i*m_stride : i + j*m_stride];
236 }
237
238 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE Packet loadPacket(Index i, Index j) const {
239 return ploadt<Packet, AlignmentType>(&operator()(i, j));
240 }
241
242 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE HalfPacket loadHalfPacket(Index i, Index j) const {
243 return ploadt<HalfPacket, AlignmentType>(&operator()(i, j));
244 }
245
246 template<typename SubPacket>
247 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE void scatterPacket(Index i, Index j, const SubPacket &p) const {
248 pscatter<Scalar, SubPacket>(&operator()(i, j), p, m_stride);
249 }
250
251 template<typename SubPacket>
252 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE SubPacket gatherPacket(Index i, Index j) const {
253 return pgather<Scalar, SubPacket>(&operator()(i, j), m_stride);
254 }
255
256 EIGEN_DEVICE_FUNC const Index stride() const { return m_stride; }
257 EIGEN_DEVICE_FUNC const Scalar* data() const { return m_data; }
258
259 EIGEN_DEVICE_FUNC Index firstAligned(Index size) const {
260 if (UIntPtr(m_data)%sizeof(Scalar)) {
261 return -1;
262 }
263 return internal::first_default_aligned(m_data, size);
264 }
265
266 protected:
267 Scalar* EIGEN_RESTRICT m_data;
268 const Index m_stride;
269};
270
271// Implementation of non-natural increment (i.e. inner-stride != 1)
272// The exposed API is not complete yet compared to the Incr==1 case
273// because some features makes less sense in this case.
274template<typename Scalar, typename Index, int AlignmentType, int Incr>
275class BlasLinearMapper
276{
277public:
278 typedef typename packet_traits<Scalar>::type Packet;
279 typedef typename packet_traits<Scalar>::half HalfPacket;
280
281 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE BlasLinearMapper(Scalar *data,Index incr) : m_data(data), m_incr(incr) {}
282
283 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE void prefetch(int i) const {
284 internal::prefetch(&operator()(i));
285 }
286
287 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE Scalar& operator()(Index i) const {
288 return m_data[i*m_incr.value()];
289 }
290
291 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE Packet loadPacket(Index i) const {
292 return pgather<Scalar,Packet>(m_data + i*m_incr.value(), m_incr.value());
293 }
294
295 template<typename PacketType>
296 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE void storePacket(Index i, const PacketType &p) const {
297 pscatter<Scalar, PacketType>(m_data + i*m_incr.value(), p, m_incr.value());
298 }
299
300protected:
301 Scalar *m_data;
302 const internal::variable_if_dynamic<Index,Incr> m_incr;
303};
304
305template<typename Scalar, typename Index, int StorageOrder, int AlignmentType,int Incr>
306class blas_data_mapper
307{
308public:
309 typedef typename packet_traits<Scalar>::type Packet;
310 typedef typename packet_traits<Scalar>::half HalfPacket;
311
312 typedef BlasLinearMapper<Scalar, Index, AlignmentType,Incr> LinearMapper;
313
314 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE blas_data_mapper(Scalar* data, Index stride, Index incr) : m_data(data), m_stride(stride), m_incr(incr) {}
315
316 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE blas_data_mapper
317 getSubMapper(Index i, Index j) const {
318 return blas_data_mapper(&operator()(i, j), m_stride, m_incr.value());
319 }
320
321 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE LinearMapper getLinearMapper(Index i, Index j) const {
322 return LinearMapper(&operator()(i, j), m_incr.value());
323 }
324
325 EIGEN_DEVICE_FUNC
326 EIGEN_ALWAYS_INLINE Scalar& operator()(Index i, Index j) const {
327 return m_data[StorageOrder==RowMajor ? j*m_incr.value() + i*m_stride : i*m_incr.value() + j*m_stride];
328 }
329
330 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE Packet loadPacket(Index i, Index j) const {
331 return pgather<Scalar,Packet>(&operator()(i, j),m_incr.value());
332 }
333
334 template <typename PacketT, int AlignmentT>
335 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE PacketT load(Index i, Index j) const {
336 return pgather<Scalar,PacketT>(&operator()(i, j),m_incr.value());
337 }
338
339 template<typename SubPacket>
340 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE void scatterPacket(Index i, Index j, const SubPacket &p) const {
341 pscatter<Scalar, SubPacket>(&operator()(i, j), p, m_stride);
342 }
343
344 template<typename SubPacket>
345 EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE SubPacket gatherPacket(Index i, Index j) const {
346 return pgather<Scalar, SubPacket>(&operator()(i, j), m_stride);
347 }
348
349protected:
350 Scalar* EIGEN_RESTRICT m_data;
351 const Index m_stride;
352 const internal::variable_if_dynamic<Index,Incr> m_incr;
353};
354
355// lightweight helper class to access matrix coefficients (const version)
356template<typename Scalar, typename Index, int StorageOrder>
357class const_blas_data_mapper : public blas_data_mapper<const Scalar, Index, StorageOrder> {
358 public:
359 EIGEN_ALWAYS_INLINE const_blas_data_mapper(const Scalar *data, Index stride) : blas_data_mapper<const Scalar, Index, StorageOrder>(data, stride) {}
360
361 EIGEN_ALWAYS_INLINE const_blas_data_mapper<Scalar, Index, StorageOrder> getSubMapper(Index i, Index j) const {
362 return const_blas_data_mapper<Scalar, Index, StorageOrder>(&(this->operator()(i, j)), this->m_stride);
363 }
364};
365
366
367/* Helper class to analyze the factors of a Product expression.
368 * In particular it allows to pop out operator-, scalar multiples,
369 * and conjugate */
370template<typename XprType> struct blas_traits
371{
372 typedef typename traits<XprType>::Scalar Scalar;
373 typedef const XprType& ExtractType;
374 typedef XprType _ExtractType;
375 enum {
376 IsComplex = NumTraits<Scalar>::IsComplex,
377 IsTransposed = false,
378 NeedToConjugate = false,
379 HasUsableDirectAccess = ( (int(XprType::Flags)&DirectAccessBit)
380 && ( bool(XprType::IsVectorAtCompileTime)
381 || int(inner_stride_at_compile_time<XprType>::ret) == 1)
382 ) ? 1 : 0
383 };
384 typedef typename conditional<bool(HasUsableDirectAccess),
385 ExtractType,
386 typename _ExtractType::PlainObject
387 >::type DirectLinearAccessType;
388 static inline ExtractType extract(const XprType& x) { return x; }
389 static inline const Scalar extractScalarFactor(const XprType&) { return Scalar(1); }
390};
391
392// pop conjugate
393template<typename Scalar, typename NestedXpr>
394struct blas_traits<CwiseUnaryOp<scalar_conjugate_op<Scalar>, NestedXpr> >
395 : blas_traits<NestedXpr>
396{
397 typedef blas_traits<NestedXpr> Base;
398 typedef CwiseUnaryOp<scalar_conjugate_op<Scalar>, NestedXpr> XprType;
399 typedef typename Base::ExtractType ExtractType;
400
401 enum {
402 IsComplex = NumTraits<Scalar>::IsComplex,
403 NeedToConjugate = Base::NeedToConjugate ? 0 : IsComplex
404 };
405 static inline ExtractType extract(const XprType& x) { return Base::extract(x.nestedExpression()); }
406 static inline Scalar extractScalarFactor(const XprType& x) { return conj(Base::extractScalarFactor(x.nestedExpression())); }
407};
408
409// pop scalar multiple
410template<typename Scalar, typename NestedXpr, typename Plain>
411struct blas_traits<CwiseBinaryOp<scalar_product_op<Scalar>, const CwiseNullaryOp<scalar_constant_op<Scalar>,Plain>, NestedXpr> >
412 : blas_traits<NestedXpr>
413{
414 typedef blas_traits<NestedXpr> Base;
415 typedef CwiseBinaryOp<scalar_product_op<Scalar>, const CwiseNullaryOp<scalar_constant_op<Scalar>,Plain>, NestedXpr> XprType;
416 typedef typename Base::ExtractType ExtractType;
417 static inline ExtractType extract(const XprType& x) { return Base::extract(x.rhs()); }
418 static inline Scalar extractScalarFactor(const XprType& x)
419 { return x.lhs().functor().m_other * Base::extractScalarFactor(x.rhs()); }
420};
421template<typename Scalar, typename NestedXpr, typename Plain>
422struct blas_traits<CwiseBinaryOp<scalar_product_op<Scalar>, NestedXpr, const CwiseNullaryOp<scalar_constant_op<Scalar>,Plain> > >
423 : blas_traits<NestedXpr>
424{
425 typedef blas_traits<NestedXpr> Base;
426 typedef CwiseBinaryOp<scalar_product_op<Scalar>, NestedXpr, const CwiseNullaryOp<scalar_constant_op<Scalar>,Plain> > XprType;
427 typedef typename Base::ExtractType ExtractType;
428 static inline ExtractType extract(const XprType& x) { return Base::extract(x.lhs()); }
429 static inline Scalar extractScalarFactor(const XprType& x)
430 { return Base::extractScalarFactor(x.lhs()) * x.rhs().functor().m_other; }
431};
432template<typename Scalar, typename Plain1, typename Plain2>
433struct blas_traits<CwiseBinaryOp<scalar_product_op<Scalar>, const CwiseNullaryOp<scalar_constant_op<Scalar>,Plain1>,
434 const CwiseNullaryOp<scalar_constant_op<Scalar>,Plain2> > >
435 : blas_traits<CwiseNullaryOp<scalar_constant_op<Scalar>,Plain1> >
436{};
437
438// pop opposite
439template<typename Scalar, typename NestedXpr>
440struct blas_traits<CwiseUnaryOp<scalar_opposite_op<Scalar>, NestedXpr> >
441 : blas_traits<NestedXpr>
442{
443 typedef blas_traits<NestedXpr> Base;
444 typedef CwiseUnaryOp<scalar_opposite_op<Scalar>, NestedXpr> XprType;
445 typedef typename Base::ExtractType ExtractType;
446 static inline ExtractType extract(const XprType& x) { return Base::extract(x.nestedExpression()); }
447 static inline Scalar extractScalarFactor(const XprType& x)
448 { return - Base::extractScalarFactor(x.nestedExpression()); }
449};
450
451// pop/push transpose
452template<typename NestedXpr>
453struct blas_traits<Transpose<NestedXpr> >
454 : blas_traits<NestedXpr>
455{
456 typedef typename NestedXpr::Scalar Scalar;
457 typedef blas_traits<NestedXpr> Base;
458 typedef Transpose<NestedXpr> XprType;
459 typedef Transpose<const typename Base::_ExtractType> ExtractType; // const to get rid of a compile error; anyway blas traits are only used on the RHS
460 typedef Transpose<const typename Base::_ExtractType> _ExtractType;
461 typedef typename conditional<bool(Base::HasUsableDirectAccess),
462 ExtractType,
463 typename ExtractType::PlainObject
464 >::type DirectLinearAccessType;
465 enum {
466 IsTransposed = Base::IsTransposed ? 0 : 1
467 };
468 static inline ExtractType extract(const XprType& x) { return ExtractType(Base::extract(x.nestedExpression())); }
469 static inline Scalar extractScalarFactor(const XprType& x) { return Base::extractScalarFactor(x.nestedExpression()); }
470};
471
472template<typename T>
473struct blas_traits<const T>
474 : blas_traits<T>
475{};
476
477template<typename T, bool HasUsableDirectAccess=blas_traits<T>::HasUsableDirectAccess>
478struct extract_data_selector {
479 static const typename T::Scalar* run(const T& m)
480 {
481 return blas_traits<T>::extract(m).data();
482 }
483};
484
485template<typename T>
486struct extract_data_selector<T,false> {
487 static typename T::Scalar* run(const T&) { return 0; }
488};
489
490template<typename T> const typename T::Scalar* extract_data(const T& m)
491{
492 return extract_data_selector<T>::run(m);
493}
494
495} // end namespace internal
496
497} // end namespace Eigen
498
499#endif // EIGEN_BLASUTIL_H
AlignmentType
Definition Constants.h:227
@ RowMajor
Definition Constants.h:322
const unsigned int DirectAccessBit
Definition Constants.h:150
Namespace containing all symbols from the Eigen library.
Definition A05_PortingFrom2To3.dox:1
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:65