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 // 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 #include "../InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
17 namespace internal {
18 
19 template<typename LhsScalar_, typename RhsScalar_> class level3_blocking;
20 
21 /* Specialization for a row-major destination matrix => simple transposition of the product */
22 template<
23  typename Index,
24  typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
25  typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs,
26  int ResInnerStride>
27 struct general_matrix_matrix_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,RowMajor,ResInnerStride>
28 {
29  typedef gebp_traits<RhsScalar,LhsScalar> Traits;
30 
31  typedef typename ScalarBinaryOpTraits<LhsScalar, RhsScalar>::ReturnType ResScalar;
32  static EIGEN_STRONG_INLINE void run(
33  Index rows, Index cols, Index depth,
34  const LhsScalar* lhs, Index lhsStride,
35  const RhsScalar* rhs, Index rhsStride,
36  ResScalar* res, Index resIncr, Index resStride,
37  ResScalar alpha,
38  level3_blocking<RhsScalar,LhsScalar>& blocking,
39  GemmParallelInfo<Index>* info = 0)
40  {
41  // transpose the product such that the result is column major
42  general_matrix_matrix_product<Index,
43  RhsScalar, RhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateRhs,
44  LhsScalar, LhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateLhs,
45  ColMajor,ResInnerStride>
46  ::run(cols,rows,depth,rhs,rhsStride,lhs,lhsStride,res,resIncr,resStride,alpha,blocking,info);
47  }
48 };
49 
50 /* Specialization for a col-major destination matrix
51  * => Blocking algorithm following Goto's paper */
52 template<
53  typename Index,
54  typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
55  typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs,
56  int ResInnerStride>
57 struct general_matrix_matrix_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,ColMajor,ResInnerStride>
58 {
59 
60 typedef gebp_traits<LhsScalar,RhsScalar> Traits;
61 
62 typedef typename ScalarBinaryOpTraits<LhsScalar, RhsScalar>::ReturnType ResScalar;
63 static void run(Index rows, Index cols, Index depth,
64  const LhsScalar* _lhs, Index lhsStride,
65  const RhsScalar* _rhs, Index rhsStride,
66  ResScalar* _res, Index resIncr, Index resStride,
67  ResScalar alpha,
68  level3_blocking<LhsScalar,RhsScalar>& blocking,
69  GemmParallelInfo<Index>* info = 0)
70 {
71  typedef const_blas_data_mapper<LhsScalar, Index, LhsStorageOrder> LhsMapper;
72  typedef const_blas_data_mapper<RhsScalar, Index, RhsStorageOrder> RhsMapper;
73  typedef blas_data_mapper<typename Traits::ResScalar, Index, ColMajor,Unaligned,ResInnerStride> ResMapper;
74  LhsMapper lhs(_lhs, lhsStride);
75  RhsMapper rhs(_rhs, rhsStride);
76  ResMapper res(_res, resStride, resIncr);
77 
78  Index kc = blocking.kc(); // cache block size along the K direction
79  Index mc = (std::min)(rows,blocking.mc()); // cache block size along the M direction
80  Index nc = (std::min)(cols,blocking.nc()); // cache block size along the N direction
81 
82  gemm_pack_lhs<LhsScalar, Index, LhsMapper, Traits::mr, Traits::LhsProgress, typename Traits::LhsPacket4Packing, LhsStorageOrder> pack_lhs;
83  gemm_pack_rhs<RhsScalar, Index, RhsMapper, Traits::nr, RhsStorageOrder> pack_rhs;
84  gebp_kernel<LhsScalar, RhsScalar, Index, ResMapper, Traits::mr, Traits::nr, ConjugateLhs, ConjugateRhs> gebp;
85 
86 #ifdef EIGEN_HAS_OPENMP
87  if(info)
88  {
89  // this is the parallel version!
90  int tid = omp_get_thread_num();
91  int threads = omp_get_num_threads();
92 
93  LhsScalar* blockA = blocking.blockA();
94  eigen_internal_assert(blockA!=0);
95 
96  std::size_t sizeB = kc*nc;
97  ei_declare_aligned_stack_constructed_variable(RhsScalar, blockB, sizeB, 0);
98 
99  // For each horizontal panel of the rhs, and corresponding vertical panel of the lhs...
100  for(Index k=0; k<depth; k+=kc)
101  {
102  const Index actual_kc = (std::min)(k+kc,depth)-k; // => rows of B', and cols of the A'
103 
104  // In order to reduce the chance that a thread has to wait for the other,
105  // let's start by packing B'.
106  pack_rhs(blockB, rhs.getSubMapper(k,0), actual_kc, nc);
107 
108  // Pack A_k to A' in a parallel fashion:
109  // each thread packs the sub block A_k,i to A'_i where i is the thread id.
110 
111  // However, before copying to A'_i, we have to make sure that no other thread is still using it,
112  // i.e., we test that info[tid].users equals 0.
113  // Then, we set info[tid].users to the number of threads to mark that all other threads are going to use it.
114  while(info[tid].users!=0) {}
115  info[tid].users = threads;
116 
117  pack_lhs(blockA+info[tid].lhs_start*actual_kc, lhs.getSubMapper(info[tid].lhs_start,k), actual_kc, info[tid].lhs_length);
118 
119  // Notify the other threads that the part A'_i is ready to go.
120  info[tid].sync = k;
121 
122  // Computes C_i += A' * B' per A'_i
123  for(int shift=0; shift<threads; ++shift)
124  {
125  int i = (tid+shift)%threads;
126 
127  // At this point we have to make sure that A'_i has been updated by the thread i,
128  // we use testAndSetOrdered to mimic a volatile access.
129  // However, no need to wait for the B' part which has been updated by the current thread!
130  if (shift>0) {
131  while(info[i].sync!=k) {
132  }
133  }
134 
135  gebp(res.getSubMapper(info[i].lhs_start, 0), blockA+info[i].lhs_start*actual_kc, blockB, info[i].lhs_length, actual_kc, nc, alpha);
136  }
137 
138  // Then keep going as usual with the remaining B'
139  for(Index j=nc; j<cols; j+=nc)
140  {
141  const Index actual_nc = (std::min)(j+nc,cols)-j;
142 
143  // pack B_k,j to B'
144  pack_rhs(blockB, rhs.getSubMapper(k,j), actual_kc, actual_nc);
145 
146  // C_j += A' * B'
147  gebp(res.getSubMapper(0, j), blockA, blockB, rows, actual_kc, actual_nc, alpha);
148  }
149 
150  // Release all the sub blocks A'_i of A' for the current thread,
151  // i.e., we simply decrement the number of users by 1
152  for(Index i=0; i<threads; ++i)
153  info[i].users -= 1;
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 
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.innerStride(), 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  std::conditional_t<StorageOrder==RowMajor,RhsScalar_,LhsScalar_>,
281  std::conditional_t<StorageOrder==RowMajor,LhsScalar_,RhsScalar_>>
282 {
283  enum {
284  Transpose = StorageOrder==RowMajor,
285  ActualRows = Transpose ? MaxCols : MaxRows,
286  ActualCols = Transpose ? MaxRows : MaxCols
287  };
288  typedef std::conditional_t<Transpose,RhsScalar_,LhsScalar_> LhsScalar;
289  typedef std::conditional_t<Transpose,LhsScalar_,RhsScalar_> RhsScalar;
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*>((std::uintptr_t(m_staticA) + (EIGEN_DEFAULT_ALIGN_BYTES-1)) & ~std::size_t(EIGEN_DEFAULT_ALIGN_BYTES-1));
315  this->m_blockB = reinterpret_cast<RhsScalar*>((std::uintptr_t(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 
327 template<int StorageOrder, typename LhsScalar_, typename RhsScalar_, int MaxRows, int MaxCols, int MaxDepth, int KcFactor>
328 class gemm_blocking_space<StorageOrder,LhsScalar_,RhsScalar_,MaxRows, MaxCols, MaxDepth, KcFactor, false>
329  : public level3_blocking<
330  std::conditional_t<StorageOrder==RowMajor,RhsScalar_,LhsScalar_>,
331  std::conditional_t<StorageOrder==RowMajor,LhsScalar_,RhsScalar_>>
332 {
333  enum {
334  Transpose = StorageOrder==RowMajor
335  };
336  typedef std::conditional_t<Transpose,RhsScalar_,LhsScalar_> LhsScalar;
337  typedef std::conditional_t<Transpose,LhsScalar_,RhsScalar_> RhsScalar;
338 
339  Index m_sizeA;
340  Index m_sizeB;
341 
342  public:
343 
344  gemm_blocking_space(Index rows, Index cols, Index depth, Index num_threads, bool l3_blocking)
345  {
346  this->m_mc = Transpose ? cols : rows;
347  this->m_nc = Transpose ? rows : cols;
348  this->m_kc = depth;
349 
350  if(l3_blocking)
351  {
352  computeProductBlockingSizes<LhsScalar,RhsScalar,KcFactor>(this->m_kc, this->m_mc, this->m_nc, num_threads);
353  }
354  else // no l3 blocking
355  {
356  Index n = this->m_nc;
357  computeProductBlockingSizes<LhsScalar,RhsScalar,KcFactor>(this->m_kc, this->m_mc, n, num_threads);
358  }
359 
360  m_sizeA = this->m_mc * this->m_kc;
361  m_sizeB = this->m_kc * this->m_nc;
362  }
363 
364  void initParallel(Index rows, Index cols, Index depth, Index num_threads)
365  {
366  this->m_mc = Transpose ? cols : rows;
367  this->m_nc = Transpose ? rows : cols;
368  this->m_kc = depth;
369 
370  eigen_internal_assert(this->m_blockA==0 && this->m_blockB==0);
371  Index m = this->m_mc;
372  computeProductBlockingSizes<LhsScalar,RhsScalar,KcFactor>(this->m_kc, m, this->m_nc, num_threads);
373  m_sizeA = this->m_mc * this->m_kc;
374  m_sizeB = this->m_kc * this->m_nc;
375  }
376 
377  void allocateA()
378  {
379  if(this->m_blockA==0)
380  this->m_blockA = aligned_new<LhsScalar>(m_sizeA);
381  }
382 
383  void allocateB()
384  {
385  if(this->m_blockB==0)
386  this->m_blockB = aligned_new<RhsScalar>(m_sizeB);
387  }
388 
389  void allocateAll()
390  {
391  allocateA();
392  allocateB();
393  }
394 
395  ~gemm_blocking_space()
396  {
397  aligned_delete(this->m_blockA, m_sizeA);
398  aligned_delete(this->m_blockB, m_sizeB);
399  }
400 };
401 
402 } // end namespace internal
403 
404 namespace internal {
405 
406 template<typename Lhs, typename Rhs>
407 struct generic_product_impl<Lhs,Rhs,DenseShape,DenseShape,GemmProduct>
408  : generic_product_impl_base<Lhs,Rhs,generic_product_impl<Lhs,Rhs,DenseShape,DenseShape,GemmProduct> >
409 {
410  typedef typename Product<Lhs,Rhs>::Scalar Scalar;
411  typedef typename Lhs::Scalar LhsScalar;
412  typedef typename Rhs::Scalar RhsScalar;
413 
414  typedef internal::blas_traits<Lhs> LhsBlasTraits;
415  typedef typename LhsBlasTraits::DirectLinearAccessType ActualLhsType;
416  typedef internal::remove_all_t<ActualLhsType> ActualLhsTypeCleaned;
417 
418  typedef internal::blas_traits<Rhs> RhsBlasTraits;
419  typedef typename RhsBlasTraits::DirectLinearAccessType ActualRhsType;
420  typedef internal::remove_all_t<ActualRhsType> ActualRhsTypeCleaned;
421 
422  enum {
423  MaxDepthAtCompileTime = min_size_prefer_fixed(Lhs::MaxColsAtCompileTime, Rhs::MaxRowsAtCompileTime)
424  };
425 
426  typedef generic_product_impl<Lhs,Rhs,DenseShape,DenseShape,CoeffBasedProductMode> lazyproduct;
427 
428  template<typename Dst>
429  static void evalTo(Dst& dst, const Lhs& lhs, const Rhs& rhs)
430  {
431  // See http://eigen.tuxfamily.org/bz/show_bug.cgi?id=404 for a discussion and helper program
432  // to determine the following heuristic.
433  // EIGEN_GEMM_TO_COEFFBASED_THRESHOLD is typically defined to 20 in GeneralProduct.h,
434  // unless it has been specialized by the user or for a given architecture.
435  // Note that the condition rhs.rows()>0 was required because lazy product is (was?) not happy with empty inputs.
436  // I'm not sure it is still required.
437  if((rhs.rows()+dst.rows()+dst.cols())<EIGEN_GEMM_TO_COEFFBASED_THRESHOLD && rhs.rows()>0)
438  lazyproduct::eval_dynamic(dst, lhs, rhs, internal::assign_op<typename Dst::Scalar,Scalar>());
439  else
440  {
441  dst.setZero();
442  scaleAndAddTo(dst, lhs, rhs, Scalar(1));
443  }
444  }
445 
446  template<typename Dst>
447  static void addTo(Dst& dst, const Lhs& lhs, const Rhs& rhs)
448  {
449  if((rhs.rows()+dst.rows()+dst.cols())<EIGEN_GEMM_TO_COEFFBASED_THRESHOLD && rhs.rows()>0)
450  lazyproduct::eval_dynamic(dst, lhs, rhs, internal::add_assign_op<typename Dst::Scalar,Scalar>());
451  else
452  scaleAndAddTo(dst,lhs, rhs, Scalar(1));
453  }
454 
455  template<typename Dst>
456  static void subTo(Dst& dst, const Lhs& lhs, const Rhs& rhs)
457  {
458  if((rhs.rows()+dst.rows()+dst.cols())<EIGEN_GEMM_TO_COEFFBASED_THRESHOLD && rhs.rows()>0)
459  lazyproduct::eval_dynamic(dst, lhs, rhs, internal::sub_assign_op<typename Dst::Scalar,Scalar>());
460  else
461  scaleAndAddTo(dst, lhs, rhs, Scalar(-1));
462  }
463 
464  template<typename Dest>
465  static void scaleAndAddTo(Dest& dst, const Lhs& a_lhs, const Rhs& a_rhs, const Scalar& alpha)
466  {
467  eigen_assert(dst.rows()==a_lhs.rows() && dst.cols()==a_rhs.cols());
468  if(a_lhs.cols()==0 || a_lhs.rows()==0 || a_rhs.cols()==0)
469  return;
470 
471  if (dst.cols() == 1)
472  {
473  // Fallback to GEMV if either the lhs or rhs is a runtime vector
474  typename Dest::ColXpr dst_vec(dst.col(0));
475  return internal::generic_product_impl<Lhs,typename Rhs::ConstColXpr,DenseShape,DenseShape,GemvProduct>
476  ::scaleAndAddTo(dst_vec, a_lhs, a_rhs.col(0), alpha);
477  }
478  else if (dst.rows() == 1)
479  {
480  // Fallback to GEMV if either the lhs or rhs is a runtime vector
481  typename Dest::RowXpr dst_vec(dst.row(0));
482  return internal::generic_product_impl<typename Lhs::ConstRowXpr,Rhs,DenseShape,DenseShape,GemvProduct>
483  ::scaleAndAddTo(dst_vec, a_lhs.row(0), a_rhs, alpha);
484  }
485 
486  add_const_on_value_type_t<ActualLhsType> lhs = LhsBlasTraits::extract(a_lhs);
487  add_const_on_value_type_t<ActualRhsType> rhs = RhsBlasTraits::extract(a_rhs);
488 
489  Scalar actualAlpha = combine_scalar_factors(alpha, a_lhs, a_rhs);
490 
491  typedef internal::gemm_blocking_space<(Dest::Flags&RowMajorBit) ? RowMajor : ColMajor,LhsScalar,RhsScalar,
492  Dest::MaxRowsAtCompileTime,Dest::MaxColsAtCompileTime,MaxDepthAtCompileTime> BlockingType;
493 
494  typedef internal::gemm_functor<
495  Scalar, Index,
496  internal::general_matrix_matrix_product<
497  Index,
498  LhsScalar, (ActualLhsTypeCleaned::Flags&RowMajorBit) ? RowMajor : ColMajor, bool(LhsBlasTraits::NeedToConjugate),
499  RhsScalar, (ActualRhsTypeCleaned::Flags&RowMajorBit) ? RowMajor : ColMajor, bool(RhsBlasTraits::NeedToConjugate),
500  (Dest::Flags&RowMajorBit) ? RowMajor : ColMajor,
501  Dest::InnerStrideAtCompileTime>,
502  ActualLhsTypeCleaned, ActualRhsTypeCleaned, Dest, BlockingType> GemmFunctor;
503 
504  BlockingType blocking(dst.rows(), dst.cols(), lhs.cols(), 1, true);
505  internal::parallelize_gemm<(Dest::MaxRowsAtCompileTime>32 || Dest::MaxRowsAtCompileTime==Dynamic)>
506  (GemmFunctor(lhs, rhs, dst, actualAlpha, blocking), a_lhs.rows(), a_rhs.cols(), a_lhs.cols(), Dest::Flags&RowMajorBit);
507  }
508 };
509 
510 } // end namespace internal
511 
512 } // end namespace Eigen
513 
514 #endif // EIGEN_GENERAL_MATRIX_MATRIX_H
Matrix3f m
int n
RowXpr row(Index i)
This is the const version of row(). *‍/.
ColXpr col(Index i)
This is the const version of col().
#define EIGEN_ALIGN_MAX
#define EIGEN_DEFAULT_ALIGN_BYTES
#define EIGEN_GEMM_TO_COEFFBASED_THRESHOLD
IndexedView_or_Block operator()(const RowIndices &rowIndices, const ColIndices &colIndices)
#define eigen_internal_assert(x)
Definition: Macros.h:908
#define EIGEN_UNUSED_VARIABLE(var)
Definition: Macros.h:957
#define eigen_assert(x)
Definition: Macros.h:902
#define ei_declare_aligned_stack_constructed_variable(TYPE, NAME, SIZE, BUFFER)
Definition: Memory.h:847
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
@ ColMajor
Definition: Constants.h:321
@ RowMajor
Definition: Constants.h:323
const unsigned int RowMajorBit
Definition: Constants.h:68
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
void parallelize_gemm(const Functor &func, Index rows, Index cols, Index depth, bool transpose)
Definition: Parallelizer.h:93
EIGEN_ALWAYS_INLINE ResScalar combine_scalar_factors(const ResScalar &alpha, const Lhs &lhs, const Rhs &rhs)
Definition: BlasUtil.h:628
void aligned_delete(T *ptr, std::size_t size)
Definition: Memory.h:447
constexpr int min_size_prefer_fixed(A a, B b)
Definition: Meta.h:553
: InteropHeaders
Definition: Core:139
@ GemmProduct
Definition: Constants.h:504
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
void initParallel()
Definition: Parallelizer.h:51
const int Dynamic
Definition: Constants.h:24
Definition: BFloat16.h:222
std::ptrdiff_t j