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