Eigen  3.3.9
 
Loading...
Searching...
No Matches
GeneralMatrixMatrix.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008-2009 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_GENERAL_MATRIX_MATRIX_H
11#define EIGEN_GENERAL_MATRIX_MATRIX_H
12
13namespace Eigen {
14
15namespace internal {
16
17template<typename _LhsScalar, typename _RhsScalar> class level3_blocking;
18
19/* Specialization for a row-major destination matrix => simple transposition of the product */
20template<
21 typename Index,
22 typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
23 typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs,
24 int ResInnerStride>
25struct general_matrix_matrix_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,RowMajor,ResInnerStride>
26{
27 typedef gebp_traits<RhsScalar,LhsScalar> Traits;
28
29 typedef typename ScalarBinaryOpTraits<LhsScalar, RhsScalar>::ReturnType ResScalar;
30 static EIGEN_STRONG_INLINE void run(
31 Index rows, Index cols, Index depth,
32 const LhsScalar* lhs, Index lhsStride,
33 const RhsScalar* rhs, Index rhsStride,
34 ResScalar* res, Index resIncr, Index resStride,
35 ResScalar alpha,
36 level3_blocking<RhsScalar,LhsScalar>& blocking,
37 GemmParallelInfo<Index>* info = 0)
38 {
39 // transpose the product such that the result is column major
40 general_matrix_matrix_product<Index,
41 RhsScalar, RhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateRhs,
42 LhsScalar, LhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateLhs,
43 ColMajor,ResInnerStride>
44 ::run(cols,rows,depth,rhs,rhsStride,lhs,lhsStride,res,resIncr,resStride,alpha,blocking,info);
45 }
46};
47
48/* Specialization for a col-major destination matrix
49 * => Blocking algorithm following Goto's paper */
50template<
51 typename Index,
52 typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
53 typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs,
54 int ResInnerStride>
55struct general_matrix_matrix_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,ColMajor,ResInnerStride>
56{
57
58typedef gebp_traits<LhsScalar,RhsScalar> Traits;
59
60typedef typename ScalarBinaryOpTraits<LhsScalar, RhsScalar>::ReturnType ResScalar;
61static void run(Index rows, Index cols, Index depth,
62 const LhsScalar* _lhs, Index lhsStride,
63 const RhsScalar* _rhs, Index rhsStride,
64 ResScalar* _res, Index resIncr, Index resStride,
65 ResScalar alpha,
66 level3_blocking<LhsScalar,RhsScalar>& blocking,
67 GemmParallelInfo<Index>* info = 0)
68{
69 typedef const_blas_data_mapper<LhsScalar, Index, LhsStorageOrder> LhsMapper;
70 typedef const_blas_data_mapper<RhsScalar, Index, RhsStorageOrder> RhsMapper;
71 typedef blas_data_mapper<typename Traits::ResScalar, Index, ColMajor,Unaligned,ResInnerStride> ResMapper;
72 LhsMapper lhs(_lhs, lhsStride);
73 RhsMapper rhs(_rhs, rhsStride);
74 ResMapper res(_res, resStride, resIncr);
75
76 Index kc = blocking.kc(); // cache block size along the K direction
77 Index mc = (std::min)(rows,blocking.mc()); // cache block size along the M direction
78 Index nc = (std::min)(cols,blocking.nc()); // cache block size along the N direction
79
80 gemm_pack_lhs<LhsScalar, Index, LhsMapper, Traits::mr, Traits::LhsProgress, LhsStorageOrder> pack_lhs;
81 gemm_pack_rhs<RhsScalar, Index, RhsMapper, Traits::nr, RhsStorageOrder> pack_rhs;
82 gebp_kernel<LhsScalar, RhsScalar, Index, ResMapper, Traits::mr, Traits::nr, ConjugateLhs, ConjugateRhs> gebp;
83
84#ifdef EIGEN_HAS_OPENMP
85 if(info)
86 {
87 // this is the parallel version!
88 int tid = omp_get_thread_num();
89 int threads = omp_get_num_threads();
90
91 LhsScalar* blockA = blocking.blockA();
92 eigen_internal_assert(blockA!=0);
93
94 std::size_t sizeB = kc*nc;
95 ei_declare_aligned_stack_constructed_variable(RhsScalar, blockB, sizeB, 0);
96
97 // For each horizontal panel of the rhs, and corresponding vertical panel of the lhs...
98 for(Index k=0; k<depth; k+=kc)
99 {
100 const Index actual_kc = (std::min)(k+kc,depth)-k; // => rows of B', and cols of the A'
101
102 // In order to reduce the chance that a thread has to wait for the other,
103 // let's start by packing B'.
104 pack_rhs(blockB, rhs.getSubMapper(k,0), actual_kc, nc);
105
106 // Pack A_k to A' in a parallel fashion:
107 // each thread packs the sub block A_k,i to A'_i where i is the thread id.
108
109 // However, before copying to A'_i, we have to make sure that no other thread is still using it,
110 // i.e., we test that info[tid].users equals 0.
111 // Then, we set info[tid].users to the number of threads to mark that all other threads are going to use it.
112 while(info[tid].users!=0) {}
113 info[tid].users += threads;
114
115 pack_lhs(blockA+info[tid].lhs_start*actual_kc, lhs.getSubMapper(info[tid].lhs_start,k), actual_kc, info[tid].lhs_length);
116
117 // Notify the other threads that the part A'_i is ready to go.
118 info[tid].sync = k;
119
120 // Computes C_i += A' * B' per A'_i
121 for(int shift=0; shift<threads; ++shift)
122 {
123 int i = (tid+shift)%threads;
124
125 // At this point we have to make sure that A'_i has been updated by the thread i,
126 // we use testAndSetOrdered to mimic a volatile access.
127 // However, no need to wait for the B' part which has been updated by the current thread!
128 if (shift>0) {
129 while(info[i].sync!=k) {
130 }
131 }
132
133 gebp(res.getSubMapper(info[i].lhs_start, 0), blockA+info[i].lhs_start*actual_kc, blockB, info[i].lhs_length, actual_kc, nc, alpha);
134 }
135
136 // Then keep going as usual with the remaining B'
137 for(Index j=nc; j<cols; j+=nc)
138 {
139 const Index actual_nc = (std::min)(j+nc,cols)-j;
140
141 // pack B_k,j to B'
142 pack_rhs(blockB, rhs.getSubMapper(k,j), actual_kc, actual_nc);
143
144 // C_j += A' * B'
145 gebp(res.getSubMapper(0, j), blockA, blockB, rows, actual_kc, actual_nc, alpha);
146 }
147
148 // Release all the sub blocks A'_i of A' for the current thread,
149 // i.e., we simply decrement the number of users by 1
150 for(Index i=0; i<threads; ++i)
151 #pragma omp atomic
152 info[i].users -= 1;
153 }
154 }
155 else
156#endif // EIGEN_HAS_OPENMP
157 {
158 EIGEN_UNUSED_VARIABLE(info);
159
160 // this is the sequential version!
161 std::size_t sizeA = kc*mc;
162 std::size_t sizeB = kc*nc;
163
164 ei_declare_aligned_stack_constructed_variable(LhsScalar, blockA, sizeA, blocking.blockA());
165 ei_declare_aligned_stack_constructed_variable(RhsScalar, blockB, sizeB, blocking.blockB());
166
167 const bool pack_rhs_once = mc!=rows && kc==depth && nc==cols;
168
169 // For each horizontal panel of the rhs, and corresponding panel of the lhs...
170 for(Index i2=0; i2<rows; i2+=mc)
171 {
172 const Index actual_mc = (std::min)(i2+mc,rows)-i2;
173
174 for(Index k2=0; k2<depth; k2+=kc)
175 {
176 const Index actual_kc = (std::min)(k2+kc,depth)-k2;
177
178 // OK, here we have selected one horizontal panel of rhs and one vertical panel of lhs.
179 // => Pack lhs's panel into a sequential chunk of memory (L2/L3 caching)
180 // Note that this panel will be read as many times as the number of blocks in the rhs's
181 // horizontal panel which is, in practice, a very low number.
182 pack_lhs(blockA, lhs.getSubMapper(i2,k2), actual_kc, actual_mc);
183
184 // For each kc x nc block of the rhs's horizontal panel...
185 for(Index j2=0; j2<cols; j2+=nc)
186 {
187 const Index actual_nc = (std::min)(j2+nc,cols)-j2;
188
189 // We pack the rhs's block into a sequential chunk of memory (L2 caching)
190 // Note that this block will be read a very high number of times, which is equal to the number of
191 // micro horizontal panel of the large rhs's panel (e.g., rows/12 times).
192 if((!pack_rhs_once) || i2==0)
193 pack_rhs(blockB, rhs.getSubMapper(k2,j2), actual_kc, actual_nc);
194
195 // Everything is packed, we can now call the panel * block kernel:
196 gebp(res.getSubMapper(i2, j2), blockA, blockB, actual_mc, actual_kc, actual_nc, alpha);
197 }
198 }
199 }
200 }
201}
202
203};
204
205/*********************************************************************************
206* Specialization of generic_product_impl for "large" GEMM, i.e.,
207* implementation of the high level wrapper to general_matrix_matrix_product
208**********************************************************************************/
209
210template<typename Scalar, typename Index, typename Gemm, typename Lhs, typename Rhs, typename Dest, typename BlockingType>
211struct gemm_functor
212{
213 gemm_functor(const Lhs& lhs, const Rhs& rhs, Dest& dest, const Scalar& actualAlpha, BlockingType& blocking)
214 : m_lhs(lhs), m_rhs(rhs), m_dest(dest), m_actualAlpha(actualAlpha), m_blocking(blocking)
215 {}
216
217 void initParallelSession(Index num_threads) const
218 {
219 m_blocking.initParallel(m_lhs.rows(), m_rhs.cols(), m_lhs.cols(), num_threads);
220 m_blocking.allocateA();
221 }
222
223 void operator() (Index row, Index rows, Index col=0, Index cols=-1, GemmParallelInfo<Index>* info=0) const
224 {
225 if(cols==-1)
226 cols = m_rhs.cols();
227
228 Gemm::run(rows, cols, m_lhs.cols(),
229 &m_lhs.coeffRef(row,0), m_lhs.outerStride(),
230 &m_rhs.coeffRef(0,col), m_rhs.outerStride(),
231 (Scalar*)&(m_dest.coeffRef(row,col)), m_dest.innerStride(), m_dest.outerStride(),
232 m_actualAlpha, m_blocking, info);
233 }
234
235 typedef typename Gemm::Traits Traits;
236
237 protected:
238 const Lhs& m_lhs;
239 const Rhs& m_rhs;
240 Dest& m_dest;
241 Scalar m_actualAlpha;
242 BlockingType& m_blocking;
243};
244
245template<int StorageOrder, typename LhsScalar, typename RhsScalar, int MaxRows, int MaxCols, int MaxDepth, int KcFactor=1,
246bool FiniteAtCompileTime = MaxRows!=Dynamic && MaxCols!=Dynamic && MaxDepth != Dynamic> class gemm_blocking_space;
247
248template<typename _LhsScalar, typename _RhsScalar>
249class level3_blocking
250{
251 typedef _LhsScalar LhsScalar;
252 typedef _RhsScalar RhsScalar;
253
254 protected:
255 LhsScalar* m_blockA;
256 RhsScalar* m_blockB;
257
258 Index m_mc;
259 Index m_nc;
260 Index m_kc;
261
262 public:
263
264 level3_blocking()
265 : m_blockA(0), m_blockB(0), m_mc(0), m_nc(0), m_kc(0)
266 {}
267
268 inline Index mc() const { return m_mc; }
269 inline Index nc() const { return m_nc; }
270 inline Index kc() const { return m_kc; }
271
272 inline LhsScalar* blockA() { return m_blockA; }
273 inline RhsScalar* blockB() { return m_blockB; }
274};
275
276template<int StorageOrder, typename _LhsScalar, typename _RhsScalar, int MaxRows, int MaxCols, int MaxDepth, int KcFactor>
277class gemm_blocking_space<StorageOrder,_LhsScalar,_RhsScalar,MaxRows, MaxCols, MaxDepth, KcFactor, true /* == FiniteAtCompileTime */>
278 : public level3_blocking<
279 typename conditional<StorageOrder==RowMajor,_RhsScalar,_LhsScalar>::type,
280 typename conditional<StorageOrder==RowMajor,_LhsScalar,_RhsScalar>::type>
281{
282 enum {
283 Transpose = StorageOrder==RowMajor,
284 ActualRows = Transpose ? MaxCols : MaxRows,
285 ActualCols = Transpose ? MaxRows : MaxCols
286 };
287 typedef typename conditional<Transpose,_RhsScalar,_LhsScalar>::type LhsScalar;
288 typedef typename conditional<Transpose,_LhsScalar,_RhsScalar>::type RhsScalar;
289 typedef gebp_traits<LhsScalar,RhsScalar> Traits;
290 enum {
291 SizeA = ActualRows * MaxDepth,
292 SizeB = ActualCols * MaxDepth
293 };
294
295#if EIGEN_MAX_STATIC_ALIGN_BYTES >= EIGEN_DEFAULT_ALIGN_BYTES
296 EIGEN_ALIGN_MAX LhsScalar m_staticA[SizeA];
297 EIGEN_ALIGN_MAX RhsScalar m_staticB[SizeB];
298#else
299 EIGEN_ALIGN_MAX char m_staticA[SizeA * sizeof(LhsScalar) + EIGEN_DEFAULT_ALIGN_BYTES-1];
300 EIGEN_ALIGN_MAX char m_staticB[SizeB * sizeof(RhsScalar) + EIGEN_DEFAULT_ALIGN_BYTES-1];
301#endif
302
303 public:
304
305 gemm_blocking_space(Index /*rows*/, Index /*cols*/, Index /*depth*/, Index /*num_threads*/, bool /*full_rows = false*/)
306 {
307 this->m_mc = ActualRows;
308 this->m_nc = ActualCols;
309 this->m_kc = MaxDepth;
310#if EIGEN_MAX_STATIC_ALIGN_BYTES >= EIGEN_DEFAULT_ALIGN_BYTES
311 this->m_blockA = m_staticA;
312 this->m_blockB = m_staticB;
313#else
314 this->m_blockA = reinterpret_cast<LhsScalar*>((internal::UIntPtr(m_staticA) + (EIGEN_DEFAULT_ALIGN_BYTES-1)) & ~std::size_t(EIGEN_DEFAULT_ALIGN_BYTES-1));
315 this->m_blockB = reinterpret_cast<RhsScalar*>((internal::UIntPtr(m_staticB) + (EIGEN_DEFAULT_ALIGN_BYTES-1)) & ~std::size_t(EIGEN_DEFAULT_ALIGN_BYTES-1));
316#endif
317 }
318
320 {}
321
322 inline void allocateA() {}
323 inline void allocateB() {}
324 inline void allocateAll() {}
325};
326
327template<int StorageOrder, typename _LhsScalar, typename _RhsScalar, int MaxRows, int MaxCols, int MaxDepth, int KcFactor>
328class gemm_blocking_space<StorageOrder,_LhsScalar,_RhsScalar,MaxRows, MaxCols, MaxDepth, KcFactor, false>
329 : public level3_blocking<
330 typename conditional<StorageOrder==RowMajor,_RhsScalar,_LhsScalar>::type,
331 typename conditional<StorageOrder==RowMajor,_LhsScalar,_RhsScalar>::type>
332{
333 enum {
334 Transpose = StorageOrder==RowMajor
335 };
336 typedef typename conditional<Transpose,_RhsScalar,_LhsScalar>::type LhsScalar;
337 typedef typename conditional<Transpose,_LhsScalar,_RhsScalar>::type RhsScalar;
338 typedef gebp_traits<LhsScalar,RhsScalar> Traits;
339
340 Index m_sizeA;
341 Index m_sizeB;
342
343 public:
344
345 gemm_blocking_space(Index rows, Index cols, Index depth, Index num_threads, bool l3_blocking)
346 {
347 this->m_mc = Transpose ? cols : rows;
348 this->m_nc = Transpose ? rows : cols;
349 this->m_kc = depth;
350
351 if(l3_blocking)
352 {
353 computeProductBlockingSizes<LhsScalar,RhsScalar,KcFactor>(this->m_kc, this->m_mc, this->m_nc, num_threads);
354 }
355 else // no l3 blocking
356 {
357 Index n = this->m_nc;
358 computeProductBlockingSizes<LhsScalar,RhsScalar,KcFactor>(this->m_kc, this->m_mc, n, num_threads);
359 }
360
361 m_sizeA = this->m_mc * this->m_kc;
362 m_sizeB = this->m_kc * this->m_nc;
363 }
364
365 void initParallel(Index rows, Index cols, Index depth, Index num_threads)
366 {
367 this->m_mc = Transpose ? cols : rows;
368 this->m_nc = Transpose ? rows : cols;
369 this->m_kc = depth;
370
371 eigen_internal_assert(this->m_blockA==0 && this->m_blockB==0);
372 Index m = this->m_mc;
373 computeProductBlockingSizes<LhsScalar,RhsScalar,KcFactor>(this->m_kc, m, this->m_nc, num_threads);
374 m_sizeA = this->m_mc * this->m_kc;
375 m_sizeB = this->m_kc * this->m_nc;
376 }
377
378 void allocateA()
379 {
380 if(this->m_blockA==0)
381 this->m_blockA = aligned_new<LhsScalar>(m_sizeA);
382 }
383
384 void allocateB()
385 {
386 if(this->m_blockB==0)
387 this->m_blockB = aligned_new<RhsScalar>(m_sizeB);
388 }
389
390 void allocateAll()
391 {
392 allocateA();
393 allocateB();
394 }
395
396 ~gemm_blocking_space()
397 {
398 aligned_delete(this->m_blockA, m_sizeA);
399 aligned_delete(this->m_blockB, m_sizeB);
400 }
401};
402
403} // end namespace internal
404
405namespace internal {
406
407template<typename Lhs, typename Rhs>
408struct generic_product_impl<Lhs,Rhs,DenseShape,DenseShape,GemmProduct>
409 : generic_product_impl_base<Lhs,Rhs,generic_product_impl<Lhs,Rhs,DenseShape,DenseShape,GemmProduct> >
410{
411 typedef typename Product<Lhs,Rhs>::Scalar Scalar;
412 typedef typename Lhs::Scalar LhsScalar;
413 typedef typename Rhs::Scalar RhsScalar;
414
415 typedef internal::blas_traits<Lhs> LhsBlasTraits;
416 typedef typename LhsBlasTraits::DirectLinearAccessType ActualLhsType;
417 typedef typename internal::remove_all<ActualLhsType>::type ActualLhsTypeCleaned;
418
419 typedef internal::blas_traits<Rhs> RhsBlasTraits;
420 typedef typename RhsBlasTraits::DirectLinearAccessType ActualRhsType;
421 typedef typename internal::remove_all<ActualRhsType>::type ActualRhsTypeCleaned;
422
423 enum {
424 MaxDepthAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(Lhs::MaxColsAtCompileTime,Rhs::MaxRowsAtCompileTime)
425 };
426
427 typedef generic_product_impl<Lhs,Rhs,DenseShape,DenseShape,CoeffBasedProductMode> lazyproduct;
428
429 template<typename Dst>
430 static void evalTo(Dst& dst, const Lhs& lhs, const Rhs& rhs)
431 {
432 if((rhs.rows()+dst.rows()+dst.cols())<20 && rhs.rows()>0)
433 lazyproduct::eval_dynamic(dst, lhs, rhs, internal::assign_op<typename Dst::Scalar,Scalar>());
434 else
435 {
436 dst.setZero();
437 scaleAndAddTo(dst, lhs, rhs, Scalar(1));
438 }
439 }
440
441 template<typename Dst>
442 static void addTo(Dst& dst, const Lhs& lhs, const Rhs& rhs)
443 {
444 if((rhs.rows()+dst.rows()+dst.cols())<20 && rhs.rows()>0)
445 lazyproduct::eval_dynamic(dst, lhs, rhs, internal::add_assign_op<typename Dst::Scalar,Scalar>());
446 else
447 scaleAndAddTo(dst,lhs, rhs, Scalar(1));
448 }
449
450 template<typename Dst>
451 static void subTo(Dst& dst, const Lhs& lhs, const Rhs& rhs)
452 {
453 if((rhs.rows()+dst.rows()+dst.cols())<20 && rhs.rows()>0)
454 lazyproduct::eval_dynamic(dst, lhs, rhs, internal::sub_assign_op<typename Dst::Scalar,Scalar>());
455 else
456 scaleAndAddTo(dst, lhs, rhs, Scalar(-1));
457 }
458
459 template<typename Dest>
460 static void scaleAndAddTo(Dest& dst, const Lhs& a_lhs, const Rhs& a_rhs, const Scalar& alpha)
461 {
462 eigen_assert(dst.rows()==a_lhs.rows() && dst.cols()==a_rhs.cols());
463 if(a_lhs.cols()==0 || a_lhs.rows()==0 || a_rhs.cols()==0)
464 return;
465
466 typename internal::add_const_on_value_type<ActualLhsType>::type lhs = LhsBlasTraits::extract(a_lhs);
467 typename internal::add_const_on_value_type<ActualRhsType>::type rhs = RhsBlasTraits::extract(a_rhs);
468
469 Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(a_lhs)
470 * RhsBlasTraits::extractScalarFactor(a_rhs);
471
472 typedef internal::gemm_blocking_space<(Dest::Flags&RowMajorBit) ? RowMajor : ColMajor,LhsScalar,RhsScalar,
473 Dest::MaxRowsAtCompileTime,Dest::MaxColsAtCompileTime,MaxDepthAtCompileTime> BlockingType;
474
475 typedef internal::gemm_functor<
476 Scalar, Index,
477 internal::general_matrix_matrix_product<
478 Index,
479 LhsScalar, (ActualLhsTypeCleaned::Flags&RowMajorBit) ? RowMajor : ColMajor, bool(LhsBlasTraits::NeedToConjugate),
480 RhsScalar, (ActualRhsTypeCleaned::Flags&RowMajorBit) ? RowMajor : ColMajor, bool(RhsBlasTraits::NeedToConjugate),
481 (Dest::Flags&RowMajorBit) ? RowMajor : ColMajor,
482 Dest::InnerStrideAtCompileTime>,
483 ActualLhsTypeCleaned, ActualRhsTypeCleaned, Dest, BlockingType> GemmFunctor;
484
485 BlockingType blocking(dst.rows(), dst.cols(), lhs.cols(), 1, true);
486 internal::parallelize_gemm<(Dest::MaxRowsAtCompileTime>32 || Dest::MaxRowsAtCompileTime==Dynamic)>
487 (GemmFunctor(lhs, rhs, dst, actualAlpha, blocking), a_lhs.rows(), a_rhs.cols(), a_lhs.cols(), Dest::Flags&RowMajorBit);
488 }
489};
490
491} // end namespace internal
492
493} // end namespace Eigen
494
495#endif // EIGEN_GENERAL_MATRIX_MATRIX_H
@ ColMajor
Definition Constants.h:320
@ RowMajor
Definition Constants.h:322
const unsigned int RowMajorBit
Definition Constants.h:61
Namespace containing all symbols from the Eigen library.
Definition A05_PortingFrom2To3.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:65
void initParallel()
Definition Parallelizer.h:49
const int Dynamic
Definition Constants.h:21