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