// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2008-2009 Gael Guennebaud // // This Source Code Form is subject to the terms of the Mozilla // Public License v. 2.0. If a copy of the MPL was not distributed // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. #ifndef EIGEN_GENERAL_BLOCK_PANEL_H #define EIGEN_GENERAL_BLOCK_PANEL_H namespace Eigen { namespace internal { enum GEBPPacketSizeType { GEBPPacketFull = 0, GEBPPacketHalf, GEBPPacketQuarter }; template class gebp_traits; /** \internal \returns b if a<=0, and returns a otherwise. */ inline std::ptrdiff_t manage_caching_sizes_helper(std::ptrdiff_t a, std::ptrdiff_t b) { return a<=0 ? b : a; } #if defined(EIGEN_DEFAULT_L1_CACHE_SIZE) #define EIGEN_SET_DEFAULT_L1_CACHE_SIZE(val) EIGEN_DEFAULT_L1_CACHE_SIZE #else #define EIGEN_SET_DEFAULT_L1_CACHE_SIZE(val) val #endif // defined(EIGEN_DEFAULT_L1_CACHE_SIZE) #if defined(EIGEN_DEFAULT_L2_CACHE_SIZE) #define EIGEN_SET_DEFAULT_L2_CACHE_SIZE(val) EIGEN_DEFAULT_L2_CACHE_SIZE #else #define EIGEN_SET_DEFAULT_L2_CACHE_SIZE(val) val #endif // defined(EIGEN_DEFAULT_L2_CACHE_SIZE) #if defined(EIGEN_DEFAULT_L3_CACHE_SIZE) #define EIGEN_SET_DEFAULT_L3_CACHE_SIZE(val) EIGEN_SET_DEFAULT_L3_CACHE_SIZE #else #define EIGEN_SET_DEFAULT_L3_CACHE_SIZE(val) val #endif // defined(EIGEN_DEFAULT_L3_CACHE_SIZE) #if EIGEN_ARCH_i386_OR_x86_64 const std::ptrdiff_t defaultL1CacheSize = EIGEN_SET_DEFAULT_L1_CACHE_SIZE(32*1024); const std::ptrdiff_t defaultL2CacheSize = EIGEN_SET_DEFAULT_L2_CACHE_SIZE(256*1024); const std::ptrdiff_t defaultL3CacheSize = EIGEN_SET_DEFAULT_L3_CACHE_SIZE(2*1024*1024); #elif EIGEN_ARCH_PPC const std::ptrdiff_t defaultL1CacheSize = EIGEN_SET_DEFAULT_L1_CACHE_SIZE(64*1024); const std::ptrdiff_t defaultL2CacheSize = EIGEN_SET_DEFAULT_L2_CACHE_SIZE(512*1024); const std::ptrdiff_t defaultL3CacheSize = EIGEN_SET_DEFAULT_L3_CACHE_SIZE(4*1024*1024); #else const std::ptrdiff_t defaultL1CacheSize = EIGEN_SET_DEFAULT_L1_CACHE_SIZE(16*1024); const std::ptrdiff_t defaultL2CacheSize = EIGEN_SET_DEFAULT_L2_CACHE_SIZE(512*1024); const std::ptrdiff_t defaultL3CacheSize = EIGEN_SET_DEFAULT_L3_CACHE_SIZE(512*1024); #endif #undef EIGEN_SET_DEFAULT_L1_CACHE_SIZE #undef EIGEN_SET_DEFAULT_L2_CACHE_SIZE #undef EIGEN_SET_DEFAULT_L3_CACHE_SIZE /** \internal */ struct CacheSizes { CacheSizes(): m_l1(-1),m_l2(-1),m_l3(-1) { int l1CacheSize, l2CacheSize, l3CacheSize; queryCacheSizes(l1CacheSize, l2CacheSize, l3CacheSize); m_l1 = manage_caching_sizes_helper(l1CacheSize, defaultL1CacheSize); m_l2 = manage_caching_sizes_helper(l2CacheSize, defaultL2CacheSize); m_l3 = manage_caching_sizes_helper(l3CacheSize, defaultL3CacheSize); } std::ptrdiff_t m_l1; std::ptrdiff_t m_l2; std::ptrdiff_t m_l3; }; /** \internal */ inline void manage_caching_sizes(Action action, std::ptrdiff_t* l1, std::ptrdiff_t* l2, std::ptrdiff_t* l3) { static CacheSizes m_cacheSizes; if(action==SetAction) { // set the cpu cache size and cache all block sizes from a global cache size in byte eigen_internal_assert(l1!=0 && l2!=0); m_cacheSizes.m_l1 = *l1; m_cacheSizes.m_l2 = *l2; m_cacheSizes.m_l3 = *l3; } else if(action==GetAction) { eigen_internal_assert(l1!=0 && l2!=0); *l1 = m_cacheSizes.m_l1; *l2 = m_cacheSizes.m_l2; *l3 = m_cacheSizes.m_l3; } else { eigen_internal_assert(false); } } /* Helper for computeProductBlockingSizes. * * Given a m x k times k x n matrix product of scalar types \c LhsScalar and \c RhsScalar, * this function computes the blocking size parameters along the respective dimensions * for matrix products and related algorithms. The blocking sizes depends on various * parameters: * - the L1 and L2 cache sizes, * - the register level blocking sizes defined by gebp_traits, * - the number of scalars that fit into a packet (when vectorization is enabled). * * \sa setCpuCacheSizes */ template void evaluateProductBlockingSizesHeuristic(Index& k, Index& m, Index& n, Index num_threads = 1) { typedef gebp_traits Traits; // Explanations: // Let's recall that the product algorithms form mc x kc vertical panels A' on the lhs and // kc x nc blocks B' on the rhs. B' has to fit into L2/L3 cache. Moreover, A' is processed // per mr x kc horizontal small panels where mr is the blocking size along the m dimension // at the register level. This small horizontal panel has to stay within L1 cache. std::ptrdiff_t l1, l2, l3; manage_caching_sizes(GetAction, &l1, &l2, &l3); #ifdef EIGEN_VECTORIZE_AVX512 // We need to find a rationale for that, but without this adjustment, // performance with AVX512 is pretty bad, like -20% slower. // One reason is that with increasing packet-size, the blocking size k // has to become pretty small if we want that 1 lhs panel fit within L1. // For instance, with the 3pX4 kernel and double, the size of the lhs+rhs panels are: // k*(3*64 + 4*8) Bytes, with l1=32kBytes, and k%8=0, we have k=144. // This is quite small for a good reuse of the accumulation registers. l1 *= 4; #endif if (num_threads > 1) { typedef typename Traits::ResScalar ResScalar; enum { kdiv = KcFactor * (Traits::mr * sizeof(LhsScalar) + Traits::nr * sizeof(RhsScalar)), ksub = Traits::mr * Traits::nr * sizeof(ResScalar), kr = 8, mr = Traits::mr, nr = Traits::nr }; // Increasing k gives us more time to prefetch the content of the "C" // registers. However once the latency is hidden there is no point in // increasing the value of k, so we'll cap it at 320 (value determined // experimentally). // To avoid that k vanishes, we make k_cache at least as big as kr const Index k_cache = numext::maxi(kr, (numext::mini)((l1-ksub)/kdiv, 320)); if (k_cache < k) { k = k_cache - (k_cache % kr); eigen_internal_assert(k > 0); } const Index n_cache = (l2-l1) / (nr * sizeof(RhsScalar) * k); const Index n_per_thread = numext::div_ceil(n, num_threads); if (n_cache <= n_per_thread) { // Don't exceed the capacity of the l2 cache. eigen_internal_assert(n_cache >= static_cast(nr)); n = n_cache - (n_cache % nr); eigen_internal_assert(n > 0); } else { n = (numext::mini)(n, (n_per_thread + nr - 1) - ((n_per_thread + nr - 1) % nr)); } if (l3 > l2) { // l3 is shared between all cores, so we'll give each thread its own chunk of l3. const Index m_cache = (l3-l2) / (sizeof(LhsScalar) * k * num_threads); const Index m_per_thread = numext::div_ceil(m, num_threads); if(m_cache < m_per_thread && m_cache >= static_cast(mr)) { m = m_cache - (m_cache % mr); eigen_internal_assert(m > 0); } else { m = (numext::mini)(m, (m_per_thread + mr - 1) - ((m_per_thread + mr - 1) % mr)); } } } else { // In unit tests we do not want to use extra large matrices, // so we reduce the cache size to check the blocking strategy is not flawed #ifdef EIGEN_DEBUG_SMALL_PRODUCT_BLOCKS l1 = 9*1024; l2 = 32*1024; l3 = 512*1024; #endif // Early return for small problems because the computation below are time consuming for small problems. // Perhaps it would make more sense to consider k*n*m?? // Note that for very tiny problem, this function should be bypassed anyway // because we use the coefficient-based implementation for them. if((numext::maxi)(k,(numext::maxi)(m,n))<48) return; typedef typename Traits::ResScalar ResScalar; enum { k_peeling = 8, k_div = KcFactor * (Traits::mr * sizeof(LhsScalar) + Traits::nr * sizeof(RhsScalar)), k_sub = Traits::mr * Traits::nr * sizeof(ResScalar) }; // ---- 1st level of blocking on L1, yields kc ---- // Blocking on the third dimension (i.e., k) is chosen so that an horizontal panel // of size mr x kc of the lhs plus a vertical panel of kc x nr of the rhs both fits within L1 cache. // We also include a register-level block of the result (mx x nr). // (In an ideal world only the lhs panel would stay in L1) // Moreover, kc has to be a multiple of 8 to be compatible with loop peeling, leading to a maximum blocking size of: const Index max_kc = numext::maxi(((l1-k_sub)/k_div) & (~(k_peeling-1)),1); const Index old_k = k; if(k>max_kc) { // We are really blocking on the third dimension: // -> reduce blocking size to make sure the last block is as large as possible // while keeping the same number of sweeps over the result. k = (k%max_kc)==0 ? max_kc : max_kc - k_peeling * ((max_kc-1-(k%max_kc))/(k_peeling*(k/max_kc+1))); eigen_internal_assert(((old_k/k) == (old_k/max_kc)) && "the number of sweeps has to remain the same"); } // ---- 2nd level of blocking on max(L2,L3), yields nc ---- // TODO find a reliable way to get the actual amount of cache per core to use for 2nd level blocking, that is: // actual_l2 = max(l2, l3/nb_core_sharing_l3) // The number below is quite conservative: it is better to underestimate the cache size rather than overestimating it) // For instance, it corresponds to 6MB of L3 shared among 4 cores. #ifdef EIGEN_DEBUG_SMALL_PRODUCT_BLOCKS const Index actual_l2 = l3; #else const Index actual_l2 = 1572864; // == 1.5 MB #endif // Here, nc is chosen such that a block of kc x nc of the rhs fit within half of L2. // The second half is implicitly reserved to access the result and lhs coefficients. // When k= Index(Traits::nr*sizeof(RhsScalar))*k) { // L1 blocking max_nc = remaining_l1 / (k*sizeof(RhsScalar)); } else { // L2 blocking max_nc = (3*actual_l2)/(2*2*max_kc*sizeof(RhsScalar)); } // WARNING Below, we assume that Traits::nr is a power of two. Index nc = numext::mini(actual_l2/(2*k*sizeof(RhsScalar)), max_nc) & (~(Traits::nr-1)); if(n>nc) { // We are really blocking over the columns: // -> reduce blocking size to make sure the last block is as large as possible // while keeping the same number of sweeps over the packed lhs. // Here we allow one more sweep if this gives us a perfect match, thus the commented "-1" n = (n%nc)==0 ? nc : (nc - Traits::nr * ((nc/*-1*/-(n%nc))/(Traits::nr*(n/nc+1)))); } else if(old_k==k) { // So far, no blocking at all, i.e., kc==k, and nc==n. // In this case, let's perform a blocking over the rows such that the packed lhs data is kept in cache L1/L2 // TODO: part of this blocking strategy is now implemented within the kernel itself, so the L1-based heuristic here should be obsolete. Index problem_size = k*n*sizeof(LhsScalar); Index actual_lm = actual_l2; Index max_mc = m; if(problem_size<=1024) { // problem is small enough to keep in L1 // Let's choose m such that lhs's block fit in 1/3 of L1 actual_lm = l1; } else if(l3!=0 && problem_size<=32768) { // we have both L2 and L3, and problem is small enough to be kept in L2 // Let's choose m such that lhs's block fit in 1/3 of L2 actual_lm = l2; max_mc = (numext::mini)(576,max_mc); } Index mc = (numext::mini)(actual_lm/(3*k*sizeof(LhsScalar)), max_mc); if (mc > Traits::mr) mc -= mc % Traits::mr; else if (mc==0) return; m = (m%mc)==0 ? mc : (mc - Traits::mr * ((mc/*-1*/-(m%mc))/(Traits::mr*(m/mc+1)))); } } } template inline bool useSpecificBlockingSizes(Index& k, Index& m, Index& n) { #ifdef EIGEN_TEST_SPECIFIC_BLOCKING_SIZES if (EIGEN_TEST_SPECIFIC_BLOCKING_SIZES) { k = numext::mini(k, EIGEN_TEST_SPECIFIC_BLOCKING_SIZE_K); m = numext::mini(m, EIGEN_TEST_SPECIFIC_BLOCKING_SIZE_M); n = numext::mini(n, EIGEN_TEST_SPECIFIC_BLOCKING_SIZE_N); return true; } #else EIGEN_UNUSED_VARIABLE(k) EIGEN_UNUSED_VARIABLE(m) EIGEN_UNUSED_VARIABLE(n) #endif return false; } /** \brief Computes the blocking parameters for a m x k times k x n matrix product * * \param[in,out] k Input: the third dimension of the product. Output: the blocking size along the same dimension. * \param[in,out] m Input: the number of rows of the left hand side. Output: the blocking size along the same dimension. * \param[in,out] n Input: the number of columns of the right hand side. Output: the blocking size along the same dimension. * * Given a m x k times k x n matrix product of scalar types \c LhsScalar and \c RhsScalar, * this function computes the blocking size parameters along the respective dimensions * for matrix products and related algorithms. * * The blocking size parameters may be evaluated: * - either by a heuristic based on cache sizes; * - or using fixed prescribed values (for testing purposes). * * \sa setCpuCacheSizes */ template void computeProductBlockingSizes(Index& k, Index& m, Index& n, Index num_threads = 1) { if (!useSpecificBlockingSizes(k, m, n)) { evaluateProductBlockingSizesHeuristic(k, m, n, num_threads); } } template inline void computeProductBlockingSizes(Index& k, Index& m, Index& n, Index num_threads = 1) { computeProductBlockingSizes(k, m, n, num_threads); } #ifdef EIGEN_HAS_SINGLE_INSTRUCTION_CJMADD #define CJMADD(CJ,A,B,C,T) C = CJ.pmadd(A,B,C); #else // FIXME (a bit overkill maybe ?) template struct gebp_madd_selector { EIGEN_ALWAYS_INLINE static void run(const CJ& cj, A& a, B& b, C& c, T& /*t*/) { c = cj.pmadd(a,b,c); } }; template struct gebp_madd_selector { EIGEN_ALWAYS_INLINE static void run(const CJ& cj, T& a, T& b, T& c, T& t) { t = b; t = cj.pmul(a,t); c = padd(c,t); } }; template EIGEN_STRONG_INLINE void gebp_madd(const CJ& cj, A& a, B& b, C& c, T& t) { gebp_madd_selector::run(cj,a,b,c,t); } #define CJMADD(CJ,A,B,C,T) gebp_madd(CJ,A,B,C,T); // #define CJMADD(CJ,A,B,C,T) T = B; T = CJ.pmul(A,T); C = padd(C,T); #endif template struct RhsPanelHelper { private: static const int remaining_registers = EIGEN_ARCH_DEFAULT_NUMBER_OF_REGISTERS - registers_taken; public: typedef typename conditional=4, RhsPacketx4, RhsPacket>::type type; }; template struct QuadPacket { Packet B_0, B1, B2, B3; const Packet& get(const FixedInt<0>&) const { return B_0; } const Packet& get(const FixedInt<1>&) const { return B1; } const Packet& get(const FixedInt<2>&) const { return B2; } const Packet& get(const FixedInt<3>&) const { return B3; } }; template struct packet_conditional { typedef T3 type; }; template struct packet_conditional { typedef T1 type; }; template struct packet_conditional { typedef T2 type; }; #define PACKET_DECL_COND_PREFIX(prefix, name, packet_size) \ typedef typename packet_conditional::type, \ typename packet_traits::half, \ typename unpacket_traits::half>::half>::type \ prefix ## name ## Packet #define PACKET_DECL_COND(name, packet_size) \ typedef typename packet_conditional::type, \ typename packet_traits::half, \ typename unpacket_traits::half>::half>::type \ name ## Packet #define PACKET_DECL_COND_SCALAR_PREFIX(prefix, packet_size) \ typedef typename packet_conditional::type, \ typename packet_traits::half, \ typename unpacket_traits::half>::half>::type \ prefix ## ScalarPacket #define PACKET_DECL_COND_SCALAR(packet_size) \ typedef typename packet_conditional::type, \ typename packet_traits::half, \ typename unpacket_traits::half>::half>::type \ ScalarPacket /* Vectorization logic * real*real: unpack rhs to constant packets, ... * * 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), * storing each res packet into two packets (2x2), * at the end combine them: swap the second and addsub them * cf*cf : same but with 2x4 blocks * cplx*real : unpack rhs to constant packets, ... * real*cplx : load lhs as (a0,a0,a1,a1), and mul as usual */ template class gebp_traits { public: typedef _LhsScalar LhsScalar; typedef _RhsScalar RhsScalar; typedef typename ScalarBinaryOpTraits::ReturnType ResScalar; PACKET_DECL_COND_PREFIX(_, Lhs, _PacketSize); PACKET_DECL_COND_PREFIX(_, Rhs, _PacketSize); PACKET_DECL_COND_PREFIX(_, Res, _PacketSize); enum { ConjLhs = _ConjLhs, ConjRhs = _ConjRhs, Vectorizable = unpacket_traits<_LhsPacket>::vectorizable && unpacket_traits<_RhsPacket>::vectorizable, LhsPacketSize = Vectorizable ? unpacket_traits<_LhsPacket>::size : 1, RhsPacketSize = Vectorizable ? unpacket_traits<_RhsPacket>::size : 1, ResPacketSize = Vectorizable ? unpacket_traits<_ResPacket>::size : 1, NumberOfRegisters = EIGEN_ARCH_DEFAULT_NUMBER_OF_REGISTERS, // register block size along the N direction must be 1 or 4 nr = 4, // register block size along the M direction (currently, this one cannot be modified) default_mr = (EIGEN_PLAIN_ENUM_MIN(16,NumberOfRegisters)/2/nr)*LhsPacketSize, #if defined(EIGEN_HAS_SINGLE_INSTRUCTION_MADD) && !defined(EIGEN_VECTORIZE_ALTIVEC) && !defined(EIGEN_VECTORIZE_VSX) \ && ((!EIGEN_COMP_MSVC) || (EIGEN_COMP_MSVC>=1914)) // we assume 16 registers or more // See bug 992, if the scalar type is not vectorizable but that EIGEN_HAS_SINGLE_INSTRUCTION_MADD is defined, // then using 3*LhsPacketSize triggers non-implemented paths in syrk. // Bug 1515: MSVC prior to v19.14 yields to register spilling. mr = Vectorizable ? 3*LhsPacketSize : default_mr, #else mr = default_mr, #endif LhsProgress = LhsPacketSize, RhsProgress = 1 }; typedef typename conditional::type LhsPacket; typedef typename conditional::type RhsPacket; typedef typename conditional::type ResPacket; typedef LhsPacket LhsPacket4Packing; typedef QuadPacket RhsPacketx4; typedef ResPacket AccPacket; EIGEN_STRONG_INLINE void initAcc(AccPacket& p) { p = pset1(ResScalar(0)); } template EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, RhsPacketType& dest) const { dest = pset1(*b); } EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, RhsPacketx4& dest) const { pbroadcast4(b, dest.B_0, dest.B1, dest.B2, dest.B3); } template EIGEN_STRONG_INLINE void updateRhs(const RhsScalar* b, RhsPacketType& dest) const { loadRhs(b, dest); } EIGEN_STRONG_INLINE void updateRhs(const RhsScalar*, RhsPacketx4&) const { } EIGEN_STRONG_INLINE void loadRhsQuad(const RhsScalar* b, RhsPacket& dest) const { dest = ploadquad(b); } template EIGEN_STRONG_INLINE void loadLhs(const LhsScalar* a, LhsPacketType& dest) const { dest = pload(a); } template EIGEN_STRONG_INLINE void loadLhsUnaligned(const LhsScalar* a, LhsPacketType& dest) const { dest = ploadu(a); } template EIGEN_STRONG_INLINE void madd(const LhsPacketType& a, const RhsPacketType& b, AccPacketType& c, RhsPacketType& tmp, const LaneIdType&) const { conj_helper cj; // It would be a lot cleaner to call pmadd all the time. Unfortunately if we // let gcc allocate the register in which to store the result of the pmul // (in the case where there is no FMA) gcc fails to figure out how to avoid // spilling register. #ifdef EIGEN_HAS_SINGLE_INSTRUCTION_MADD EIGEN_UNUSED_VARIABLE(tmp); c = cj.pmadd(a,b,c); #else tmp = b; tmp = cj.pmul(a,tmp); c = padd(c,tmp); #endif } template EIGEN_STRONG_INLINE void madd(const LhsPacketType& a, const RhsPacketx4& b, AccPacketType& c, RhsPacket& tmp, const LaneIdType& lane) const { madd(a, b.get(lane), c, tmp, lane); } EIGEN_STRONG_INLINE void acc(const AccPacket& c, const ResPacket& alpha, ResPacket& r) const { r = pmadd(c,alpha,r); } template EIGEN_STRONG_INLINE void acc(const ResPacketHalf& c, const ResPacketHalf& alpha, ResPacketHalf& r) const { r = pmadd(c,alpha,r); } }; template class gebp_traits, RealScalar, _ConjLhs, false, Arch, _PacketSize> { public: typedef std::complex LhsScalar; typedef RealScalar RhsScalar; typedef typename ScalarBinaryOpTraits::ReturnType ResScalar; PACKET_DECL_COND_PREFIX(_, Lhs, _PacketSize); PACKET_DECL_COND_PREFIX(_, Rhs, _PacketSize); PACKET_DECL_COND_PREFIX(_, Res, _PacketSize); enum { ConjLhs = _ConjLhs, ConjRhs = false, Vectorizable = unpacket_traits<_LhsPacket>::vectorizable && unpacket_traits<_RhsPacket>::vectorizable, LhsPacketSize = Vectorizable ? unpacket_traits<_LhsPacket>::size : 1, RhsPacketSize = Vectorizable ? unpacket_traits<_RhsPacket>::size : 1, ResPacketSize = Vectorizable ? unpacket_traits<_ResPacket>::size : 1, NumberOfRegisters = EIGEN_ARCH_DEFAULT_NUMBER_OF_REGISTERS, nr = 4, #if defined(EIGEN_HAS_SINGLE_INSTRUCTION_MADD) && !defined(EIGEN_VECTORIZE_ALTIVEC) && !defined(EIGEN_VECTORIZE_VSX) // we assume 16 registers mr = 3*LhsPacketSize, #else mr = (EIGEN_PLAIN_ENUM_MIN(16,NumberOfRegisters)/2/nr)*LhsPacketSize, #endif LhsProgress = LhsPacketSize, RhsProgress = 1 }; typedef typename conditional::type LhsPacket; typedef typename conditional::type RhsPacket; typedef typename conditional::type ResPacket; typedef LhsPacket LhsPacket4Packing; typedef QuadPacket RhsPacketx4; typedef ResPacket AccPacket; EIGEN_STRONG_INLINE void initAcc(AccPacket& p) { p = pset1(ResScalar(0)); } template EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, RhsPacketType& dest) const { dest = pset1(*b); } EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, RhsPacketx4& dest) const { pbroadcast4(b, dest.B_0, dest.B1, dest.B2, dest.B3); } template EIGEN_STRONG_INLINE void updateRhs(const RhsScalar* b, RhsPacketType& dest) const { loadRhs(b, dest); } EIGEN_STRONG_INLINE void updateRhs(const RhsScalar*, RhsPacketx4&) const {} EIGEN_STRONG_INLINE void loadRhsQuad(const RhsScalar* b, RhsPacket& dest) const { loadRhsQuad_impl(b,dest, typename conditional::type()); } EIGEN_STRONG_INLINE void loadRhsQuad_impl(const RhsScalar* b, RhsPacket& dest, const true_type&) const { // FIXME we can do better! // what we want here is a ploadheight RhsScalar tmp[4] = {b[0],b[0],b[1],b[1]}; dest = ploadquad(tmp); } EIGEN_STRONG_INLINE void loadRhsQuad_impl(const RhsScalar* b, RhsPacket& dest, const false_type&) const { eigen_internal_assert(RhsPacketSize<=8); dest = pset1(*b); } EIGEN_STRONG_INLINE void loadLhs(const LhsScalar* a, LhsPacket& dest) const { dest = pload(a); } template EIGEN_STRONG_INLINE void loadLhsUnaligned(const LhsScalar* a, LhsPacketType& dest) const { dest = ploadu(a); } template EIGEN_STRONG_INLINE void madd(const LhsPacketType& a, const RhsPacketType& b, AccPacketType& c, RhsPacketType& tmp, const LaneIdType&) const { madd_impl(a, b, c, tmp, typename conditional::type()); } template EIGEN_STRONG_INLINE void madd_impl(const LhsPacketType& a, const RhsPacketType& b, AccPacketType& c, RhsPacketType& tmp, const true_type&) const { #ifdef EIGEN_HAS_SINGLE_INSTRUCTION_MADD EIGEN_UNUSED_VARIABLE(tmp); c.v = pmadd(a.v,b,c.v); #else tmp = b; tmp = pmul(a.v,tmp); c.v = padd(c.v,tmp); #endif } EIGEN_STRONG_INLINE void madd_impl(const LhsScalar& a, const RhsScalar& b, ResScalar& c, RhsScalar& /*tmp*/, const false_type&) const { c += a * b; } template EIGEN_STRONG_INLINE void madd(const LhsPacketType& a, const RhsPacketx4& b, AccPacketType& c, RhsPacket& tmp, const LaneIdType& lane) const { madd(a, b.get(lane), c, tmp, lane); } template EIGEN_STRONG_INLINE void acc(const AccPacketType& c, const ResPacketType& alpha, ResPacketType& r) const { conj_helper cj; r = cj.pmadd(c,alpha,r); } protected: }; template struct DoublePacket { Packet first; Packet second; }; template DoublePacket padd(const DoublePacket &a, const DoublePacket &b) { DoublePacket res; res.first = padd(a.first, b.first); res.second = padd(a.second,b.second); return res; } // note that for DoublePacket the "4" in "downto4" // corresponds to the number of complexes, so it means "8" // it terms of real coefficients. template const DoublePacket& predux_half_dowto4(const DoublePacket &a, typename enable_if::size<=8>::type* = 0) { return a; } template DoublePacket::half> predux_half_dowto4(const DoublePacket &a, typename enable_if::size==16>::type* = 0) { // yes, that's pretty hackish :( DoublePacket::half> res; typedef std::complex::type> Cplx; typedef typename packet_traits::type CplxPacket; res.first = predux_half_dowto4(CplxPacket(a.first)).v; res.second = predux_half_dowto4(CplxPacket(a.second)).v; return res; } // same here, "quad" actually means "8" in terms of real coefficients template void loadQuadToDoublePacket(const Scalar* b, DoublePacket& dest, typename enable_if::size<=8>::type* = 0) { dest.first = pset1(numext::real(*b)); dest.second = pset1(numext::imag(*b)); } template void loadQuadToDoublePacket(const Scalar* b, DoublePacket& dest, typename enable_if::size==16>::type* = 0) { // yes, that's pretty hackish too :( typedef typename NumTraits::Real RealScalar; RealScalar r[4] = {numext::real(b[0]), numext::real(b[0]), numext::real(b[1]), numext::real(b[1])}; RealScalar i[4] = {numext::imag(b[0]), numext::imag(b[0]), numext::imag(b[1]), numext::imag(b[1])}; dest.first = ploadquad(r); dest.second = ploadquad(i); } template struct unpacket_traits > { typedef DoublePacket::half> half; }; // template // DoublePacket pmadd(const DoublePacket &a, const DoublePacket &b) // { // DoublePacket res; // res.first = padd(a.first, b.first); // res.second = padd(a.second,b.second); // return res; // } template class gebp_traits, std::complex, _ConjLhs, _ConjRhs, Arch, _PacketSize > { public: typedef std::complex Scalar; typedef std::complex LhsScalar; typedef std::complex RhsScalar; typedef std::complex ResScalar; PACKET_DECL_COND_PREFIX(_, Lhs, _PacketSize); PACKET_DECL_COND_PREFIX(_, Rhs, _PacketSize); PACKET_DECL_COND_PREFIX(_, Res, _PacketSize); PACKET_DECL_COND(Real, _PacketSize); PACKET_DECL_COND_SCALAR(_PacketSize); enum { ConjLhs = _ConjLhs, ConjRhs = _ConjRhs, Vectorizable = unpacket_traits::vectorizable && unpacket_traits::vectorizable, ResPacketSize = Vectorizable ? unpacket_traits<_ResPacket>::size : 1, LhsPacketSize = Vectorizable ? unpacket_traits<_LhsPacket>::size : 1, RhsPacketSize = Vectorizable ? unpacket_traits::size : 1, RealPacketSize = Vectorizable ? unpacket_traits::size : 1, // FIXME: should depend on NumberOfRegisters nr = 4, mr = ResPacketSize, LhsProgress = ResPacketSize, RhsProgress = 1 }; typedef DoublePacket DoublePacketType; typedef typename conditional::type LhsPacket4Packing; typedef typename conditional::type LhsPacket; typedef typename conditional::type RhsPacket; typedef typename conditional::type ResPacket; typedef typename conditional::type AccPacket; // this actualy holds 8 packets! typedef QuadPacket RhsPacketx4; EIGEN_STRONG_INLINE void initAcc(Scalar& p) { p = Scalar(0); } EIGEN_STRONG_INLINE void initAcc(DoublePacketType& p) { p.first = pset1(RealScalar(0)); p.second = pset1(RealScalar(0)); } // Scalar path EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, ScalarPacket& dest) const { dest = pset1(*b); } // Vectorized path template EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, DoublePacket& dest) const { dest.first = pset1(numext::real(*b)); dest.second = pset1(numext::imag(*b)); } EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, RhsPacketx4& dest) const { loadRhs(b, dest.B_0); loadRhs(b + 1, dest.B1); loadRhs(b + 2, dest.B2); loadRhs(b + 3, dest.B3); } // Scalar path EIGEN_STRONG_INLINE void updateRhs(const RhsScalar* b, ScalarPacket& dest) const { loadRhs(b, dest); } // Vectorized path template EIGEN_STRONG_INLINE void updateRhs(const RhsScalar* b, DoublePacket& dest) const { loadRhs(b, dest); } EIGEN_STRONG_INLINE void updateRhs(const RhsScalar*, RhsPacketx4&) const {} EIGEN_STRONG_INLINE void loadRhsQuad(const RhsScalar* b, ResPacket& dest) const { loadRhs(b,dest); } EIGEN_STRONG_INLINE void loadRhsQuad(const RhsScalar* b, DoublePacketType& dest) const { loadQuadToDoublePacket(b,dest); } // nothing special here EIGEN_STRONG_INLINE void loadLhs(const LhsScalar* a, LhsPacket& dest) const { dest = pload((const typename unpacket_traits::type*)(a)); } template EIGEN_STRONG_INLINE void loadLhsUnaligned(const LhsScalar* a, LhsPacketType& dest) const { dest = ploadu((const typename unpacket_traits::type*)(a)); } template EIGEN_STRONG_INLINE typename enable_if::value>::type madd(const LhsPacketType& a, const RhsPacketType& b, DoublePacket& c, TmpType& /*tmp*/, const LaneIdType&) const { c.first = padd(pmul(a,b.first), c.first); c.second = padd(pmul(a,b.second),c.second); } template EIGEN_STRONG_INLINE void madd(const LhsPacket& a, const RhsPacket& b, ResPacket& c, RhsPacket& /*tmp*/, const LaneIdType&) const { c = cj.pmadd(a,b,c); } template EIGEN_STRONG_INLINE void madd(const LhsPacketType& a, const RhsPacketx4& b, AccPacketType& c, RhsPacket& tmp, const LaneIdType& lane) const { madd(a, b.get(lane), c, tmp, lane); } EIGEN_STRONG_INLINE void acc(const Scalar& c, const Scalar& alpha, Scalar& r) const { r += alpha * c; } template EIGEN_STRONG_INLINE void acc(const DoublePacket& c, const ResPacketType& alpha, ResPacketType& r) const { // assemble c ResPacketType tmp; if((!ConjLhs)&&(!ConjRhs)) { tmp = pcplxflip(pconj(ResPacketType(c.second))); tmp = padd(ResPacketType(c.first),tmp); } else if((!ConjLhs)&&(ConjRhs)) { tmp = pconj(pcplxflip(ResPacketType(c.second))); tmp = padd(ResPacketType(c.first),tmp); } else if((ConjLhs)&&(!ConjRhs)) { tmp = pcplxflip(ResPacketType(c.second)); tmp = padd(pconj(ResPacketType(c.first)),tmp); } else if((ConjLhs)&&(ConjRhs)) { tmp = pcplxflip(ResPacketType(c.second)); tmp = psub(pconj(ResPacketType(c.first)),tmp); } r = pmadd(tmp,alpha,r); } protected: conj_helper cj; }; template class gebp_traits, false, _ConjRhs, Arch, _PacketSize > { public: typedef std::complex Scalar; typedef RealScalar LhsScalar; typedef Scalar RhsScalar; typedef Scalar ResScalar; PACKET_DECL_COND_PREFIX(_, Lhs, _PacketSize); PACKET_DECL_COND_PREFIX(_, Rhs, _PacketSize); PACKET_DECL_COND_PREFIX(_, Res, _PacketSize); PACKET_DECL_COND_PREFIX(_, Real, _PacketSize); PACKET_DECL_COND_SCALAR_PREFIX(_, _PacketSize); #undef PACKET_DECL_COND_SCALAR_PREFIX #undef PACKET_DECL_COND_PREFIX #undef PACKET_DECL_COND_SCALAR #undef PACKET_DECL_COND enum { ConjLhs = false, ConjRhs = _ConjRhs, Vectorizable = unpacket_traits<_RealPacket>::vectorizable && unpacket_traits<_ScalarPacket>::vectorizable, LhsPacketSize = Vectorizable ? unpacket_traits<_LhsPacket>::size : 1, RhsPacketSize = Vectorizable ? unpacket_traits<_RhsPacket>::size : 1, ResPacketSize = Vectorizable ? unpacket_traits<_ResPacket>::size : 1, NumberOfRegisters = EIGEN_ARCH_DEFAULT_NUMBER_OF_REGISTERS, // FIXME: should depend on NumberOfRegisters nr = 4, mr = (EIGEN_PLAIN_ENUM_MIN(16,NumberOfRegisters)/2/nr)*ResPacketSize, LhsProgress = ResPacketSize, RhsProgress = 1 }; typedef typename conditional::type LhsPacket; typedef typename conditional::type RhsPacket; typedef typename conditional::type ResPacket; typedef LhsPacket LhsPacket4Packing; typedef QuadPacket RhsPacketx4; typedef ResPacket AccPacket; EIGEN_STRONG_INLINE void initAcc(AccPacket& p) { p = pset1(ResScalar(0)); } template EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, RhsPacketType& dest) const { dest = pset1(*b); } EIGEN_STRONG_INLINE void loadRhs(const RhsScalar* b, RhsPacketx4& dest) const { pbroadcast4(b, dest.B_0, dest.B1, dest.B2, dest.B3); } template EIGEN_STRONG_INLINE void updateRhs(const RhsScalar* b, RhsPacketType& dest) const { loadRhs(b, dest); } EIGEN_STRONG_INLINE void updateRhs(const RhsScalar*, RhsPacketx4&) const {} EIGEN_STRONG_INLINE void loadLhs(const LhsScalar* a, LhsPacket& dest) const { dest = ploaddup(a); } EIGEN_STRONG_INLINE void loadRhsQuad(const RhsScalar* b, RhsPacket& dest) const { dest = ploadquad(b); } template EIGEN_STRONG_INLINE void loadLhsUnaligned(const LhsScalar* a, LhsPacketType& dest) const { dest = ploaddup(a); } template EIGEN_STRONG_INLINE void madd(const LhsPacketType& a, const RhsPacketType& b, AccPacketType& c, RhsPacketType& tmp, const LaneIdType&) const { madd_impl(a, b, c, tmp, typename conditional::type()); } template EIGEN_STRONG_INLINE void madd_impl(const LhsPacketType& a, const RhsPacketType& b, AccPacketType& c, RhsPacketType& tmp, const true_type&) const { #ifdef EIGEN_HAS_SINGLE_INSTRUCTION_MADD EIGEN_UNUSED_VARIABLE(tmp); c.v = pmadd(a,b.v,c.v); #else tmp = b; tmp.v = pmul(a,tmp.v); c = padd(c,tmp); #endif } EIGEN_STRONG_INLINE void madd_impl(const LhsScalar& a, const RhsScalar& b, ResScalar& c, RhsScalar& /*tmp*/, const false_type&) const { c += a * b; } template EIGEN_STRONG_INLINE void madd(const LhsPacketType& a, const RhsPacketx4& b, AccPacketType& c, RhsPacket& tmp, const LaneIdType& lane) const { madd(a, b.get(lane), c, tmp, lane); } template EIGEN_STRONG_INLINE void acc(const AccPacketType& c, const ResPacketType& alpha, ResPacketType& r) const { conj_helper cj; r = cj.pmadd(alpha,c,r); } protected: }; /* optimized General packed Block * packed Panel product kernel * * Mixing type logic: C += A * B * | A | B | comments * |real |cplx | no vectorization yet, would require to pack A with duplication * |cplx |real | easy vectorization */ template struct gebp_kernel { typedef gebp_traits Traits; typedef gebp_traits HalfTraits; typedef gebp_traits QuarterTraits; typedef typename Traits::ResScalar ResScalar; typedef typename Traits::LhsPacket LhsPacket; typedef typename Traits::RhsPacket RhsPacket; typedef typename Traits::ResPacket ResPacket; typedef typename Traits::AccPacket AccPacket; typedef typename Traits::RhsPacketx4 RhsPacketx4; typedef typename RhsPanelHelper::type RhsPanel15; typedef gebp_traits SwappedTraits; typedef typename SwappedTraits::ResScalar SResScalar; typedef typename SwappedTraits::LhsPacket SLhsPacket; typedef typename SwappedTraits::RhsPacket SRhsPacket; typedef typename SwappedTraits::ResPacket SResPacket; typedef typename SwappedTraits::AccPacket SAccPacket; typedef typename HalfTraits::LhsPacket LhsPacketHalf; typedef typename HalfTraits::RhsPacket RhsPacketHalf; typedef typename HalfTraits::ResPacket ResPacketHalf; typedef typename HalfTraits::AccPacket AccPacketHalf; typedef typename QuarterTraits::LhsPacket LhsPacketQuarter; typedef typename QuarterTraits::RhsPacket RhsPacketQuarter; typedef typename QuarterTraits::ResPacket ResPacketQuarter; typedef typename QuarterTraits::AccPacket AccPacketQuarter; typedef typename DataMapper::LinearMapper LinearMapper; enum { Vectorizable = Traits::Vectorizable, LhsProgress = Traits::LhsProgress, LhsProgressHalf = HalfTraits::LhsProgress, LhsProgressQuarter = QuarterTraits::LhsProgress, RhsProgress = Traits::RhsProgress, RhsProgressHalf = HalfTraits::RhsProgress, RhsProgressQuarter = QuarterTraits::RhsProgress, ResPacketSize = Traits::ResPacketSize }; EIGEN_DONT_INLINE void operator()(const DataMapper& res, const LhsScalar* blockA, const RhsScalar* blockB, Index rows, Index depth, Index cols, ResScalar alpha, Index strideA=-1, Index strideB=-1, Index offsetA=0, Index offsetB=0); }; template::LhsProgress> struct last_row_process_16_packets { typedef gebp_traits Traits; typedef gebp_traits SwappedTraits; typedef typename Traits::ResScalar ResScalar; typedef typename SwappedTraits::LhsPacket SLhsPacket; typedef typename SwappedTraits::RhsPacket SRhsPacket; typedef typename SwappedTraits::ResPacket SResPacket; typedef typename SwappedTraits::AccPacket SAccPacket; EIGEN_STRONG_INLINE void operator()(const DataMapper& res, SwappedTraits &straits, const LhsScalar* blA, const RhsScalar* blB, Index depth, const Index endk, Index i, Index j2, ResScalar alpha, SAccPacket &C0) { EIGEN_UNUSED_VARIABLE(res); EIGEN_UNUSED_VARIABLE(straits); EIGEN_UNUSED_VARIABLE(blA); EIGEN_UNUSED_VARIABLE(blB); EIGEN_UNUSED_VARIABLE(depth); EIGEN_UNUSED_VARIABLE(endk); EIGEN_UNUSED_VARIABLE(i); EIGEN_UNUSED_VARIABLE(j2); EIGEN_UNUSED_VARIABLE(alpha); EIGEN_UNUSED_VARIABLE(C0); } }; template struct last_row_process_16_packets { typedef gebp_traits Traits; typedef gebp_traits SwappedTraits; typedef typename Traits::ResScalar ResScalar; typedef typename SwappedTraits::LhsPacket SLhsPacket; typedef typename SwappedTraits::RhsPacket SRhsPacket; typedef typename SwappedTraits::ResPacket SResPacket; typedef typename SwappedTraits::AccPacket SAccPacket; EIGEN_STRONG_INLINE void operator()(const DataMapper& res, SwappedTraits &straits, const LhsScalar* blA, const RhsScalar* blB, Index depth, const Index endk, Index i, Index j2, ResScalar alpha, SAccPacket &C0) { typedef typename unpacket_traits::half>::half SResPacketQuarter; typedef typename unpacket_traits::half>::half SLhsPacketQuarter; typedef typename unpacket_traits::half>::half SRhsPacketQuarter; typedef typename unpacket_traits::half>::half SAccPacketQuarter; SResPacketQuarter R = res.template gatherPacket(i, j2); SResPacketQuarter alphav = pset1(alpha); if (depth - endk > 0) { // We have to handle the last row(s) of the rhs, which // correspond to a half-packet SAccPacketQuarter c0 = predux_half_dowto4(predux_half_dowto4(C0)); for (Index kk = endk; kk < depth; kk++) { SLhsPacketQuarter a0; SRhsPacketQuarter b0; straits.loadLhsUnaligned(blB, a0); straits.loadRhs(blA, b0); straits.madd(a0,b0,c0,b0, fix<0>); blB += SwappedTraits::LhsProgress/4; blA += 1; } straits.acc(c0, alphav, R); } else { straits.acc(predux_half_dowto4(predux_half_dowto4(C0)), alphav, R); } res.scatterPacket(i, j2, R); } }; template struct lhs_process_one_packet { typedef typename GEBPTraits::RhsPacketx4 RhsPacketx4; EIGEN_STRONG_INLINE void peeled_kc_onestep(Index K, const LhsScalar* blA, const RhsScalar* blB, GEBPTraits traits, LhsPacket *A0, RhsPacketx4 *rhs_panel, RhsPacket *T0, AccPacket *C0, AccPacket *C1, AccPacket *C2, AccPacket *C3) { EIGEN_ASM_COMMENT("begin step of gebp micro kernel 1X4"); EIGEN_ASM_COMMENT("Note: these asm comments work around bug 935!"); traits.loadLhs(&blA[(0+1*K)*LhsProgress], *A0); traits.loadRhs(&blB[(0+4*K)*RhsProgress], *rhs_panel); traits.madd(*A0, *rhs_panel, *C0, *T0, fix<0>); traits.madd(*A0, *rhs_panel, *C1, *T0, fix<1>); traits.madd(*A0, *rhs_panel, *C2, *T0, fix<2>); traits.madd(*A0, *rhs_panel, *C3, *T0, fix<3>); #if EIGEN_GNUC_AT_LEAST(6,0) && defined(EIGEN_VECTORIZE_SSE) __asm__ ("" : "+x,m" (*A0)); #endif EIGEN_ASM_COMMENT("end step of gebp micro kernel 1X4"); } EIGEN_STRONG_INLINE void operator()( const DataMapper& res, const LhsScalar* blockA, const RhsScalar* blockB, ResScalar alpha, Index peelStart, Index peelEnd, Index strideA, Index strideB, Index offsetA, Index offsetB, int prefetch_res_offset, Index peeled_kc, Index pk, Index cols, Index depth, Index packet_cols4) { GEBPTraits traits; // loops on each largest micro horizontal panel of lhs // (LhsProgress x depth) for(Index i=peelStart; i(alpha); R0 = r0.template loadPacket(0); R1 = r1.template loadPacket(0); traits.acc(C0, alphav, R0); traits.acc(C1, alphav, R1); r0.storePacket(0, R0); r1.storePacket(0, R1); R0 = r2.template loadPacket(0); R1 = r3.template loadPacket(0); traits.acc(C2, alphav, R0); traits.acc(C3, alphav, R1); r2.storePacket(0, R0); r3.storePacket(0, R1); } // Deal with remaining columns of the rhs for(Index j2=packet_cols4; j2); \ EIGEN_ASM_COMMENT("end step of gebp micro kernel 1/half/quarterX1"); \ } while(false); EIGEN_GEBGP_ONESTEP(0); EIGEN_GEBGP_ONESTEP(1); EIGEN_GEBGP_ONESTEP(2); EIGEN_GEBGP_ONESTEP(3); EIGEN_GEBGP_ONESTEP(4); EIGEN_GEBGP_ONESTEP(5); EIGEN_GEBGP_ONESTEP(6); EIGEN_GEBGP_ONESTEP(7); blB += pk*RhsProgress; blA += pk*LhsProgress; EIGEN_ASM_COMMENT("end gebp micro kernel 1/half/quarterX1"); } // process remaining peeled loop for(Index k=peeled_kc; k(alpha); R0 = r0.template loadPacket(0); traits.acc(C0, alphav, R0); r0.storePacket(0, R0); } } } }; template struct lhs_process_fraction_of_packet : lhs_process_one_packet { EIGEN_STRONG_INLINE void peeled_kc_onestep(Index K, const LhsScalar* blA, const RhsScalar* blB, GEBPTraits traits, LhsPacket *A0, RhsPacket *B_0, RhsPacket *B1, RhsPacket *B2, RhsPacket *B3, AccPacket *C0, AccPacket *C1, AccPacket *C2, AccPacket *C3) { EIGEN_ASM_COMMENT("begin step of gebp micro kernel 1X4"); EIGEN_ASM_COMMENT("Note: these asm comments work around bug 935!"); traits.loadLhsUnaligned(&blA[(0+1*K)*(LhsProgress)], *A0); traits.broadcastRhs(&blB[(0+4*K)*RhsProgress], *B_0, *B1, *B2, *B3); traits.madd(*A0, *B_0, *C0, *B_0); traits.madd(*A0, *B1, *C1, *B1); traits.madd(*A0, *B2, *C2, *B2); traits.madd(*A0, *B3, *C3, *B3); EIGEN_ASM_COMMENT("end step of gebp micro kernel 1X4"); } }; template EIGEN_DONT_INLINE void gebp_kernel ::operator()(const DataMapper& res, const LhsScalar* blockA, const RhsScalar* blockB, Index rows, Index depth, Index cols, ResScalar alpha, Index strideA, Index strideB, Index offsetA, Index offsetB) { Traits traits; SwappedTraits straits; if(strideA==-1) strideA = depth; if(strideB==-1) strideB = depth; conj_helper cj; Index packet_cols4 = nr>=4 ? (cols/4) * 4 : 0; const Index peeled_mc3 = mr>=3*Traits::LhsProgress ? (rows/(3*LhsProgress))*(3*LhsProgress) : 0; const Index peeled_mc2 = mr>=2*Traits::LhsProgress ? peeled_mc3+((rows-peeled_mc3)/(2*LhsProgress))*(2*LhsProgress) : 0; const Index peeled_mc1 = mr>=1*Traits::LhsProgress ? peeled_mc2+((rows-peeled_mc2)/(1*LhsProgress))*(1*LhsProgress) : 0; const Index peeled_mc_half = mr>=LhsProgressHalf ? peeled_mc1+((rows-peeled_mc1)/(LhsProgressHalf))*(LhsProgressHalf) : 0; const Index peeled_mc_quarter = mr>=LhsProgressQuarter ? peeled_mc_half+((rows-peeled_mc_half)/(LhsProgressQuarter))*(LhsProgressQuarter) : 0; enum { pk = 8 }; // NOTE Such a large peeling factor is important for large matrices (~ +5% when >1000 on Haswell) const Index peeled_kc = depth & ~(pk-1); const int prefetch_res_offset = 32/sizeof(ResScalar); // const Index depth2 = depth & ~1; //---------- Process 3 * LhsProgress rows at once ---------- // This corresponds to 3*LhsProgress x nr register blocks. // Usually, make sense only with FMA if(mr>=3*Traits::LhsProgress) { // Here, the general idea is to loop on each largest micro horizontal panel of the lhs (3*Traits::LhsProgress x depth) // and on each largest micro vertical panel of the rhs (depth * nr). // Blocking sizes, i.e., 'depth' has been computed so that the micro horizontal panel of the lhs fit in L1. // However, if depth is too small, we can extend the number of rows of these horizontal panels. // This actual number of rows is computed as follow: const Index l1 = defaultL1CacheSize; // in Bytes, TODO, l1 should be passed to this function. // The max(1, ...) here is needed because we may be using blocking params larger than what our known l1 cache size // suggests we should be using: either because our known l1 cache size is inaccurate (e.g. on Android, we can only guess), // or because we are testing specific blocking sizes. const Index actual_panel_rows = (3*LhsProgress) * std::max(1,( (l1 - sizeof(ResScalar)*mr*nr - depth*nr*sizeof(RhsScalar)) / (depth * sizeof(LhsScalar) * 3*LhsProgress) )); for(Index i1=0; i1); \ traits.madd(A1, rhs_panel, C4, T0, fix<0>); \ traits.madd(A2, rhs_panel, C8, T0, fix<0>); \ traits.updateRhs(blB + (1+4*K) * Traits::RhsProgress, rhs_panel); \ traits.madd(A0, rhs_panel, C1, T0, fix<1>); \ traits.madd(A1, rhs_panel, C5, T0, fix<1>); \ traits.madd(A2, rhs_panel, C9, T0, fix<1>); \ traits.updateRhs(blB + (2+4*K) * Traits::RhsProgress, rhs_panel); \ traits.madd(A0, rhs_panel, C2, T0, fix<2>); \ traits.madd(A1, rhs_panel, C6, T0, fix<2>); \ traits.madd(A2, rhs_panel, C10, T0, fix<2>); \ traits.updateRhs(blB + (3+4*K) * Traits::RhsProgress, rhs_panel); \ traits.madd(A0, rhs_panel, C3, T0, fix<3>); \ traits.madd(A1, rhs_panel, C7, T0, fix<3>); \ traits.madd(A2, rhs_panel, C11, T0, fix<3>); \ EIGEN_ASM_COMMENT("end step of gebp micro kernel 3pX4"); \ } while (false) internal::prefetch(blB); EIGEN_GEBP_ONESTEP(0); EIGEN_GEBP_ONESTEP(1); EIGEN_GEBP_ONESTEP(2); EIGEN_GEBP_ONESTEP(3); EIGEN_GEBP_ONESTEP(4); EIGEN_GEBP_ONESTEP(5); EIGEN_GEBP_ONESTEP(6); EIGEN_GEBP_ONESTEP(7); blB += pk*4*RhsProgress; blA += pk*3*Traits::LhsProgress; EIGEN_ASM_COMMENT("end gebp micro kernel 3pX4"); } // process remaining peeled loop for(Index k=peeled_kc; k(alpha); R0 = r0.template loadPacket(0 * Traits::ResPacketSize); R1 = r0.template loadPacket(1 * Traits::ResPacketSize); R2 = r0.template loadPacket(2 * Traits::ResPacketSize); traits.acc(C0, alphav, R0); traits.acc(C4, alphav, R1); traits.acc(C8, alphav, R2); r0.storePacket(0 * Traits::ResPacketSize, R0); r0.storePacket(1 * Traits::ResPacketSize, R1); r0.storePacket(2 * Traits::ResPacketSize, R2); R0 = r1.template loadPacket(0 * Traits::ResPacketSize); R1 = r1.template loadPacket(1 * Traits::ResPacketSize); R2 = r1.template loadPacket(2 * Traits::ResPacketSize); traits.acc(C1, alphav, R0); traits.acc(C5, alphav, R1); traits.acc(C9, alphav, R2); r1.storePacket(0 * Traits::ResPacketSize, R0); r1.storePacket(1 * Traits::ResPacketSize, R1); r1.storePacket(2 * Traits::ResPacketSize, R2); R0 = r2.template loadPacket(0 * Traits::ResPacketSize); R1 = r2.template loadPacket(1 * Traits::ResPacketSize); R2 = r2.template loadPacket(2 * Traits::ResPacketSize); traits.acc(C2, alphav, R0); traits.acc(C6, alphav, R1); traits.acc(C10, alphav, R2); r2.storePacket(0 * Traits::ResPacketSize, R0); r2.storePacket(1 * Traits::ResPacketSize, R1); r2.storePacket(2 * Traits::ResPacketSize, R2); R0 = r3.template loadPacket(0 * Traits::ResPacketSize); R1 = r3.template loadPacket(1 * Traits::ResPacketSize); R2 = r3.template loadPacket(2 * Traits::ResPacketSize); traits.acc(C3, alphav, R0); traits.acc(C7, alphav, R1); traits.acc(C11, alphav, R2); r3.storePacket(0 * Traits::ResPacketSize, R0); r3.storePacket(1 * Traits::ResPacketSize, R1); r3.storePacket(2 * Traits::ResPacketSize, R2); } } // Deal with remaining columns of the rhs for(Index j2=packet_cols4; j2); \ traits.madd(A1, B_0, C4, B_0, fix<0>); \ traits.madd(A2, B_0, C8, B_0, fix<0>); \ EIGEN_ASM_COMMENT("end step of gebp micro kernel 3pX1"); \ } while (false) EIGEN_GEBGP_ONESTEP(0); EIGEN_GEBGP_ONESTEP(1); EIGEN_GEBGP_ONESTEP(2); EIGEN_GEBGP_ONESTEP(3); EIGEN_GEBGP_ONESTEP(4); EIGEN_GEBGP_ONESTEP(5); EIGEN_GEBGP_ONESTEP(6); EIGEN_GEBGP_ONESTEP(7); blB += int(pk) * int(RhsProgress); blA += int(pk) * 3 * int(Traits::LhsProgress); EIGEN_ASM_COMMENT("end gebp micro kernel 3pX1"); } // process remaining peeled loop for(Index k=peeled_kc; k(alpha); R0 = r0.template loadPacket(0 * Traits::ResPacketSize); R1 = r0.template loadPacket(1 * Traits::ResPacketSize); R2 = r0.template loadPacket(2 * Traits::ResPacketSize); traits.acc(C0, alphav, R0); traits.acc(C4, alphav, R1); traits.acc(C8, alphav, R2); r0.storePacket(0 * Traits::ResPacketSize, R0); r0.storePacket(1 * Traits::ResPacketSize, R1); r0.storePacket(2 * Traits::ResPacketSize, R2); } } } } //---------- Process 2 * LhsProgress rows at once ---------- if(mr>=2*Traits::LhsProgress) { const Index l1 = defaultL1CacheSize; // in Bytes, TODO, l1 should be passed to this function. // The max(1, ...) here is needed because we may be using blocking params larger than what our known l1 cache size // suggests we should be using: either because our known l1 cache size is inaccurate (e.g. on Android, we can only guess), // or because we are testing specific blocking sizes. Index actual_panel_rows = (2*LhsProgress) * std::max(1,( (l1 - sizeof(ResScalar)*mr*nr - depth*nr*sizeof(RhsScalar)) / (depth * sizeof(LhsScalar) * 2*LhsProgress) )); for(Index i1=peeled_mc3; i1=6 without FMA (bug 1637) #if EIGEN_GNUC_AT_LEAST(6,0) && defined(EIGEN_VECTORIZE_SSE) #define EIGEN_GEBP_2PX4_SPILLING_WORKAROUND __asm__ ("" : [a0] "+x,m" (A0),[a1] "+x,m" (A1)); #else #define EIGEN_GEBP_2PX4_SPILLING_WORKAROUND #endif #define EIGEN_GEBGP_ONESTEP(K) \ do { \ EIGEN_ASM_COMMENT("begin step of gebp micro kernel 2pX4"); \ traits.loadLhs(&blA[(0 + 2 * K) * LhsProgress], A0); \ traits.loadLhs(&blA[(1 + 2 * K) * LhsProgress], A1); \ traits.loadRhs(&blB[(0 + 4 * K) * RhsProgress], rhs_panel); \ traits.madd(A0, rhs_panel, C0, T0, fix<0>); \ traits.madd(A1, rhs_panel, C4, T0, fix<0>); \ traits.madd(A0, rhs_panel, C1, T0, fix<1>); \ traits.madd(A1, rhs_panel, C5, T0, fix<1>); \ traits.madd(A0, rhs_panel, C2, T0, fix<2>); \ traits.madd(A1, rhs_panel, C6, T0, fix<2>); \ traits.madd(A0, rhs_panel, C3, T0, fix<3>); \ traits.madd(A1, rhs_panel, C7, T0, fix<3>); \ EIGEN_GEBP_2PX4_SPILLING_WORKAROUND \ EIGEN_ASM_COMMENT("end step of gebp micro kernel 2pX4"); \ } while (false) internal::prefetch(blB+(48+0)); EIGEN_GEBGP_ONESTEP(0); EIGEN_GEBGP_ONESTEP(1); EIGEN_GEBGP_ONESTEP(2); EIGEN_GEBGP_ONESTEP(3); internal::prefetch(blB+(48+16)); EIGEN_GEBGP_ONESTEP(4); EIGEN_GEBGP_ONESTEP(5); EIGEN_GEBGP_ONESTEP(6); EIGEN_GEBGP_ONESTEP(7); blB += pk*4*RhsProgress; blA += pk*(2*Traits::LhsProgress); EIGEN_ASM_COMMENT("end gebp micro kernel 2pX4"); } // process remaining peeled loop for(Index k=peeled_kc; k(alpha); R0 = r0.template loadPacket(0 * Traits::ResPacketSize); R1 = r0.template loadPacket(1 * Traits::ResPacketSize); R2 = r1.template loadPacket(0 * Traits::ResPacketSize); R3 = r1.template loadPacket(1 * Traits::ResPacketSize); traits.acc(C0, alphav, R0); traits.acc(C4, alphav, R1); traits.acc(C1, alphav, R2); traits.acc(C5, alphav, R3); r0.storePacket(0 * Traits::ResPacketSize, R0); r0.storePacket(1 * Traits::ResPacketSize, R1); r1.storePacket(0 * Traits::ResPacketSize, R2); r1.storePacket(1 * Traits::ResPacketSize, R3); R0 = r2.template loadPacket(0 * Traits::ResPacketSize); R1 = r2.template loadPacket(1 * Traits::ResPacketSize); R2 = r3.template loadPacket(0 * Traits::ResPacketSize); R3 = r3.template loadPacket(1 * Traits::ResPacketSize); traits.acc(C2, alphav, R0); traits.acc(C6, alphav, R1); traits.acc(C3, alphav, R2); traits.acc(C7, alphav, R3); r2.storePacket(0 * Traits::ResPacketSize, R0); r2.storePacket(1 * Traits::ResPacketSize, R1); r3.storePacket(0 * Traits::ResPacketSize, R2); r3.storePacket(1 * Traits::ResPacketSize, R3); } } // Deal with remaining columns of the rhs for(Index j2=packet_cols4; j2); \ traits.madd(A1, B_0, C4, B_0, fix<0>); \ EIGEN_ASM_COMMENT("end step of gebp micro kernel 2pX1"); \ } while(false) EIGEN_GEBGP_ONESTEP(0); EIGEN_GEBGP_ONESTEP(1); EIGEN_GEBGP_ONESTEP(2); EIGEN_GEBGP_ONESTEP(3); EIGEN_GEBGP_ONESTEP(4); EIGEN_GEBGP_ONESTEP(5); EIGEN_GEBGP_ONESTEP(6); EIGEN_GEBGP_ONESTEP(7); blB += int(pk) * int(RhsProgress); blA += int(pk) * 2 * int(Traits::LhsProgress); EIGEN_ASM_COMMENT("end gebp micro kernel 2pX1"); } // process remaining peeled loop for(Index k=peeled_kc; k(alpha); R0 = r0.template loadPacket(0 * Traits::ResPacketSize); R1 = r0.template loadPacket(1 * Traits::ResPacketSize); traits.acc(C0, alphav, R0); traits.acc(C4, alphav, R1); r0.storePacket(0 * Traits::ResPacketSize, R0); r0.storePacket(1 * Traits::ResPacketSize, R1); } } } } //---------- Process 1 * LhsProgress rows at once ---------- if(mr>=1*Traits::LhsProgress) { lhs_process_one_packet p; p(res, blockA, blockB, alpha, peeled_mc2, peeled_mc1, strideA, strideB, offsetA, offsetB, prefetch_res_offset, peeled_kc, pk, cols, depth, packet_cols4); } //---------- Process LhsProgressHalf rows at once ---------- if((LhsProgressHalf < LhsProgress) && mr>=LhsProgressHalf) { lhs_process_fraction_of_packet p; p(res, blockA, blockB, alpha, peeled_mc1, peeled_mc_half, strideA, strideB, offsetA, offsetB, prefetch_res_offset, peeled_kc, pk, cols, depth, packet_cols4); } //---------- Process LhsProgressQuarter rows at once ---------- if((LhsProgressQuarter < LhsProgressHalf) && mr>=LhsProgressQuarter) { lhs_process_fraction_of_packet p; p(res, blockA, blockB, alpha, peeled_mc_half, peeled_mc_quarter, strideA, strideB, offsetA, offsetB, prefetch_res_offset, peeled_kc, pk, cols, depth, packet_cols4); } //---------- Process remaining rows, 1 at once ---------- if(peeled_mc_quarter::half>::size; const int SResPacketQuarterSize = unpacket_traits::half>::half>::size; if ((SwappedTraits::LhsProgress % 4) == 0 && (SwappedTraits::LhsProgress<=16) && (SwappedTraits::LhsProgress!=8 || SResPacketHalfSize==nr) && (SwappedTraits::LhsProgress!=16 || SResPacketQuarterSize==nr)) { SAccPacket C0, C1, C2, C3; straits.initAcc(C0); straits.initAcc(C1); straits.initAcc(C2); straits.initAcc(C3); const Index spk = (std::max)(1,SwappedTraits::LhsProgress/4); const Index endk = (depth/spk)*spk; const Index endk4 = (depth/(spk*4))*(spk*4); Index k=0; for(; k); straits.madd(A1,B_1,C1,B_1, fix<0>); straits.loadLhsUnaligned(blB+2*SwappedTraits::LhsProgress, A0); straits.loadLhsUnaligned(blB+3*SwappedTraits::LhsProgress, A1); straits.loadRhsQuad(blA+2*spk, B_0); straits.loadRhsQuad(blA+3*spk, B_1); straits.madd(A0,B_0,C2,B_0, fix<0>); straits.madd(A1,B_1,C3,B_1, fix<0>); blB += 4*SwappedTraits::LhsProgress; blA += 4*spk; } C0 = padd(padd(C0,C1),padd(C2,C3)); for(; k); blB += SwappedTraits::LhsProgress; blA += spk; } if(SwappedTraits::LhsProgress==8) { // Special case where we have to first reduce the accumulation register C0 typedef typename conditional=8,typename unpacket_traits::half,SResPacket>::type SResPacketHalf; typedef typename conditional=8,typename unpacket_traits::half,SLhsPacket>::type SLhsPacketHalf; typedef typename conditional=8,typename unpacket_traits::half,SRhsPacket>::type SRhsPacketHalf; typedef typename conditional=8,typename unpacket_traits::half,SAccPacket>::type SAccPacketHalf; SResPacketHalf R = res.template gatherPacket(i, j2); SResPacketHalf alphav = pset1(alpha); if(depth-endk>0) { // We have to handle the last row of the rhs which corresponds to a half-packet SLhsPacketHalf a0; SRhsPacketHalf b0; straits.loadLhsUnaligned(blB, a0); straits.loadRhs(blA, b0); SAccPacketHalf c0 = predux_half_dowto4(C0); straits.madd(a0,b0,c0,b0, fix<0>); straits.acc(c0, alphav, R); } else { straits.acc(predux_half_dowto4(C0), alphav, R); } res.scatterPacket(i, j2, R); } else if (SwappedTraits::LhsProgress==16) { // Special case where we have to first reduce the // accumulation register C0. We specialize the block in // template form, so that LhsProgress < 16 paths don't // fail to compile last_row_process_16_packets p; p(res, straits, blA, blB, depth, endk, i, j2,alpha, C0); } else { SResPacket R = res.template gatherPacket(i, j2); SResPacket alphav = pset1(alpha); straits.acc(C0, alphav, R); res.scatterPacket(i, j2, R); } } else // scalar path { // get a 1 x 4 res block as registers ResScalar C0(0), C1(0), C2(0), C3(0); for(Index k=0; k struct gemm_pack_lhs { typedef typename DataMapper::LinearMapper LinearMapper; EIGEN_DONT_INLINE void operator()(Scalar* blockA, const DataMapper& lhs, Index depth, Index rows, Index stride=0, Index offset=0); }; template EIGEN_DONT_INLINE void gemm_pack_lhs ::operator()(Scalar* blockA, const DataMapper& lhs, Index depth, Index rows, Index stride, Index offset) { typedef typename unpacket_traits::half HalfPacket; typedef typename unpacket_traits::half>::half QuarterPacket; enum { PacketSize = unpacket_traits::size, HalfPacketSize = unpacket_traits::size, QuarterPacketSize = unpacket_traits::size, HasHalf = (int)HalfPacketSize < (int)PacketSize, HasQuarter = (int)QuarterPacketSize < (int)HalfPacketSize}; EIGEN_ASM_COMMENT("EIGEN PRODUCT PACK LHS"); EIGEN_UNUSED_VARIABLE(stride); EIGEN_UNUSED_VARIABLE(offset); eigen_assert(((!PanelMode) && stride==0 && offset==0) || (PanelMode && stride>=depth && offset<=stride)); eigen_assert( ((Pack1%PacketSize)==0 && Pack1<=4*PacketSize) || (Pack1<=4) ); conj_if::IsComplex && Conjugate> cj; Index count = 0; const Index peeled_mc3 = Pack1>=3*PacketSize ? (rows/(3*PacketSize))*(3*PacketSize) : 0; const Index peeled_mc2 = Pack1>=2*PacketSize ? peeled_mc3+((rows-peeled_mc3)/(2*PacketSize))*(2*PacketSize) : 0; const Index peeled_mc1 = Pack1>=1*PacketSize ? peeled_mc2+((rows-peeled_mc2)/(1*PacketSize))*(1*PacketSize) : 0; const Index peeled_mc_half = Pack1>=HalfPacketSize ? peeled_mc1+((rows-peeled_mc1)/(HalfPacketSize))*(HalfPacketSize) : 0; const Index peeled_mc_quarter = Pack1>=QuarterPacketSize ? (rows/(QuarterPacketSize))*(QuarterPacketSize) : 0; const Index last_lhs_progress = rows > peeled_mc_quarter ? (rows - peeled_mc_quarter) & ~1 : 0; const Index peeled_mc0 = Pack2>=PacketSize ? peeled_mc_quarter : Pack2>1 && last_lhs_progress ? (rows/last_lhs_progress)*last_lhs_progress : 0; Index i=0; // Pack 3 packets if(Pack1>=3*PacketSize) { for(; i(i+0*PacketSize, k); B = lhs.template loadPacket(i+1*PacketSize, k); C = lhs.template loadPacket(i+2*PacketSize, k); pstore(blockA+count, cj.pconj(A)); count+=PacketSize; pstore(blockA+count, cj.pconj(B)); count+=PacketSize; pstore(blockA+count, cj.pconj(C)); count+=PacketSize; } if(PanelMode) count += (3*PacketSize) * (stride-offset-depth); } } // Pack 2 packets if(Pack1>=2*PacketSize) { for(; i(i+0*PacketSize, k); B = lhs.template loadPacket(i+1*PacketSize, k); pstore(blockA+count, cj.pconj(A)); count+=PacketSize; pstore(blockA+count, cj.pconj(B)); count+=PacketSize; } if(PanelMode) count += (2*PacketSize) * (stride-offset-depth); } } // Pack 1 packets if(Pack1>=1*PacketSize) { for(; i(i+0*PacketSize, k); pstore(blockA+count, cj.pconj(A)); count+=PacketSize; } if(PanelMode) count += (1*PacketSize) * (stride-offset-depth); } } // Pack half packets if(HasHalf && Pack1>=HalfPacketSize) { for(; i(i+0*(HalfPacketSize), k); pstoreu(blockA+count, cj.pconj(A)); count+=HalfPacketSize; } if(PanelMode) count += (HalfPacketSize) * (stride-offset-depth); } } // Pack quarter packets if(HasQuarter && Pack1>=QuarterPacketSize) { for(; i(i+0*(QuarterPacketSize), k); pstoreu(blockA+count, cj.pconj(A)); count+=QuarterPacketSize; } if(PanelMode) count += (QuarterPacketSize) * (stride-offset-depth); } } // Pack2 may be *smaller* than PacketSize—that happens for // products like real * complex, where we have to go half the // progress on the lhs in order to duplicate those operands to // address both real & imaginary parts on the rhs. This portion will // pack those half ones until they match the number expected on the // last peeling loop at this point (for the rhs). if(Pack21) { for(; i struct gemm_pack_lhs { typedef typename DataMapper::LinearMapper LinearMapper; EIGEN_DONT_INLINE void operator()(Scalar* blockA, const DataMapper& lhs, Index depth, Index rows, Index stride=0, Index offset=0); }; template EIGEN_DONT_INLINE void gemm_pack_lhs ::operator()(Scalar* blockA, const DataMapper& lhs, Index depth, Index rows, Index stride, Index offset) { typedef typename unpacket_traits::half HalfPacket; typedef typename unpacket_traits::half>::half QuarterPacket; enum { PacketSize = unpacket_traits::size, HalfPacketSize = unpacket_traits::size, QuarterPacketSize = unpacket_traits::size, HasHalf = (int)HalfPacketSize < (int)PacketSize, HasQuarter = (int)QuarterPacketSize < (int)HalfPacketSize}; EIGEN_ASM_COMMENT("EIGEN PRODUCT PACK LHS"); EIGEN_UNUSED_VARIABLE(stride); EIGEN_UNUSED_VARIABLE(offset); eigen_assert(((!PanelMode) && stride==0 && offset==0) || (PanelMode && stride>=depth && offset<=stride)); conj_if::IsComplex && Conjugate> cj; Index count = 0; bool gone_half = false, gone_quarter = false, gone_last = false; Index i = 0; int pack = Pack1; int psize = PacketSize; while(pack>0) { Index remaining_rows = rows-i; Index peeled_mc = gone_last ? Pack2>1 ? (rows/pack)*pack : 0 : i+(remaining_rows/pack)*pack; Index starting_pos = i; for(; i=psize && psize >= QuarterPacketSize) { const Index peeled_k = (depth/psize)*psize; for(; k kernel; for (int p = 0; p < psize; ++p) kernel.packet[p] = lhs.template loadPacket(i+p+m, k); ptranspose(kernel); for (int p = 0; p < psize; ++p) pstore(blockA+count+m+(pack)*p, cj.pconj(kernel.packet[p])); } else if (HasHalf && psize == HalfPacketSize) { gone_half = true; PacketBlock kernel_half; for (int p = 0; p < psize; ++p) kernel_half.packet[p] = lhs.template loadPacket(i+p+m, k); ptranspose(kernel_half); for (int p = 0; p < psize; ++p) pstore(blockA+count+m+(pack)*p, cj.pconj(kernel_half.packet[p])); } else if (HasQuarter && psize == QuarterPacketSize) { gone_quarter = true; PacketBlock kernel_quarter; for (int p = 0; p < psize; ++p) kernel_quarter.packet[p] = lhs.template loadPacket(i+p+m, k); ptranspose(kernel_quarter); for (int p = 0; p < psize; ++p) pstore(blockA+count+m+(pack)*p, cj.pconj(kernel_quarter.packet[p])); } } count += psize*pack; } } for(; k= psize/2 || left >= psize/4) && ((psize/2 == HalfPacketSize && HasHalf && !gone_half) || (psize/2 == QuarterPacketSize && HasQuarter && !gone_quarter))) { psize /= 2; pack = psize; continue; } // Pack2 may be *smaller* than PacketSize—that happens for // products like real * complex, where we have to go half the // progress on the lhs in order to duplicate those operands to // address both real & imaginary parts on the rhs. This portion will // pack those half ones until they match the number expected on the // last peeling loop at this point (for the rhs). if (Pack2 < PacketSize && !gone_last) { gone_last = true; psize = pack = left & ~1; } } } for(; i struct gemm_pack_rhs { typedef typename packet_traits::type Packet; typedef typename DataMapper::LinearMapper LinearMapper; enum { PacketSize = packet_traits::size }; EIGEN_DONT_INLINE void operator()(Scalar* blockB, const DataMapper& rhs, Index depth, Index cols, Index stride=0, Index offset=0); }; template EIGEN_DONT_INLINE void gemm_pack_rhs ::operator()(Scalar* blockB, const DataMapper& rhs, Index depth, Index cols, Index stride, Index offset) { EIGEN_ASM_COMMENT("EIGEN PRODUCT PACK RHS COLMAJOR"); EIGEN_UNUSED_VARIABLE(stride); EIGEN_UNUSED_VARIABLE(offset); eigen_assert(((!PanelMode) && stride==0 && offset==0) || (PanelMode && stride>=depth && offset<=stride)); conj_if::IsComplex && Conjugate> cj; Index packet_cols8 = nr>=8 ? (cols/8) * 8 : 0; Index packet_cols4 = nr>=4 ? (cols/4) * 4 : 0; Index count = 0; const Index peeled_k = (depth/PacketSize)*PacketSize; // if(nr>=8) // { // for(Index j2=0; j2 kernel; // for (int p = 0; p < PacketSize; ++p) { // kernel.packet[p] = ploadu(&rhs[(j2+p)*rhsStride+k]); // } // ptranspose(kernel); // for (int p = 0; p < PacketSize; ++p) { // pstoreu(blockB+count, cj.pconj(kernel.packet[p])); // count+=PacketSize; // } // } // } // for(; k=4) { for(Index j2=packet_cols8; j2 kernel; kernel.packet[0 ] = dm0.template loadPacket(k); kernel.packet[1%PacketSize] = dm1.template loadPacket(k); kernel.packet[2%PacketSize] = dm2.template loadPacket(k); kernel.packet[3%PacketSize] = dm3.template loadPacket(k); ptranspose(kernel); pstoreu(blockB+count+0*PacketSize, cj.pconj(kernel.packet[0])); pstoreu(blockB+count+1*PacketSize, cj.pconj(kernel.packet[1%PacketSize])); pstoreu(blockB+count+2*PacketSize, cj.pconj(kernel.packet[2%PacketSize])); pstoreu(blockB+count+3*PacketSize, cj.pconj(kernel.packet[3%PacketSize])); count+=4*PacketSize; } } for(; k struct gemm_pack_rhs { typedef typename packet_traits::type Packet; typedef typename unpacket_traits::half HalfPacket; typedef typename unpacket_traits::half>::half QuarterPacket; typedef typename DataMapper::LinearMapper LinearMapper; enum { PacketSize = packet_traits::size, HalfPacketSize = unpacket_traits::size, QuarterPacketSize = unpacket_traits::size}; EIGEN_DONT_INLINE void operator()(Scalar* blockB, const DataMapper& rhs, Index depth, Index cols, Index stride=0, Index offset=0) { EIGEN_ASM_COMMENT("EIGEN PRODUCT PACK RHS ROWMAJOR"); EIGEN_UNUSED_VARIABLE(stride); EIGEN_UNUSED_VARIABLE(offset); eigen_assert(((!PanelMode) && stride==0 && offset==0) || (PanelMode && stride>=depth && offset<=stride)); const bool HasHalf = (int)HalfPacketSize < (int)PacketSize; const bool HasQuarter = (int)QuarterPacketSize < (int)HalfPacketSize; conj_if::IsComplex && Conjugate> cj; Index packet_cols8 = nr>=8 ? (cols/8) * 8 : 0; Index packet_cols4 = nr>=4 ? (cols/4) * 4 : 0; Index count = 0; // if(nr>=8) // { // for(Index j2=0; j2(&rhs[k*rhsStride + j2]); // pstoreu(blockB+count, cj.pconj(A)); // } else if (PacketSize==4) { // Packet A = ploadu(&rhs[k*rhsStride + j2]); // Packet B = ploadu(&rhs[k*rhsStride + j2 + PacketSize]); // pstoreu(blockB+count, cj.pconj(A)); // pstoreu(blockB+count+PacketSize, cj.pconj(B)); // } else { // const Scalar* b0 = &rhs[k*rhsStride + j2]; // blockB[count+0] = cj(b0[0]); // blockB[count+1] = cj(b0[1]); // blockB[count+2] = cj(b0[2]); // blockB[count+3] = cj(b0[3]); // blockB[count+4] = cj(b0[4]); // blockB[count+5] = cj(b0[5]); // blockB[count+6] = cj(b0[6]); // blockB[count+7] = cj(b0[7]); // } // count += 8; // } // // skip what we have after // if(PanelMode) count += 8 * (stride-offset-depth); // } // } if(nr>=4) { for(Index j2=packet_cols8; j2(k, j2); pstoreu(blockB+count, cj.pconj(A)); count += PacketSize; } else if (HasHalf && HalfPacketSize==4) { HalfPacket A = rhs.template loadPacket(k, j2); pstoreu(blockB+count, cj.pconj(A)); count += HalfPacketSize; } else if (HasQuarter && QuarterPacketSize==4) { QuarterPacket A = rhs.template loadPacket(k, j2); pstoreu(blockB+count, cj.pconj(A)); count += QuarterPacketSize; } else { const LinearMapper dm0 = rhs.getLinearMapper(k, j2); blockB[count+0] = cj(dm0(0)); blockB[count+1] = cj(dm0(1)); blockB[count+2] = cj(dm0(2)); blockB[count+3] = cj(dm0(3)); count += 4; } } // skip what we have after if(PanelMode) count += 4 * (stride-offset-depth); } } // copy the remaining columns one at a time (nr==1) for(Index j2=packet_cols4; j2