GeneralMatrixMatrix.h
Go to the documentation of this file.
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 // Eigen is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 3 of the License, or (at your option) any later version.
10 //
11 // Alternatively, you can redistribute it and/or
12 // modify it under the terms of the GNU General Public License as
13 // published by the Free Software Foundation; either version 2 of
14 // the License, or (at your option) any later version.
15 //
16 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
17 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
19 // GNU General Public License for more details.
20 //
21 // You should have received a copy of the GNU Lesser General Public
22 // License and a copy of the GNU General Public License along with
23 // Eigen. If not, see <http://www.gnu.org/licenses/>.
24 
25 #ifndef EIGEN_GENERAL_MATRIX_MATRIX_H
26 #define EIGEN_GENERAL_MATRIX_MATRIX_H
27 
28 namespace Eigen {
29 
30 namespace internal {
31 
32 template<typename _LhsScalar, typename _RhsScalar> class level3_blocking;
33 
34 /* Specialization for a row-major destination matrix => simple transposition of the product */
35 template<
36  typename Index,
37  typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
38  typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs>
39 struct general_matrix_matrix_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,RowMajor>
40 {
41  typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
42  static EIGEN_STRONG_INLINE void run(
43  Index rows, Index cols, Index depth,
44  const LhsScalar* lhs, Index lhsStride,
45  const RhsScalar* rhs, Index rhsStride,
46  ResScalar* res, Index resStride,
47  ResScalar alpha,
48  level3_blocking<RhsScalar,LhsScalar>& blocking,
49  GemmParallelInfo<Index>* info = 0)
50  {
51  // transpose the product such that the result is column major
52  general_matrix_matrix_product<Index,
53  RhsScalar, RhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateRhs,
54  LhsScalar, LhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateLhs,
55  ColMajor>
56  ::run(cols,rows,depth,rhs,rhsStride,lhs,lhsStride,res,resStride,alpha,blocking,info);
57  }
58 };
59 
60 /* Specialization for a col-major destination matrix
61  * => Blocking algorithm following Goto's paper */
62 template<
63  typename Index,
64  typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
65  typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs>
66 struct general_matrix_matrix_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,ColMajor>
67 {
68 typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
69 static void run(Index rows, Index cols, Index depth,
70  const LhsScalar* _lhs, Index lhsStride,
71  const RhsScalar* _rhs, Index rhsStride,
72  ResScalar* res, Index resStride,
73  ResScalar alpha,
74  level3_blocking<LhsScalar,RhsScalar>& blocking,
75  GemmParallelInfo<Index>* info = 0)
76 {
77  const_blas_data_mapper<LhsScalar, Index, LhsStorageOrder> lhs(_lhs,lhsStride);
78  const_blas_data_mapper<RhsScalar, Index, RhsStorageOrder> rhs(_rhs,rhsStride);
79 
80  typedef gebp_traits<LhsScalar,RhsScalar> Traits;
81 
82  Index kc = blocking.kc(); // cache block size along the K direction
83  Index mc = (std::min)(rows,blocking.mc()); // cache block size along the M direction
84  //Index nc = blocking.nc(); // cache block size along the N direction
85 
86  gemm_pack_lhs<LhsScalar, Index, Traits::mr, Traits::LhsProgress, LhsStorageOrder> pack_lhs;
87  gemm_pack_rhs<RhsScalar, Index, Traits::nr, RhsStorageOrder> pack_rhs;
88  gebp_kernel<LhsScalar, RhsScalar, Index, Traits::mr, Traits::nr, ConjugateLhs, ConjugateRhs> gebp;
89 
90 #ifdef EIGEN_HAS_OPENMP
91  if(info)
92  {
93  // this is the parallel version!
94  Index tid = omp_get_thread_num();
95  Index threads = omp_get_num_threads();
96 
97  std::size_t sizeA = kc*mc;
98  std::size_t sizeW = kc*Traits::WorkSpaceFactor;
99  ei_declare_aligned_stack_constructed_variable(LhsScalar, blockA, sizeA, 0);
100  ei_declare_aligned_stack_constructed_variable(RhsScalar, w, sizeW, 0);
101 
102  RhsScalar* blockB = blocking.blockB();
103  eigen_internal_assert(blockB!=0);
104 
105  // For each horizontal panel of the rhs, and corresponding vertical panel of the lhs...
106  for(Index k=0; k<depth; k+=kc)
107  {
108  const Index actual_kc = (std::min)(k+kc,depth)-k; // => rows of B', and cols of the A'
109 
110  // In order to reduce the chance that a thread has to wait for the other,
111  // let's start by packing A'.
112  pack_lhs(blockA, &lhs(0,k), lhsStride, actual_kc, mc);
113 
114  // Pack B_k to B' in a parallel fashion:
115  // each thread packs the sub block B_k,j to B'_j where j is the thread id.
116 
117  // However, before copying to B'_j, we have to make sure that no other thread is still using it,
118  // i.e., we test that info[tid].users equals 0.
119  // Then, we set info[tid].users to the number of threads to mark that all other threads are going to use it.
120  while(info[tid].users!=0) {}
121  info[tid].users += threads;
122 
123  pack_rhs(blockB+info[tid].rhs_start*actual_kc, &rhs(k,info[tid].rhs_start), rhsStride, actual_kc, info[tid].rhs_length);
124 
125  // Notify the other threads that the part B'_j is ready to go.
126  info[tid].sync = k;
127 
128  // Computes C_i += A' * B' per B'_j
129  for(Index shift=0; shift<threads; ++shift)
130  {
131  Index j = (tid+shift)%threads;
132 
133  // At this point we have to make sure that B'_j has been updated by the thread j,
134  // we use testAndSetOrdered to mimic a volatile access.
135  // However, no need to wait for the B' part which has been updated by the current thread!
136  if(shift>0)
137  while(info[j].sync!=k) {}
138 
139  gebp(res+info[j].rhs_start*resStride, resStride, blockA, blockB+info[j].rhs_start*actual_kc, mc, actual_kc, info[j].rhs_length, alpha, -1,-1,0,0, w);
140  }
141 
142  // Then keep going as usual with the remaining A'
143  for(Index i=mc; i<rows; i+=mc)
144  {
145  const Index actual_mc = (std::min)(i+mc,rows)-i;
146 
147  // pack A_i,k to A'
148  pack_lhs(blockA, &lhs(i,k), lhsStride, actual_kc, actual_mc);
149 
150  // C_i += A' * B'
151  gebp(res+i, resStride, blockA, blockB, actual_mc, actual_kc, cols, alpha, -1,-1,0,0, w);
152  }
153 
154  // Release all the sub blocks B'_j of B' for the current thread,
155  // i.e., we simply decrement the number of users by 1
156  for(Index j=0; j<threads; ++j)
157  #pragma omp atomic
158  --(info[j].users);
159  }
160  }
161  else
162 #endif // EIGEN_HAS_OPENMP
163  {
164  EIGEN_UNUSED_VARIABLE(info);
165 
166  // this is the sequential version!
167  std::size_t sizeA = kc*mc;
168  std::size_t sizeB = kc*cols;
169  std::size_t sizeW = kc*Traits::WorkSpaceFactor;
170 
171  ei_declare_aligned_stack_constructed_variable(LhsScalar, blockA, sizeA, blocking.blockA());
172  ei_declare_aligned_stack_constructed_variable(RhsScalar, blockB, sizeB, blocking.blockB());
173  ei_declare_aligned_stack_constructed_variable(RhsScalar, blockW, sizeW, blocking.blockW());
174 
175  // For each horizontal panel of the rhs, and corresponding panel of the lhs...
176  // (==GEMM_VAR1)
177  for(Index k2=0; k2<depth; k2+=kc)
178  {
179  const Index actual_kc = (std::min)(k2+kc,depth)-k2;
180 
181  // OK, here we have selected one horizontal panel of rhs and one vertical panel of lhs.
182  // => Pack rhs's panel into a sequential chunk of memory (L2 caching)
183  // Note that this panel will be read as many times as the number of blocks in the lhs's
184  // vertical panel which is, in practice, a very low number.
185  pack_rhs(blockB, &rhs(k2,0), rhsStride, actual_kc, cols);
186 
187 
188  // For each mc x kc block of the lhs's vertical panel...
189  // (==GEPP_VAR1)
190  for(Index i2=0; i2<rows; i2+=mc)
191  {
192  const Index actual_mc = (std::min)(i2+mc,rows)-i2;
193 
194  // We pack the lhs's block into a sequential chunk of memory (L1 caching)
195  // Note that this block will be read a very high number of times, which is equal to the number of
196  // micro vertical panel of the large rhs's panel (e.g., cols/4 times).
197  pack_lhs(blockA, &lhs(i2,k2), lhsStride, actual_kc, actual_mc);
198 
199  // Everything is packed, we can now call the block * panel kernel:
200  gebp(res+i2, resStride, blockA, blockB, actual_mc, actual_kc, cols, alpha, -1, -1, 0, 0, blockW);
201 
202  }
203  }
204  }
205 }
206 
207 };
208 
209 /*********************************************************************************
210 * Specialization of GeneralProduct<> for "large" GEMM, i.e.,
211 * implementation of the high level wrapper to general_matrix_matrix_product
212 **********************************************************************************/
213 
214 template<typename Lhs, typename Rhs>
215 struct traits<GeneralProduct<Lhs,Rhs,GemmProduct> >
216  : traits<ProductBase<GeneralProduct<Lhs,Rhs,GemmProduct>, Lhs, Rhs> >
217 {};
218 
219 template<typename Scalar, typename Index, typename Gemm, typename Lhs, typename Rhs, typename Dest, typename BlockingType>
220 struct gemm_functor
221 {
222  gemm_functor(const Lhs& lhs, const Rhs& rhs, Dest& dest, Scalar actualAlpha,
223  BlockingType& blocking)
224  : m_lhs(lhs), m_rhs(rhs), m_dest(dest), m_actualAlpha(actualAlpha), m_blocking(blocking)
225  {}
226 
227  void initParallelSession() const
228  {
229  m_blocking.allocateB();
230  }
231 
232  void operator() (Index row, Index rows, Index col=0, Index cols=-1, GemmParallelInfo<Index>* info=0) const
233  {
234  if(cols==-1)
235  cols = m_rhs.cols();
236 
237  Gemm::run(rows, cols, m_lhs.cols(),
238  /*(const Scalar*)*/&m_lhs.coeffRef(row,0), m_lhs.outerStride(),
239  /*(const Scalar*)*/&m_rhs.coeffRef(0,col), m_rhs.outerStride(),
240  (Scalar*)&(m_dest.coeffRef(row,col)), m_dest.outerStride(),
241  m_actualAlpha, m_blocking, info);
242  }
243 
244  protected:
245  const Lhs& m_lhs;
246  const Rhs& m_rhs;
247  Dest& m_dest;
248  Scalar m_actualAlpha;
249  BlockingType& m_blocking;
250 };
251 
252 template<int StorageOrder, typename LhsScalar, typename RhsScalar, int MaxRows, int MaxCols, int MaxDepth,
253 bool FiniteAtCompileTime = MaxRows!=Dynamic && MaxCols!=Dynamic && MaxDepth != Dynamic> class gemm_blocking_space;
254 
255 template<typename _LhsScalar, typename _RhsScalar>
256 class level3_blocking
257 {
258  typedef _LhsScalar LhsScalar;
259  typedef _RhsScalar RhsScalar;
260 
261  protected:
262  LhsScalar* m_blockA;
263  RhsScalar* m_blockB;
264  RhsScalar* m_blockW;
265 
266  DenseIndex m_mc;
267  DenseIndex m_nc;
268  DenseIndex m_kc;
269 
270  public:
271 
272  level3_blocking()
273  : m_blockA(0), m_blockB(0), m_blockW(0), m_mc(0), m_nc(0), m_kc(0)
274  {}
275 
276  inline DenseIndex mc() const { return m_mc; }
277  inline DenseIndex nc() const { return m_nc; }
278  inline DenseIndex kc() const { return m_kc; }
279 
280  inline LhsScalar* blockA() { return m_blockA; }
281  inline RhsScalar* blockB() { return m_blockB; }
282  inline RhsScalar* blockW() { return m_blockW; }
283 };
284 
285 template<int StorageOrder, typename _LhsScalar, typename _RhsScalar, int MaxRows, int MaxCols, int MaxDepth>
286 class gemm_blocking_space<StorageOrder,_LhsScalar,_RhsScalar,MaxRows, MaxCols, MaxDepth, true>
287  : public level3_blocking<
288  typename conditional<StorageOrder==RowMajor,_RhsScalar,_LhsScalar>::type,
289  typename conditional<StorageOrder==RowMajor,_LhsScalar,_RhsScalar>::type>
290 {
291  enum {
292  Transpose = StorageOrder==RowMajor,
293  ActualRows = Transpose ? MaxCols : MaxRows,
294  ActualCols = Transpose ? MaxRows : MaxCols
295  };
296  typedef typename conditional<Transpose,_RhsScalar,_LhsScalar>::type LhsScalar;
297  typedef typename conditional<Transpose,_LhsScalar,_RhsScalar>::type RhsScalar;
298  typedef gebp_traits<LhsScalar,RhsScalar> Traits;
299  enum {
300  SizeA = ActualRows * MaxDepth,
301  SizeB = ActualCols * MaxDepth,
302  SizeW = MaxDepth * Traits::WorkSpaceFactor
303  };
304 
305  EIGEN_ALIGN16 LhsScalar m_staticA[SizeA];
306  EIGEN_ALIGN16 RhsScalar m_staticB[SizeB];
307  EIGEN_ALIGN16 RhsScalar m_staticW[SizeW];
308 
309  public:
310 
311  gemm_blocking_space(DenseIndex /*rows*/, DenseIndex /*cols*/, DenseIndex /*depth*/)
312  {
313  this->m_mc = ActualRows;
314  this->m_nc = ActualCols;
315  this->m_kc = MaxDepth;
316  this->m_blockA = m_staticA;
317  this->m_blockB = m_staticB;
318  this->m_blockW = m_staticW;
319  }
320 
321  inline void allocateA() {}
322  inline void allocateB() {}
323  inline void allocateW() {}
324  inline void allocateAll() {}
325 };
326 
327 template<int StorageOrder, typename _LhsScalar, typename _RhsScalar, int MaxRows, int MaxCols, int MaxDepth>
328 class gemm_blocking_space<StorageOrder,_LhsScalar,_RhsScalar,MaxRows, MaxCols, MaxDepth, 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  DenseIndex m_sizeA;
341  DenseIndex m_sizeB;
342  DenseIndex m_sizeW;
343 
344  public:
345 
346  gemm_blocking_space(DenseIndex rows, DenseIndex cols, DenseIndex depth)
347  {
348  this->m_mc = Transpose ? cols : rows;
349  this->m_nc = Transpose ? rows : cols;
350  this->m_kc = depth;
351 
352  computeProductBlockingSizes<LhsScalar,RhsScalar>(this->m_kc, this->m_mc, this->m_nc);
353  m_sizeA = this->m_mc * this->m_kc;
354  m_sizeB = this->m_kc * this->m_nc;
355  m_sizeW = this->m_kc*Traits::WorkSpaceFactor;
356  }
357 
358  void allocateA()
359  {
360  if(this->m_blockA==0)
361  this->m_blockA = aligned_new<LhsScalar>(m_sizeA);
362  }
363 
364  void allocateB()
365  {
366  if(this->m_blockB==0)
367  this->m_blockB = aligned_new<RhsScalar>(m_sizeB);
368  }
369 
370  void allocateW()
371  {
372  if(this->m_blockW==0)
373  this->m_blockW = aligned_new<RhsScalar>(m_sizeW);
374  }
375 
376  void allocateAll()
377  {
378  allocateA();
379  allocateB();
380  allocateW();
381  }
382 
383  ~gemm_blocking_space()
384  {
385  aligned_delete(this->m_blockA, m_sizeA);
386  aligned_delete(this->m_blockB, m_sizeB);
387  aligned_delete(this->m_blockW, m_sizeW);
388  }
389 };
390 
391 } // end namespace internal
392 
393 template<typename Lhs, typename Rhs>
394 class GeneralProduct<Lhs, Rhs, GemmProduct>
395  : public ProductBase<GeneralProduct<Lhs,Rhs,GemmProduct>, Lhs, Rhs>
396 {
397  enum {
398  MaxDepthAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(Lhs::MaxColsAtCompileTime,Rhs::MaxRowsAtCompileTime)
399  };
400  public:
402 
403  typedef typename Lhs::Scalar LhsScalar;
404  typedef typename Rhs::Scalar RhsScalar;
405  typedef Scalar ResScalar;
406 
407  GeneralProduct(const Lhs& lhs, const Rhs& rhs) : Base(lhs,rhs)
408  {
409  typedef internal::scalar_product_op<LhsScalar,RhsScalar> BinOp;
410  EIGEN_CHECK_BINARY_COMPATIBILIY(BinOp,LhsScalar,RhsScalar);
411  }
412 
413  template<typename Dest> void scaleAndAddTo(Dest& dst, Scalar alpha) const
414  {
415  eigen_assert(dst.rows()==m_lhs.rows() && dst.cols()==m_rhs.cols());
416 
417  typename internal::add_const_on_value_type<ActualLhsType>::type lhs = LhsBlasTraits::extract(m_lhs);
418  typename internal::add_const_on_value_type<ActualRhsType>::type rhs = RhsBlasTraits::extract(m_rhs);
419 
420  Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(m_lhs)
421  * RhsBlasTraits::extractScalarFactor(m_rhs);
422 
423  typedef internal::gemm_blocking_space<(Dest::Flags&RowMajorBit) ? RowMajor : ColMajor,LhsScalar,RhsScalar,
424  Dest::MaxRowsAtCompileTime,Dest::MaxColsAtCompileTime,MaxDepthAtCompileTime> BlockingType;
425 
426  typedef internal::gemm_functor<
427  Scalar, Index,
428  internal::general_matrix_matrix_product<
429  Index,
430  LhsScalar, (_ActualLhsType::Flags&RowMajorBit) ? RowMajor : ColMajor, bool(LhsBlasTraits::NeedToConjugate),
431  RhsScalar, (_ActualRhsType::Flags&RowMajorBit) ? RowMajor : ColMajor, bool(RhsBlasTraits::NeedToConjugate),
432  (Dest::Flags&RowMajorBit) ? RowMajor : ColMajor>,
433  _ActualLhsType, _ActualRhsType, Dest, BlockingType> GemmFunctor;
434 
435  BlockingType blocking(dst.rows(), dst.cols(), lhs.cols());
436 
437  internal::parallelize_gemm<(Dest::MaxRowsAtCompileTime>32 || Dest::MaxRowsAtCompileTime==Dynamic)>(GemmFunctor(lhs, rhs, dst, actualAlpha, blocking), this->rows(), this->cols(), Dest::Flags&RowMajorBit);
438  }
439 };
440 
441 } // end namespace Eigen
442 
443 #endif // EIGEN_GENERAL_MATRIX_MATRIX_H