GeneralBlockPanelKernel.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_BLOCK_PANEL_H
26 #define EIGEN_GENERAL_BLOCK_PANEL_H
27 
28 namespace Eigen {
29 
30 namespace internal {
31 
32 template<typename _LhsScalar, typename _RhsScalar, bool _ConjLhs=false, bool _ConjRhs=false>
33 class gebp_traits;
34 
35 
37 inline std::ptrdiff_t manage_caching_sizes_helper(std::ptrdiff_t a, std::ptrdiff_t b)
38 {
39  return a<=0 ? b : a;
40 }
41 
43 inline void manage_caching_sizes(Action action, std::ptrdiff_t* l1=0, std::ptrdiff_t* l2=0)
44 {
45  static std::ptrdiff_t m_l1CacheSize = manage_caching_sizes_helper(queryL1CacheSize(),8 * 1024);
46  static std::ptrdiff_t m_l2CacheSize = manage_caching_sizes_helper(queryTopLevelCacheSize(),1*1024*1024);
47 
48  if(action==SetAction)
49  {
50  // set the cpu cache size and cache all block sizes from a global cache size in byte
51  eigen_internal_assert(l1!=0 && l2!=0);
52  m_l1CacheSize = *l1;
53  m_l2CacheSize = *l2;
54  }
55  else if(action==GetAction)
56  {
57  eigen_internal_assert(l1!=0 && l2!=0);
58  *l1 = m_l1CacheSize;
59  *l2 = m_l2CacheSize;
60  }
61  else
62  {
63  eigen_internal_assert(false);
64  }
65 }
66 
82 template<typename LhsScalar, typename RhsScalar, int KcFactor>
83 void computeProductBlockingSizes(std::ptrdiff_t& k, std::ptrdiff_t& m, std::ptrdiff_t& n)
84 {
86  // Explanations:
87  // Let's recall the product algorithms form kc x nc horizontal panels B' on the rhs and
88  // mc x kc blocks A' on the lhs. A' has to fit into L2 cache. Moreover, B' is processed
89  // per kc x nr vertical small panels where nr is the blocking size along the n dimension
90  // at the register level. For vectorization purpose, these small vertical panels are unpacked,
91  // e.g., each coefficient is replicated to fit a packet. This small vertical panel has to
92  // stay in L1 cache.
93  std::ptrdiff_t l1, l2;
94 
95  typedef gebp_traits<LhsScalar,RhsScalar> Traits;
96  enum {
97  kdiv = KcFactor * 2 * Traits::nr
98  * Traits::RhsProgress * sizeof(RhsScalar),
99  mr = gebp_traits<LhsScalar,RhsScalar>::mr,
100  mr_mask = (0xffffffff/mr)*mr
101  };
102 
103  manage_caching_sizes(GetAction, &l1, &l2);
104  k = std::min<std::ptrdiff_t>(k, l1/kdiv);
105  std::ptrdiff_t _m = k>0 ? l2/(4 * sizeof(LhsScalar) * k) : 0;
106  if(_m<m) m = _m & mr_mask;
107 }
108 
109 template<typename LhsScalar, typename RhsScalar>
110 inline void computeProductBlockingSizes(std::ptrdiff_t& k, std::ptrdiff_t& m, std::ptrdiff_t& n)
111 {
112  computeProductBlockingSizes<LhsScalar,RhsScalar,1>(k, m, n);
113 }
114 
115 #ifdef EIGEN_HAS_FUSE_CJMADD
116  #define MADD(CJ,A,B,C,T) C = CJ.pmadd(A,B,C);
117 #else
118 
119  // FIXME (a bit overkill maybe ?)
120 
121  template<typename CJ, typename A, typename B, typename C, typename T> struct gebp_madd_selector {
122  EIGEN_ALWAYS_INLINE static void run(const CJ& cj, A& a, B& b, C& c, T& /*t*/)
123  {
124  c = cj.pmadd(a,b,c);
125  }
126  };
127 
128  template<typename CJ, typename T> struct gebp_madd_selector<CJ,T,T,T,T> {
129  EIGEN_ALWAYS_INLINE static void run(const CJ& cj, T& a, T& b, T& c, T& t)
130  {
131  t = b; t = cj.pmul(a,t); c = padd(c,t);
132  }
133  };
134 
135  template<typename CJ, typename A, typename B, typename C, typename T>
136  EIGEN_STRONG_INLINE void gebp_madd(const CJ& cj, A& a, B& b, C& c, T& t)
137  {
138  gebp_madd_selector<CJ,A,B,C,T>::run(cj,a,b,c,t);
139  }
140 
141  #define MADD(CJ,A,B,C,T) gebp_madd(CJ,A,B,C,T);
142 // #define MADD(CJ,A,B,C,T) T = B; T = CJ.pmul(A,T); C = padd(C,T);
143 #endif
144 
145 /* Vectorization logic
146  * real*real: unpack rhs to constant packets, ...
147  *
148  * cd*cd : unpack rhs to (b_r,b_r), (b_i,b_i), mul to get (a_r b_r,a_i b_r) (a_r b_i,a_i b_i),
149  * storing each res packet into two packets (2x2),
150  * at the end combine them: swap the second and addsub them
151  * cf*cf : same but with 2x4 blocks
152  * cplx*real : unpack rhs to constant packets, ...
153  * real*cplx : load lhs as (a0,a0,a1,a1), and mul as usual
154  */
155 template<typename _LhsScalar, typename _RhsScalar, bool _ConjLhs, bool _ConjRhs>
156 class gebp_traits
157 {
158 public:
159  typedef _LhsScalar LhsScalar;
160  typedef _RhsScalar RhsScalar;
161  typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
162 
163  enum {
164  ConjLhs = _ConjLhs,
165  ConjRhs = _ConjRhs,
166  Vectorizable = packet_traits<LhsScalar>::Vectorizable && packet_traits<RhsScalar>::Vectorizable,
167  LhsPacketSize = Vectorizable ? packet_traits<LhsScalar>::size : 1,
168  RhsPacketSize = Vectorizable ? packet_traits<RhsScalar>::size : 1,
169  ResPacketSize = Vectorizable ? packet_traits<ResScalar>::size : 1,
170 
171  NumberOfRegisters = EIGEN_ARCH_DEFAULT_NUMBER_OF_REGISTERS,
172 
173  // register block size along the N direction (must be either 2 or 4)
174  nr = NumberOfRegisters/4,
175 
176  // register block size along the M direction (currently, this one cannot be modified)
177  mr = 2 * LhsPacketSize,
178 
179  WorkSpaceFactor = nr * RhsPacketSize,
180 
181  LhsProgress = LhsPacketSize,
182  RhsProgress = RhsPacketSize
183  };
184 
185  typedef typename packet_traits<LhsScalar>::type _LhsPacket;
186  typedef typename packet_traits<RhsScalar>::type _RhsPacket;
187  typedef typename packet_traits<ResScalar>::type _ResPacket;
188 
189  typedef typename conditional<Vectorizable,_LhsPacket,LhsScalar>::type LhsPacket;
190  typedef typename conditional<Vectorizable,_RhsPacket,RhsScalar>::type RhsPacket;
191  typedef typename conditional<Vectorizable,_ResPacket,ResScalar>::type ResPacket;
192 
193  typedef ResPacket AccPacket;
194 
195  EIGEN_STRONG_INLINE void initAcc(AccPacket& p)
196  {
197  p = pset1<ResPacket>(ResScalar(0));
198  }
199 
200  EIGEN_STRONG_INLINE void unpackRhs(DenseIndex n, const RhsScalar* rhs, RhsScalar* b)
201  {
202  for(DenseIndex k=0; k<n; k++)
203  pstore1<RhsPacket>(&b[k*RhsPacketSize], rhs[k]);
204  }
205 
206  EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, RhsPacket& dest) const
207  {
208  dest = pload<RhsPacket>(b);
209  }
210 
211  EIGEN_STRONG_INLINE void loadLhs(const LhsScalar* a, LhsPacket& dest) const
212  {
213  dest = pload<LhsPacket>(a);
214  }
215 
216  EIGEN_STRONG_INLINE void madd(const LhsPacket& a, const RhsPacket& b, AccPacket& c, AccPacket& tmp) const
217  {
218  tmp = b; tmp = pmul(a,tmp); c = padd(c,tmp);
219  }
220 
221  EIGEN_STRONG_INLINE void acc(const AccPacket& c, const ResPacket& alpha, ResPacket& r) const
222  {
223  r = pmadd(c,alpha,r);
224  }
225 
226 protected:
227 // conj_helper<LhsScalar,RhsScalar,ConjLhs,ConjRhs> cj;
228 // conj_helper<LhsPacket,RhsPacket,ConjLhs,ConjRhs> pcj;
229 };
230 
231 template<typename RealScalar, bool _ConjLhs>
232 class gebp_traits<std::complex<RealScalar>, RealScalar, _ConjLhs, false>
233 {
234 public:
235  typedef std::complex<RealScalar> LhsScalar;
236  typedef RealScalar RhsScalar;
237  typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
238 
239  enum {
240  ConjLhs = _ConjLhs,
241  ConjRhs = false,
242  Vectorizable = packet_traits<LhsScalar>::Vectorizable && packet_traits<RhsScalar>::Vectorizable,
243  LhsPacketSize = Vectorizable ? packet_traits<LhsScalar>::size : 1,
244  RhsPacketSize = Vectorizable ? packet_traits<RhsScalar>::size : 1,
245  ResPacketSize = Vectorizable ? packet_traits<ResScalar>::size : 1,
246 
247  NumberOfRegisters = EIGEN_ARCH_DEFAULT_NUMBER_OF_REGISTERS,
248  nr = NumberOfRegisters/4,
249  mr = 2 * LhsPacketSize,
250  WorkSpaceFactor = nr*RhsPacketSize,
251 
252  LhsProgress = LhsPacketSize,
253  RhsProgress = RhsPacketSize
254  };
255 
256  typedef typename packet_traits<LhsScalar>::type _LhsPacket;
257  typedef typename packet_traits<RhsScalar>::type _RhsPacket;
258  typedef typename packet_traits<ResScalar>::type _ResPacket;
259 
260  typedef typename conditional<Vectorizable,_LhsPacket,LhsScalar>::type LhsPacket;
261  typedef typename conditional<Vectorizable,_RhsPacket,RhsScalar>::type RhsPacket;
262  typedef typename conditional<Vectorizable,_ResPacket,ResScalar>::type ResPacket;
263 
264  typedef ResPacket AccPacket;
265 
266  EIGEN_STRONG_INLINE void initAcc(AccPacket& p)
267  {
268  p = pset1<ResPacket>(ResScalar(0));
269  }
270 
271  EIGEN_STRONG_INLINE void unpackRhs(DenseIndex n, const RhsScalar* rhs, RhsScalar* b)
272  {
273  for(DenseIndex k=0; k<n; k++)
274  pstore1<RhsPacket>(&b[k*RhsPacketSize], rhs[k]);
275  }
276 
277  EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, RhsPacket& dest) const
278  {
279  dest = pload<RhsPacket>(b);
280  }
281 
282  EIGEN_STRONG_INLINE void loadLhs(const LhsScalar* a, LhsPacket& dest) const
283  {
284  dest = pload<LhsPacket>(a);
285  }
286 
287  EIGEN_STRONG_INLINE void madd(const LhsPacket& a, const RhsPacket& b, AccPacket& c, RhsPacket& tmp) const
288  {
289  madd_impl(a, b, c, tmp, typename conditional<Vectorizable,true_type,false_type>::type());
290  }
291 
292  EIGEN_STRONG_INLINE void madd_impl(const LhsPacket& a, const RhsPacket& b, AccPacket& c, RhsPacket& tmp, const true_type&) const
293  {
294  tmp = b; tmp = pmul(a.v,tmp); c.v = padd(c.v,tmp);
295  }
296 
297  EIGEN_STRONG_INLINE void madd_impl(const LhsScalar& a, const RhsScalar& b, ResScalar& c, RhsScalar& /*tmp*/, const false_type&) const
298  {
299  c += a * b;
300  }
301 
302  EIGEN_STRONG_INLINE void acc(const AccPacket& c, const ResPacket& alpha, ResPacket& r) const
303  {
304  r = cj.pmadd(c,alpha,r);
305  }
306 
307 protected:
308  conj_helper<ResPacket,ResPacket,ConjLhs,false> cj;
309 };
310 
311 template<typename RealScalar, bool _ConjLhs, bool _ConjRhs>
312 class gebp_traits<std::complex<RealScalar>, std::complex<RealScalar>, _ConjLhs, _ConjRhs >
313 {
314 public:
315  typedef std::complex<RealScalar> Scalar;
316  typedef std::complex<RealScalar> LhsScalar;
317  typedef std::complex<RealScalar> RhsScalar;
318  typedef std::complex<RealScalar> ResScalar;
319 
320  enum {
321  ConjLhs = _ConjLhs,
322  ConjRhs = _ConjRhs,
323  Vectorizable = packet_traits<RealScalar>::Vectorizable
324  && packet_traits<Scalar>::Vectorizable,
325  RealPacketSize = Vectorizable ? packet_traits<RealScalar>::size : 1,
326  ResPacketSize = Vectorizable ? packet_traits<ResScalar>::size : 1,
327 
328  nr = 2,
329  mr = 2 * ResPacketSize,
330  WorkSpaceFactor = Vectorizable ? 2*nr*RealPacketSize : nr,
331 
332  LhsProgress = ResPacketSize,
333  RhsProgress = Vectorizable ? 2*ResPacketSize : 1
334  };
335 
336  typedef typename packet_traits<RealScalar>::type RealPacket;
337  typedef typename packet_traits<Scalar>::type ScalarPacket;
338  struct DoublePacket
339  {
340  RealPacket first;
341  RealPacket second;
342  };
343 
344  typedef typename conditional<Vectorizable,RealPacket, Scalar>::type LhsPacket;
345  typedef typename conditional<Vectorizable,DoublePacket,Scalar>::type RhsPacket;
346  typedef typename conditional<Vectorizable,ScalarPacket,Scalar>::type ResPacket;
347  typedef typename conditional<Vectorizable,DoublePacket,Scalar>::type AccPacket;
348 
349  EIGEN_STRONG_INLINE void initAcc(Scalar& p) { p = Scalar(0); }
350 
351  EIGEN_STRONG_INLINE void initAcc(DoublePacket& p)
352  {
353  p.first = pset1<RealPacket>(RealScalar(0));
354  p.second = pset1<RealPacket>(RealScalar(0));
355  }
356 
357  /* Unpack the rhs coeff such that each complex coefficient is spread into
358  * two packects containing respectively the real and imaginary coefficient
359  * duplicated as many time as needed: (x+iy) => [x, ..., x] [y, ..., y]
360  */
361  EIGEN_STRONG_INLINE void unpackRhs(DenseIndex n, const Scalar* rhs, Scalar* b)
362  {
363  for(DenseIndex k=0; k<n; k++)
364  {
365  if(Vectorizable)
366  {
367  pstore1<RealPacket>((RealScalar*)&b[k*ResPacketSize*2+0], real(rhs[k]));
368  pstore1<RealPacket>((RealScalar*)&b[k*ResPacketSize*2+ResPacketSize], imag(rhs[k]));
369  }
370  else
371  b[k] = rhs[k];
372  }
373  }
374 
375  EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, ResPacket& dest) const { dest = *b; }
376 
377  EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, DoublePacket& dest) const
378  {
379  dest.first = pload<RealPacket>((const RealScalar*)b);
380  dest.second = pload<RealPacket>((const RealScalar*)(b+ResPacketSize));
381  }
382 
383  // nothing special here
384  EIGEN_STRONG_INLINE void loadLhs(const LhsScalar* a, LhsPacket& dest) const
385  {
386  dest = pload<LhsPacket>((const typename unpacket_traits<LhsPacket>::type*)(a));
387  }
388 
389  EIGEN_STRONG_INLINE void madd(const LhsPacket& a, const RhsPacket& b, DoublePacket& c, RhsPacket& /*tmp*/) const
390  {
391  c.first = padd(pmul(a,b.first), c.first);
392  c.second = padd(pmul(a,b.second),c.second);
393  }
394 
395  EIGEN_STRONG_INLINE void madd(const LhsPacket& a, const RhsPacket& b, ResPacket& c, RhsPacket& /*tmp*/) const
396  {
397  c = cj.pmadd(a,b,c);
398  }
399 
400  EIGEN_STRONG_INLINE void acc(const Scalar& c, const Scalar& alpha, Scalar& r) const { r += alpha * c; }
401 
402  EIGEN_STRONG_INLINE void acc(const DoublePacket& c, const ResPacket& alpha, ResPacket& r) const
403  {
404  // assemble c
405  ResPacket tmp;
406  if((!ConjLhs)&&(!ConjRhs))
407  {
408  tmp = pcplxflip(pconj(ResPacket(c.second)));
409  tmp = padd(ResPacket(c.first),tmp);
410  }
411  else if((!ConjLhs)&&(ConjRhs))
412  {
413  tmp = pconj(pcplxflip(ResPacket(c.second)));
414  tmp = padd(ResPacket(c.first),tmp);
415  }
416  else if((ConjLhs)&&(!ConjRhs))
417  {
418  tmp = pcplxflip(ResPacket(c.second));
419  tmp = padd(pconj(ResPacket(c.first)),tmp);
420  }
421  else if((ConjLhs)&&(ConjRhs))
422  {
423  tmp = pcplxflip(ResPacket(c.second));
424  tmp = psub(pconj(ResPacket(c.first)),tmp);
425  }
426 
427  r = pmadd(tmp,alpha,r);
428  }
429 
430 protected:
431  conj_helper<LhsScalar,RhsScalar,ConjLhs,ConjRhs> cj;
432 };
433 
434 template<typename RealScalar, bool _ConjRhs>
435 class gebp_traits<RealScalar, std::complex<RealScalar>, false, _ConjRhs >
436 {
437 public:
438  typedef std::complex<RealScalar> Scalar;
439  typedef RealScalar LhsScalar;
440  typedef Scalar RhsScalar;
441  typedef Scalar ResScalar;
442 
443  enum {
444  ConjLhs = false,
445  ConjRhs = _ConjRhs,
446  Vectorizable = packet_traits<RealScalar>::Vectorizable
447  && packet_traits<Scalar>::Vectorizable,
448  LhsPacketSize = Vectorizable ? packet_traits<LhsScalar>::size : 1,
449  RhsPacketSize = Vectorizable ? packet_traits<RhsScalar>::size : 1,
450  ResPacketSize = Vectorizable ? packet_traits<ResScalar>::size : 1,
451 
452  NumberOfRegisters = EIGEN_ARCH_DEFAULT_NUMBER_OF_REGISTERS,
453  nr = 4,
454  mr = 2*ResPacketSize,
455  WorkSpaceFactor = nr*RhsPacketSize,
456 
457  LhsProgress = ResPacketSize,
458  RhsProgress = ResPacketSize
459  };
460 
461  typedef typename packet_traits<LhsScalar>::type _LhsPacket;
462  typedef typename packet_traits<RhsScalar>::type _RhsPacket;
463  typedef typename packet_traits<ResScalar>::type _ResPacket;
464 
465  typedef typename conditional<Vectorizable,_LhsPacket,LhsScalar>::type LhsPacket;
466  typedef typename conditional<Vectorizable,_RhsPacket,RhsScalar>::type RhsPacket;
467  typedef typename conditional<Vectorizable,_ResPacket,ResScalar>::type ResPacket;
468 
469  typedef ResPacket AccPacket;
470 
471  EIGEN_STRONG_INLINE void initAcc(AccPacket& p)
472  {
473  p = pset1<ResPacket>(ResScalar(0));
474  }
475 
476  EIGEN_STRONG_INLINE void unpackRhs(DenseIndex n, const RhsScalar* rhs, RhsScalar* b)
477  {
478  for(DenseIndex k=0; k<n; k++)
479  pstore1<RhsPacket>(&b[k*RhsPacketSize], rhs[k]);
480  }
481 
482  EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, RhsPacket& dest) const
483  {
484  dest = pload<RhsPacket>(b);
485  }
486 
487  EIGEN_STRONG_INLINE void loadLhs(const LhsScalar* a, LhsPacket& dest) const
488  {
489  dest = ploaddup<LhsPacket>(a);
490  }
491 
492  EIGEN_STRONG_INLINE void madd(const LhsPacket& a, const RhsPacket& b, AccPacket& c, RhsPacket& tmp) const
493  {
494  madd_impl(a, b, c, tmp, typename conditional<Vectorizable,true_type,false_type>::type());
495  }
496 
497  EIGEN_STRONG_INLINE void madd_impl(const LhsPacket& a, const RhsPacket& b, AccPacket& c, RhsPacket& tmp, const true_type&) const
498  {
499  tmp = b; tmp.v = pmul(a,tmp.v); c = padd(c,tmp);
500  }
501 
502  EIGEN_STRONG_INLINE void madd_impl(const LhsScalar& a, const RhsScalar& b, ResScalar& c, RhsScalar& /*tmp*/, const false_type&) const
503  {
504  c += a * b;
505  }
506 
507  EIGEN_STRONG_INLINE void acc(const AccPacket& c, const ResPacket& alpha, ResPacket& r) const
508  {
509  r = cj.pmadd(alpha,c,r);
510  }
511 
512 protected:
513  conj_helper<ResPacket,ResPacket,false,ConjRhs> cj;
514 };
515 
516 /* optimized GEneral packed Block * packed Panel product kernel
517  *
518  * Mixing type logic: C += A * B
519  * | A | B | comments
520  * |real |cplx | no vectorization yet, would require to pack A with duplication
521  * |cplx |real | easy vectorization
522  */
523 template<typename LhsScalar, typename RhsScalar, typename Index, int mr, int nr, bool ConjugateLhs, bool ConjugateRhs>
524 struct gebp_kernel
525 {
526  typedef gebp_traits<LhsScalar,RhsScalar,ConjugateLhs,ConjugateRhs> Traits;
527  typedef typename Traits::ResScalar ResScalar;
528  typedef typename Traits::LhsPacket LhsPacket;
529  typedef typename Traits::RhsPacket RhsPacket;
530  typedef typename Traits::ResPacket ResPacket;
531  typedef typename Traits::AccPacket AccPacket;
532 
533  enum {
534  Vectorizable = Traits::Vectorizable,
535  LhsProgress = Traits::LhsProgress,
536  RhsProgress = Traits::RhsProgress,
537  ResPacketSize = Traits::ResPacketSize
538  };
539 
541  void operator()(ResScalar* res, Index resStride, const LhsScalar* blockA, const RhsScalar* blockB, Index rows, Index depth, Index cols, ResScalar alpha,
542  Index strideA=-1, Index strideB=-1, Index offsetA=0, Index offsetB=0, RhsScalar* unpackedB = 0)
543  {
544  Traits traits;
545 
546  if(strideA==-1) strideA = depth;
547  if(strideB==-1) strideB = depth;
548  conj_helper<LhsScalar,RhsScalar,ConjugateLhs,ConjugateRhs> cj;
549 // conj_helper<LhsPacket,RhsPacket,ConjugateLhs,ConjugateRhs> pcj;
550  Index packet_cols = (cols/nr) * nr;
551  const Index peeled_mc = (rows/mr)*mr;
552  // FIXME:
553  const Index peeled_mc2 = peeled_mc + (rows-peeled_mc >= LhsProgress ? LhsProgress : 0);
554  const Index peeled_kc = (depth/4)*4;
555 
556  if(unpackedB==0)
557  unpackedB = const_cast<RhsScalar*>(blockB - strideB * nr * RhsProgress);
558 
559  // loops on each micro vertical panel of rhs (depth x nr)
560  for(Index j2=0; j2<packet_cols; j2+=nr)
561  {
562  traits.unpackRhs(depth*nr,&blockB[j2*strideB+offsetB*nr],unpackedB);
563 
564  // loops on each largest micro horizontal panel of lhs (mr x depth)
565  // => we select a mr x nr micro block of res which is entirely
566  // stored into mr/packet_size x nr registers.
567  for(Index i=0; i<peeled_mc; i+=mr)
568  {
569  const LhsScalar* blA = &blockA[i*strideA+offsetA*mr];
570  prefetch(&blA[0]);
571 
572  // gets res block as register
573  AccPacket C0, C1, C2, C3, C4, C5, C6, C7;
574  traits.initAcc(C0);
575  traits.initAcc(C1);
576  if(nr==4) traits.initAcc(C2);
577  if(nr==4) traits.initAcc(C3);
578  traits.initAcc(C4);
579  traits.initAcc(C5);
580  if(nr==4) traits.initAcc(C6);
581  if(nr==4) traits.initAcc(C7);
582 
583  ResScalar* r0 = &res[(j2+0)*resStride + i];
584  ResScalar* r1 = r0 + resStride;
585  ResScalar* r2 = r1 + resStride;
586  ResScalar* r3 = r2 + resStride;
587 
588  prefetch(r0+16);
589  prefetch(r1+16);
590  prefetch(r2+16);
591  prefetch(r3+16);
592 
593  // performs "inner" product
594  // TODO let's check wether the folowing peeled loop could not be
595  // optimized via optimal prefetching from one loop to the other
596  const RhsScalar* blB = unpackedB;
597  for(Index k=0; k<peeled_kc; k+=4)
598  {
599  if(nr==2)
600  {
601  LhsPacket A0, A1;
602  RhsPacket B_0;
603  RhsPacket T0;
604 
605 EIGEN_ASM_COMMENT("mybegin2");
606  traits.loadLhs(&blA[0*LhsProgress], A0);
607  traits.loadLhs(&blA[1*LhsProgress], A1);
608  traits.loadRhs(&blB[0*RhsProgress], B_0);
609  traits.madd(A0,B_0,C0,T0);
610  traits.madd(A1,B_0,C4,B_0);
611  traits.loadRhs(&blB[1*RhsProgress], B_0);
612  traits.madd(A0,B_0,C1,T0);
613  traits.madd(A1,B_0,C5,B_0);
614 
615  traits.loadLhs(&blA[2*LhsProgress], A0);
616  traits.loadLhs(&blA[3*LhsProgress], A1);
617  traits.loadRhs(&blB[2*RhsProgress], B_0);
618  traits.madd(A0,B_0,C0,T0);
619  traits.madd(A1,B_0,C4,B_0);
620  traits.loadRhs(&blB[3*RhsProgress], B_0);
621  traits.madd(A0,B_0,C1,T0);
622  traits.madd(A1,B_0,C5,B_0);
623 
624  traits.loadLhs(&blA[4*LhsProgress], A0);
625  traits.loadLhs(&blA[5*LhsProgress], A1);
626  traits.loadRhs(&blB[4*RhsProgress], B_0);
627  traits.madd(A0,B_0,C0,T0);
628  traits.madd(A1,B_0,C4,B_0);
629  traits.loadRhs(&blB[5*RhsProgress], B_0);
630  traits.madd(A0,B_0,C1,T0);
631  traits.madd(A1,B_0,C5,B_0);
632 
633  traits.loadLhs(&blA[6*LhsProgress], A0);
634  traits.loadLhs(&blA[7*LhsProgress], A1);
635  traits.loadRhs(&blB[6*RhsProgress], B_0);
636  traits.madd(A0,B_0,C0,T0);
637  traits.madd(A1,B_0,C4,B_0);
638  traits.loadRhs(&blB[7*RhsProgress], B_0);
639  traits.madd(A0,B_0,C1,T0);
640  traits.madd(A1,B_0,C5,B_0);
641 EIGEN_ASM_COMMENT("myend");
642  }
643  else
644  {
645 EIGEN_ASM_COMMENT("mybegin4");
646  LhsPacket A0, A1;
647  RhsPacket B_0, B1, B2, B3;
648  RhsPacket T0;
649 
650  traits.loadLhs(&blA[0*LhsProgress], A0);
651  traits.loadLhs(&blA[1*LhsProgress], A1);
652  traits.loadRhs(&blB[0*RhsProgress], B_0);
653  traits.loadRhs(&blB[1*RhsProgress], B1);
654 
655  traits.madd(A0,B_0,C0,T0);
656  traits.loadRhs(&blB[2*RhsProgress], B2);
657  traits.madd(A1,B_0,C4,B_0);
658  traits.loadRhs(&blB[3*RhsProgress], B3);
659  traits.loadRhs(&blB[4*RhsProgress], B_0);
660  traits.madd(A0,B1,C1,T0);
661  traits.madd(A1,B1,C5,B1);
662  traits.loadRhs(&blB[5*RhsProgress], B1);
663  traits.madd(A0,B2,C2,T0);
664  traits.madd(A1,B2,C6,B2);
665  traits.loadRhs(&blB[6*RhsProgress], B2);
666  traits.madd(A0,B3,C3,T0);
667  traits.loadLhs(&blA[2*LhsProgress], A0);
668  traits.madd(A1,B3,C7,B3);
669  traits.loadLhs(&blA[3*LhsProgress], A1);
670  traits.loadRhs(&blB[7*RhsProgress], B3);
671  traits.madd(A0,B_0,C0,T0);
672  traits.madd(A1,B_0,C4,B_0);
673  traits.loadRhs(&blB[8*RhsProgress], B_0);
674  traits.madd(A0,B1,C1,T0);
675  traits.madd(A1,B1,C5,B1);
676  traits.loadRhs(&blB[9*RhsProgress], B1);
677  traits.madd(A0,B2,C2,T0);
678  traits.madd(A1,B2,C6,B2);
679  traits.loadRhs(&blB[10*RhsProgress], B2);
680  traits.madd(A0,B3,C3,T0);
681  traits.loadLhs(&blA[4*LhsProgress], A0);
682  traits.madd(A1,B3,C7,B3);
683  traits.loadLhs(&blA[5*LhsProgress], A1);
684  traits.loadRhs(&blB[11*RhsProgress], B3);
685 
686  traits.madd(A0,B_0,C0,T0);
687  traits.madd(A1,B_0,C4,B_0);
688  traits.loadRhs(&blB[12*RhsProgress], B_0);
689  traits.madd(A0,B1,C1,T0);
690  traits.madd(A1,B1,C5,B1);
691  traits.loadRhs(&blB[13*RhsProgress], B1);
692  traits.madd(A0,B2,C2,T0);
693  traits.madd(A1,B2,C6,B2);
694  traits.loadRhs(&blB[14*RhsProgress], B2);
695  traits.madd(A0,B3,C3,T0);
696  traits.loadLhs(&blA[6*LhsProgress], A0);
697  traits.madd(A1,B3,C7,B3);
698  traits.loadLhs(&blA[7*LhsProgress], A1);
699  traits.loadRhs(&blB[15*RhsProgress], B3);
700  traits.madd(A0,B_0,C0,T0);
701  traits.madd(A1,B_0,C4,B_0);
702  traits.madd(A0,B1,C1,T0);
703  traits.madd(A1,B1,C5,B1);
704  traits.madd(A0,B2,C2,T0);
705  traits.madd(A1,B2,C6,B2);
706  traits.madd(A0,B3,C3,T0);
707  traits.madd(A1,B3,C7,B3);
708  }
709 
710  blB += 4*nr*RhsProgress;
711  blA += 4*mr;
712  }
713  // process remaining peeled loop
714  for(Index k=peeled_kc; k<depth; k++)
715  {
716  if(nr==2)
717  {
718  LhsPacket A0, A1;
719  RhsPacket B_0;
720  RhsPacket T0;
721 
722  traits.loadLhs(&blA[0*LhsProgress], A0);
723  traits.loadLhs(&blA[1*LhsProgress], A1);
724  traits.loadRhs(&blB[0*RhsProgress], B_0);
725  traits.madd(A0,B_0,C0,T0);
726  traits.madd(A1,B_0,C4,B_0);
727  traits.loadRhs(&blB[1*RhsProgress], B_0);
728  traits.madd(A0,B_0,C1,T0);
729  traits.madd(A1,B_0,C5,B_0);
730  }
731  else
732  {
733  LhsPacket A0, A1;
734  RhsPacket B_0, B1, B2, B3;
735  RhsPacket T0;
736 
737  traits.loadLhs(&blA[0*LhsProgress], A0);
738  traits.loadLhs(&blA[1*LhsProgress], A1);
739  traits.loadRhs(&blB[0*RhsProgress], B_0);
740  traits.loadRhs(&blB[1*RhsProgress], B1);
741 
742  traits.madd(A0,B_0,C0,T0);
743  traits.loadRhs(&blB[2*RhsProgress], B2);
744  traits.madd(A1,B_0,C4,B_0);
745  traits.loadRhs(&blB[3*RhsProgress], B3);
746  traits.madd(A0,B1,C1,T0);
747  traits.madd(A1,B1,C5,B1);
748  traits.madd(A0,B2,C2,T0);
749  traits.madd(A1,B2,C6,B2);
750  traits.madd(A0,B3,C3,T0);
751  traits.madd(A1,B3,C7,B3);
752  }
753 
754  blB += nr*RhsProgress;
755  blA += mr;
756  }
757 
758  if(nr==4)
759  {
760  ResPacket R0, R1, R2, R3, R4, R5, R6;
761  ResPacket alphav = pset1<ResPacket>(alpha);
762 
763  R0 = ploadu<ResPacket>(r0);
764  R1 = ploadu<ResPacket>(r1);
765  R2 = ploadu<ResPacket>(r2);
766  R3 = ploadu<ResPacket>(r3);
767  R4 = ploadu<ResPacket>(r0 + ResPacketSize);
768  R5 = ploadu<ResPacket>(r1 + ResPacketSize);
769  R6 = ploadu<ResPacket>(r2 + ResPacketSize);
770  traits.acc(C0, alphav, R0);
771  pstoreu(r0, R0);
772  R0 = ploadu<ResPacket>(r3 + ResPacketSize);
773 
774  traits.acc(C1, alphav, R1);
775  traits.acc(C2, alphav, R2);
776  traits.acc(C3, alphav, R3);
777  traits.acc(C4, alphav, R4);
778  traits.acc(C5, alphav, R5);
779  traits.acc(C6, alphav, R6);
780  traits.acc(C7, alphav, R0);
781 
782  pstoreu(r1, R1);
783  pstoreu(r2, R2);
784  pstoreu(r3, R3);
785  pstoreu(r0 + ResPacketSize, R4);
786  pstoreu(r1 + ResPacketSize, R5);
787  pstoreu(r2 + ResPacketSize, R6);
788  pstoreu(r3 + ResPacketSize, R0);
789  }
790  else
791  {
792  ResPacket R0, R1, R4;
793  ResPacket alphav = pset1<ResPacket>(alpha);
794 
795  R0 = ploadu<ResPacket>(r0);
796  R1 = ploadu<ResPacket>(r1);
797  R4 = ploadu<ResPacket>(r0 + ResPacketSize);
798  traits.acc(C0, alphav, R0);
799  pstoreu(r0, R0);
800  R0 = ploadu<ResPacket>(r1 + ResPacketSize);
801  traits.acc(C1, alphav, R1);
802  traits.acc(C4, alphav, R4);
803  traits.acc(C5, alphav, R0);
804  pstoreu(r1, R1);
805  pstoreu(r0 + ResPacketSize, R4);
806  pstoreu(r1 + ResPacketSize, R0);
807  }
808 
809  }
810 
811  if(rows-peeled_mc>=LhsProgress)
812  {
813  Index i = peeled_mc;
814  const LhsScalar* blA = &blockA[i*strideA+offsetA*LhsProgress];
815  prefetch(&blA[0]);
816 
817  // gets res block as register
818  AccPacket C0, C1, C2, C3;
819  traits.initAcc(C0);
820  traits.initAcc(C1);
821  if(nr==4) traits.initAcc(C2);
822  if(nr==4) traits.initAcc(C3);
823 
824  // performs "inner" product
825  const RhsScalar* blB = unpackedB;
826  for(Index k=0; k<peeled_kc; k+=4)
827  {
828  if(nr==2)
829  {
830  LhsPacket A0;
831  RhsPacket B_0, B1;
832 
833  traits.loadLhs(&blA[0*LhsProgress], A0);
834  traits.loadRhs(&blB[0*RhsProgress], B_0);
835  traits.loadRhs(&blB[1*RhsProgress], B1);
836  traits.madd(A0,B_0,C0,B_0);
837  traits.loadRhs(&blB[2*RhsProgress], B_0);
838  traits.madd(A0,B1,C1,B1);
839  traits.loadLhs(&blA[1*LhsProgress], A0);
840  traits.loadRhs(&blB[3*RhsProgress], B1);
841  traits.madd(A0,B_0,C0,B_0);
842  traits.loadRhs(&blB[4*RhsProgress], B_0);
843  traits.madd(A0,B1,C1,B1);
844  traits.loadLhs(&blA[2*LhsProgress], A0);
845  traits.loadRhs(&blB[5*RhsProgress], B1);
846  traits.madd(A0,B_0,C0,B_0);
847  traits.loadRhs(&blB[6*RhsProgress], B_0);
848  traits.madd(A0,B1,C1,B1);
849  traits.loadLhs(&blA[3*LhsProgress], A0);
850  traits.loadRhs(&blB[7*RhsProgress], B1);
851  traits.madd(A0,B_0,C0,B_0);
852  traits.madd(A0,B1,C1,B1);
853  }
854  else
855  {
856  LhsPacket A0;
857  RhsPacket B_0, B1, B2, B3;
858 
859  traits.loadLhs(&blA[0*LhsProgress], A0);
860  traits.loadRhs(&blB[0*RhsProgress], B_0);
861  traits.loadRhs(&blB[1*RhsProgress], B1);
862 
863  traits.madd(A0,B_0,C0,B_0);
864  traits.loadRhs(&blB[2*RhsProgress], B2);
865  traits.loadRhs(&blB[3*RhsProgress], B3);
866  traits.loadRhs(&blB[4*RhsProgress], B_0);
867  traits.madd(A0,B1,C1,B1);
868  traits.loadRhs(&blB[5*RhsProgress], B1);
869  traits.madd(A0,B2,C2,B2);
870  traits.loadRhs(&blB[6*RhsProgress], B2);
871  traits.madd(A0,B3,C3,B3);
872  traits.loadLhs(&blA[1*LhsProgress], A0);
873  traits.loadRhs(&blB[7*RhsProgress], B3);
874  traits.madd(A0,B_0,C0,B_0);
875  traits.loadRhs(&blB[8*RhsProgress], B_0);
876  traits.madd(A0,B1,C1,B1);
877  traits.loadRhs(&blB[9*RhsProgress], B1);
878  traits.madd(A0,B2,C2,B2);
879  traits.loadRhs(&blB[10*RhsProgress], B2);
880  traits.madd(A0,B3,C3,B3);
881  traits.loadLhs(&blA[2*LhsProgress], A0);
882  traits.loadRhs(&blB[11*RhsProgress], B3);
883 
884  traits.madd(A0,B_0,C0,B_0);
885  traits.loadRhs(&blB[12*RhsProgress], B_0);
886  traits.madd(A0,B1,C1,B1);
887  traits.loadRhs(&blB[13*RhsProgress], B1);
888  traits.madd(A0,B2,C2,B2);
889  traits.loadRhs(&blB[14*RhsProgress], B2);
890  traits.madd(A0,B3,C3,B3);
891 
892  traits.loadLhs(&blA[3*LhsProgress], A0);
893  traits.loadRhs(&blB[15*RhsProgress], B3);
894  traits.madd(A0,B_0,C0,B_0);
895  traits.madd(A0,B1,C1,B1);
896  traits.madd(A0,B2,C2,B2);
897  traits.madd(A0,B3,C3,B3);
898  }
899 
900  blB += nr*4*RhsProgress;
901  blA += 4*LhsProgress;
902  }
903  // process remaining peeled loop
904  for(Index k=peeled_kc; k<depth; k++)
905  {
906  if(nr==2)
907  {
908  LhsPacket A0;
909  RhsPacket B_0, B1;
910 
911  traits.loadLhs(&blA[0*LhsProgress], A0);
912  traits.loadRhs(&blB[0*RhsProgress], B_0);
913  traits.loadRhs(&blB[1*RhsProgress], B1);
914  traits.madd(A0,B_0,C0,B_0);
915  traits.madd(A0,B1,C1,B1);
916  }
917  else
918  {
919  LhsPacket A0;
920  RhsPacket B_0, B1, B2, B3;
921 
922  traits.loadLhs(&blA[0*LhsProgress], A0);
923  traits.loadRhs(&blB[0*RhsProgress], B_0);
924  traits.loadRhs(&blB[1*RhsProgress], B1);
925  traits.loadRhs(&blB[2*RhsProgress], B2);
926  traits.loadRhs(&blB[3*RhsProgress], B3);
927 
928  traits.madd(A0,B_0,C0,B_0);
929  traits.madd(A0,B1,C1,B1);
930  traits.madd(A0,B2,C2,B2);
931  traits.madd(A0,B3,C3,B3);
932  }
933 
934  blB += nr*RhsProgress;
935  blA += LhsProgress;
936  }
937 
938  ResPacket R0, R1, R2, R3;
939  ResPacket alphav = pset1<ResPacket>(alpha);
940 
941  ResScalar* r0 = &res[(j2+0)*resStride + i];
942  ResScalar* r1 = r0 + resStride;
943  ResScalar* r2 = r1 + resStride;
944  ResScalar* r3 = r2 + resStride;
945 
946  R0 = ploadu<ResPacket>(r0);
947  R1 = ploadu<ResPacket>(r1);
948  if(nr==4) R2 = ploadu<ResPacket>(r2);
949  if(nr==4) R3 = ploadu<ResPacket>(r3);
950 
951  traits.acc(C0, alphav, R0);
952  traits.acc(C1, alphav, R1);
953  if(nr==4) traits.acc(C2, alphav, R2);
954  if(nr==4) traits.acc(C3, alphav, R3);
955 
956  pstoreu(r0, R0);
957  pstoreu(r1, R1);
958  if(nr==4) pstoreu(r2, R2);
959  if(nr==4) pstoreu(r3, R3);
960  }
961  for(Index i=peeled_mc2; i<rows; i++)
962  {
963  const LhsScalar* blA = &blockA[i*strideA+offsetA];
964  prefetch(&blA[0]);
965 
966  // gets a 1 x nr res block as registers
967  ResScalar C0(0), C1(0), C2(0), C3(0);
968  // TODO directly use blockB ???
969  const RhsScalar* blB = &blockB[j2*strideB+offsetB*nr];
970  for(Index k=0; k<depth; k++)
971  {
972  if(nr==2)
973  {
974  LhsScalar A0;
975  RhsScalar B_0, B1;
976 
977  A0 = blA[k];
978  B_0 = blB[0];
979  B1 = blB[1];
980  MADD(cj,A0,B_0,C0,B_0);
981  MADD(cj,A0,B1,C1,B1);
982  }
983  else
984  {
985  LhsScalar A0;
986  RhsScalar B_0, B1, B2, B3;
987 
988  A0 = blA[k];
989  B_0 = blB[0];
990  B1 = blB[1];
991  B2 = blB[2];
992  B3 = blB[3];
993 
994  MADD(cj,A0,B_0,C0,B_0);
995  MADD(cj,A0,B1,C1,B1);
996  MADD(cj,A0,B2,C2,B2);
997  MADD(cj,A0,B3,C3,B3);
998  }
999 
1000  blB += nr;
1001  }
1002  res[(j2+0)*resStride + i] += alpha*C0;
1003  res[(j2+1)*resStride + i] += alpha*C1;
1004  if(nr==4) res[(j2+2)*resStride + i] += alpha*C2;
1005  if(nr==4) res[(j2+3)*resStride + i] += alpha*C3;
1006  }
1007  }
1008  // process remaining rhs/res columns one at a time
1009  // => do the same but with nr==1
1010  for(Index j2=packet_cols; j2<cols; j2++)
1011  {
1012  // unpack B
1013  traits.unpackRhs(depth, &blockB[j2*strideB+offsetB], unpackedB);
1014 
1015  for(Index i=0; i<peeled_mc; i+=mr)
1016  {
1017  const LhsScalar* blA = &blockA[i*strideA+offsetA*mr];
1018  prefetch(&blA[0]);
1019 
1020  // TODO move the res loads to the stores
1021 
1022  // get res block as registers
1023  AccPacket C0, C4;
1024  traits.initAcc(C0);
1025  traits.initAcc(C4);
1026 
1027  const RhsScalar* blB = unpackedB;
1028  for(Index k=0; k<depth; k++)
1029  {
1030  LhsPacket A0, A1;
1031  RhsPacket B_0;
1032  RhsPacket T0;
1033 
1034  traits.loadLhs(&blA[0*LhsProgress], A0);
1035  traits.loadLhs(&blA[1*LhsProgress], A1);
1036  traits.loadRhs(&blB[0*RhsProgress], B_0);
1037  traits.madd(A0,B_0,C0,T0);
1038  traits.madd(A1,B_0,C4,B_0);
1039 
1040  blB += RhsProgress;
1041  blA += 2*LhsProgress;
1042  }
1043  ResPacket R0, R4;
1044  ResPacket alphav = pset1<ResPacket>(alpha);
1045 
1046  ResScalar* r0 = &res[(j2+0)*resStride + i];
1047 
1048  R0 = ploadu<ResPacket>(r0);
1049  R4 = ploadu<ResPacket>(r0+ResPacketSize);
1050 
1051  traits.acc(C0, alphav, R0);
1052  traits.acc(C4, alphav, R4);
1053 
1054  pstoreu(r0, R0);
1055  pstoreu(r0+ResPacketSize, R4);
1056  }
1057  if(rows-peeled_mc>=LhsProgress)
1058  {
1059  Index i = peeled_mc;
1060  const LhsScalar* blA = &blockA[i*strideA+offsetA*LhsProgress];
1061  prefetch(&blA[0]);
1062 
1063  AccPacket C0;
1064  traits.initAcc(C0);
1065 
1066  const RhsScalar* blB = unpackedB;
1067  for(Index k=0; k<depth; k++)
1068  {
1069  LhsPacket A0;
1070  RhsPacket B_0;
1071  traits.loadLhs(blA, A0);
1072  traits.loadRhs(blB, B_0);
1073  traits.madd(A0, B_0, C0, B_0);
1074  blB += RhsProgress;
1075  blA += LhsProgress;
1076  }
1077 
1078  ResPacket alphav = pset1<ResPacket>(alpha);
1079  ResPacket R0 = ploadu<ResPacket>(&res[(j2+0)*resStride + i]);
1080  traits.acc(C0, alphav, R0);
1081  pstoreu(&res[(j2+0)*resStride + i], R0);
1082  }
1083  for(Index i=peeled_mc2; i<rows; i++)
1084  {
1085  const LhsScalar* blA = &blockA[i*strideA+offsetA];
1086  prefetch(&blA[0]);
1087 
1088  // gets a 1 x 1 res block as registers
1089  ResScalar C0(0);
1090  // FIXME directly use blockB ??
1091  const RhsScalar* blB = &blockB[j2*strideB+offsetB];
1092  for(Index k=0; k<depth; k++)
1093  {
1094  LhsScalar A0 = blA[k];
1095  RhsScalar B_0 = blB[k];
1096  MADD(cj, A0, B_0, C0, B_0);
1097  }
1098  res[(j2+0)*resStride + i] += alpha*C0;
1099  }
1100  }
1101  }
1102 };
1103 
1104 #undef CJMADD
1105 
1106 // pack a block of the lhs
1107 // The traversal is as follow (mr==4):
1108 // 0 4 8 12 ...
1109 // 1 5 9 13 ...
1110 // 2 6 10 14 ...
1111 // 3 7 11 15 ...
1112 //
1113 // 16 20 24 28 ...
1114 // 17 21 25 29 ...
1115 // 18 22 26 30 ...
1116 // 19 23 27 31 ...
1117 //
1118 // 32 33 34 35 ...
1119 // 36 36 38 39 ...
1120 template<typename Scalar, typename Index, int Pack1, int Pack2, int StorageOrder, bool Conjugate, bool PanelMode>
1121 struct gemm_pack_lhs
1122 {
1123  EIGEN_DONT_INLINE void operator()(Scalar* blockA, const Scalar* EIGEN_RESTRICT _lhs, Index lhsStride, Index depth, Index rows,
1124  Index stride=0, Index offset=0)
1125  {
1126  typedef typename packet_traits<Scalar>::type Packet;
1127  enum { PacketSize = packet_traits<Scalar>::size };
1128 
1129  EIGEN_ASM_COMMENT("EIGEN PRODUCT PACK LHS");
1130  eigen_assert(((!PanelMode) && stride==0 && offset==0) || (PanelMode && stride>=depth && offset<=stride));
1131  eigen_assert( (StorageOrder==RowMajor) || ((Pack1%PacketSize)==0 && Pack1<=4*PacketSize) );
1132  conj_if<NumTraits<Scalar>::IsComplex && Conjugate> cj;
1133  const_blas_data_mapper<Scalar, Index, StorageOrder> lhs(_lhs,lhsStride);
1134  Index count = 0;
1135  Index peeled_mc = (rows/Pack1)*Pack1;
1136  for(Index i=0; i<peeled_mc; i+=Pack1)
1137  {
1138  if(PanelMode) count += Pack1 * offset;
1139 
1140  if(StorageOrder==ColMajor)
1141  {
1142  for(Index k=0; k<depth; k++)
1143  {
1144  Packet A, B, C, D;
1145  if(Pack1>=1*PacketSize) A = ploadu<Packet>(&lhs(i+0*PacketSize, k));
1146  if(Pack1>=2*PacketSize) B = ploadu<Packet>(&lhs(i+1*PacketSize, k));
1147  if(Pack1>=3*PacketSize) C = ploadu<Packet>(&lhs(i+2*PacketSize, k));
1148  if(Pack1>=4*PacketSize) D = ploadu<Packet>(&lhs(i+3*PacketSize, k));
1149  if(Pack1>=1*PacketSize) { pstore(blockA+count, cj.pconj(A)); count+=PacketSize; }
1150  if(Pack1>=2*PacketSize) { pstore(blockA+count, cj.pconj(B)); count+=PacketSize; }
1151  if(Pack1>=3*PacketSize) { pstore(blockA+count, cj.pconj(C)); count+=PacketSize; }
1152  if(Pack1>=4*PacketSize) { pstore(blockA+count, cj.pconj(D)); count+=PacketSize; }
1153  }
1154  }
1155  else
1156  {
1157  for(Index k=0; k<depth; k++)
1158  {
1159  // TODO add a vectorized transpose here
1160  Index w=0;
1161  for(; w<Pack1-3; w+=4)
1162  {
1163  Scalar a(cj(lhs(i+w+0, k))),
1164  b(cj(lhs(i+w+1, k))),
1165  c(cj(lhs(i+w+2, k))),
1166  d(cj(lhs(i+w+3, k)));
1167  blockA[count++] = a;
1168  blockA[count++] = b;
1169  blockA[count++] = c;
1170  blockA[count++] = d;
1171  }
1172  if(Pack1%4)
1173  for(;w<Pack1;++w)
1174  blockA[count++] = cj(lhs(i+w, k));
1175  }
1176  }
1177  if(PanelMode) count += Pack1 * (stride-offset-depth);
1178  }
1179  if(rows-peeled_mc>=Pack2)
1180  {
1181  if(PanelMode) count += Pack2*offset;
1182  for(Index k=0; k<depth; k++)
1183  for(Index w=0; w<Pack2; w++)
1184  blockA[count++] = cj(lhs(peeled_mc+w, k));
1185  if(PanelMode) count += Pack2 * (stride-offset-depth);
1186  peeled_mc += Pack2;
1187  }
1188  for(Index i=peeled_mc; i<rows; i++)
1189  {
1190  if(PanelMode) count += offset;
1191  for(Index k=0; k<depth; k++)
1192  blockA[count++] = cj(lhs(i, k));
1193  if(PanelMode) count += (stride-offset-depth);
1194  }
1195  }
1196 };
1197 
1198 // copy a complete panel of the rhs
1199 // this version is optimized for column major matrices
1200 // The traversal order is as follow: (nr==4):
1201 // 0 1 2 3 12 13 14 15 24 27
1202 // 4 5 6 7 16 17 18 19 25 28
1203 // 8 9 10 11 20 21 22 23 26 29
1204 // . . . . . . . . . .
1205 template<typename Scalar, typename Index, int nr, bool Conjugate, bool PanelMode>
1206 struct gemm_pack_rhs<Scalar, Index, nr, ColMajor, Conjugate, PanelMode>
1207 {
1208  typedef typename packet_traits<Scalar>::type Packet;
1209  enum { PacketSize = packet_traits<Scalar>::size };
1210  EIGEN_DONT_INLINE void operator()(Scalar* blockB, const Scalar* rhs, Index rhsStride, Index depth, Index cols,
1211  Index stride=0, Index offset=0)
1212  {
1213  EIGEN_ASM_COMMENT("EIGEN PRODUCT PACK RHS COLMAJOR");
1214  eigen_assert(((!PanelMode) && stride==0 && offset==0) || (PanelMode && stride>=depth && offset<=stride));
1215  conj_if<NumTraits<Scalar>::IsComplex && Conjugate> cj;
1216  Index packet_cols = (cols/nr) * nr;
1217  Index count = 0;
1218  for(Index j2=0; j2<packet_cols; j2+=nr)
1219  {
1220  // skip what we have before
1221  if(PanelMode) count += nr * offset;
1222  const Scalar* b0 = &rhs[(j2+0)*rhsStride];
1223  const Scalar* b1 = &rhs[(j2+1)*rhsStride];
1224  const Scalar* b2 = &rhs[(j2+2)*rhsStride];
1225  const Scalar* b3 = &rhs[(j2+3)*rhsStride];
1226  for(Index k=0; k<depth; k++)
1227  {
1228  blockB[count+0] = cj(b0[k]);
1229  blockB[count+1] = cj(b1[k]);
1230  if(nr==4) blockB[count+2] = cj(b2[k]);
1231  if(nr==4) blockB[count+3] = cj(b3[k]);
1232  count += nr;
1233  }
1234  // skip what we have after
1235  if(PanelMode) count += nr * (stride-offset-depth);
1236  }
1237 
1238  // copy the remaining columns one at a time (nr==1)
1239  for(Index j2=packet_cols; j2<cols; ++j2)
1240  {
1241  if(PanelMode) count += offset;
1242  const Scalar* b0 = &rhs[(j2+0)*rhsStride];
1243  for(Index k=0; k<depth; k++)
1244  {
1245  blockB[count] = cj(b0[k]);
1246  count += 1;
1247  }
1248  if(PanelMode) count += (stride-offset-depth);
1249  }
1250  }
1251 };
1252 
1253 // this version is optimized for row major matrices
1254 template<typename Scalar, typename Index, int nr, bool Conjugate, bool PanelMode>
1255 struct gemm_pack_rhs<Scalar, Index, nr, RowMajor, Conjugate, PanelMode>
1256 {
1257  enum { PacketSize = packet_traits<Scalar>::size };
1258  EIGEN_DONT_INLINE void operator()(Scalar* blockB, const Scalar* rhs, Index rhsStride, Index depth, Index cols,
1259  Index stride=0, Index offset=0)
1260  {
1261  EIGEN_ASM_COMMENT("EIGEN PRODUCT PACK RHS ROWMAJOR");
1262  eigen_assert(((!PanelMode) && stride==0 && offset==0) || (PanelMode && stride>=depth && offset<=stride));
1263  conj_if<NumTraits<Scalar>::IsComplex && Conjugate> cj;
1264  Index packet_cols = (cols/nr) * nr;
1265  Index count = 0;
1266  for(Index j2=0; j2<packet_cols; j2+=nr)
1267  {
1268  // skip what we have before
1269  if(PanelMode) count += nr * offset;
1270  for(Index k=0; k<depth; k++)
1271  {
1272  const Scalar* b0 = &rhs[k*rhsStride + j2];
1273  blockB[count+0] = cj(b0[0]);
1274  blockB[count+1] = cj(b0[1]);
1275  if(nr==4) blockB[count+2] = cj(b0[2]);
1276  if(nr==4) blockB[count+3] = cj(b0[3]);
1277  count += nr;
1278  }
1279  // skip what we have after
1280  if(PanelMode) count += nr * (stride-offset-depth);
1281  }
1282  // copy the remaining columns one at a time (nr==1)
1283  for(Index j2=packet_cols; j2<cols; ++j2)
1284  {
1285  if(PanelMode) count += offset;
1286  const Scalar* b0 = &rhs[j2];
1287  for(Index k=0; k<depth; k++)
1288  {
1289  blockB[count] = cj(b0[k*rhsStride]);
1290  count += 1;
1291  }
1292  if(PanelMode) count += stride-offset-depth;
1293  }
1294  }
1295 };
1296 
1297 } // end namespace internal
1298 
1301 inline std::ptrdiff_t l1CacheSize()
1302 {
1303  std::ptrdiff_t l1, l2;
1305  return l1;
1306 }
1307 
1310 inline std::ptrdiff_t l2CacheSize()
1311 {
1312  std::ptrdiff_t l1, l2;
1314  return l2;
1315 }
1316 
1322 inline void setCpuCacheSizes(std::ptrdiff_t l1, std::ptrdiff_t l2)
1323 {
1325 }
1326 
1327 } // end namespace Eigen
1328 
1329 #endif // EIGEN_GENERAL_BLOCK_PANEL_H