00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #ifndef EIGEN_GENERAL_MATRIX_MATRIX_H
00011 #define EIGEN_GENERAL_MATRIX_MATRIX_H
00012
00013 namespace Eigen {
00014
00015 namespace internal {
00016
00017 template<typename _LhsScalar, typename _RhsScalar> class level3_blocking;
00018
00019
00020 template<
00021 typename Index,
00022 typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
00023 typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs>
00024 struct general_matrix_matrix_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,RowMajor>
00025 {
00026 typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
00027 static EIGEN_STRONG_INLINE void run(
00028 Index rows, Index cols, Index depth,
00029 const LhsScalar* lhs, Index lhsStride,
00030 const RhsScalar* rhs, Index rhsStride,
00031 ResScalar* res, Index resStride,
00032 ResScalar alpha,
00033 level3_blocking<RhsScalar,LhsScalar>& blocking,
00034 GemmParallelInfo<Index>* info = 0)
00035 {
00036
00037 general_matrix_matrix_product<Index,
00038 RhsScalar, RhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateRhs,
00039 LhsScalar, LhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateLhs,
00040 ColMajor>
00041 ::run(cols,rows,depth,rhs,rhsStride,lhs,lhsStride,res,resStride,alpha,blocking,info);
00042 }
00043 };
00044
00045
00046
00047 template<
00048 typename Index,
00049 typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
00050 typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs>
00051 struct general_matrix_matrix_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,ColMajor>
00052 {
00053 typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
00054 static void run(Index rows, Index cols, Index depth,
00055 const LhsScalar* _lhs, Index lhsStride,
00056 const RhsScalar* _rhs, Index rhsStride,
00057 ResScalar* res, Index resStride,
00058 ResScalar alpha,
00059 level3_blocking<LhsScalar,RhsScalar>& blocking,
00060 GemmParallelInfo<Index>* info = 0)
00061 {
00062 const_blas_data_mapper<LhsScalar, Index, LhsStorageOrder> lhs(_lhs,lhsStride);
00063 const_blas_data_mapper<RhsScalar, Index, RhsStorageOrder> rhs(_rhs,rhsStride);
00064
00065 typedef gebp_traits<LhsScalar,RhsScalar> Traits;
00066
00067 Index kc = blocking.kc();
00068 Index mc = (std::min)(rows,blocking.mc());
00069
00070
00071 gemm_pack_lhs<LhsScalar, Index, Traits::mr, Traits::LhsProgress, LhsStorageOrder> pack_lhs;
00072 gemm_pack_rhs<RhsScalar, Index, Traits::nr, RhsStorageOrder> pack_rhs;
00073 gebp_kernel<LhsScalar, RhsScalar, Index, Traits::mr, Traits::nr, ConjugateLhs, ConjugateRhs> gebp;
00074
00075 #ifdef EIGEN_HAS_OPENMP
00076 if(info)
00077 {
00078
00079 Index tid = omp_get_thread_num();
00080 Index threads = omp_get_num_threads();
00081
00082 std::size_t sizeA = kc*mc;
00083 std::size_t sizeW = kc*Traits::WorkSpaceFactor;
00084 ei_declare_aligned_stack_constructed_variable(LhsScalar, blockA, sizeA, 0);
00085 ei_declare_aligned_stack_constructed_variable(RhsScalar, w, sizeW, 0);
00086
00087 RhsScalar* blockB = blocking.blockB();
00088 eigen_internal_assert(blockB!=0);
00089
00090
00091 for(Index k=0; k<depth; k+=kc)
00092 {
00093 const Index actual_kc = (std::min)(k+kc,depth)-k;
00094
00095
00096
00097 pack_lhs(blockA, &lhs(0,k), lhsStride, actual_kc, mc);
00098
00099
00100
00101
00102
00103
00104
00105 while(info[tid].users!=0) {}
00106 info[tid].users += threads;
00107
00108 pack_rhs(blockB+info[tid].rhs_start*actual_kc, &rhs(k,info[tid].rhs_start), rhsStride, actual_kc, info[tid].rhs_length);
00109
00110
00111 info[tid].sync = k;
00112
00113
00114 for(Index shift=0; shift<threads; ++shift)
00115 {
00116 Index j = (tid+shift)%threads;
00117
00118
00119
00120
00121 if(shift>0)
00122 while(info[j].sync!=k) {}
00123
00124 gebp(res+info[j].rhs_start*resStride, resStride, blockA, blockB+info[j].rhs_start*actual_kc, mc, actual_kc, info[j].rhs_length, alpha, -1,-1,0,0, w);
00125 }
00126
00127
00128 for(Index i=mc; i<rows; i+=mc)
00129 {
00130 const Index actual_mc = (std::min)(i+mc,rows)-i;
00131
00132
00133 pack_lhs(blockA, &lhs(i,k), lhsStride, actual_kc, actual_mc);
00134
00135
00136 gebp(res+i, resStride, blockA, blockB, actual_mc, actual_kc, cols, alpha, -1,-1,0,0, w);
00137 }
00138
00139
00140
00141 for(Index j=0; j<threads; ++j)
00142 #pragma omp atomic
00143 --(info[j].users);
00144 }
00145 }
00146 else
00147 #endif // EIGEN_HAS_OPENMP
00148 {
00149 EIGEN_UNUSED_VARIABLE(info);
00150
00151
00152 std::size_t sizeA = kc*mc;
00153 std::size_t sizeB = kc*cols;
00154 std::size_t sizeW = kc*Traits::WorkSpaceFactor;
00155
00156 ei_declare_aligned_stack_constructed_variable(LhsScalar, blockA, sizeA, blocking.blockA());
00157 ei_declare_aligned_stack_constructed_variable(RhsScalar, blockB, sizeB, blocking.blockB());
00158 ei_declare_aligned_stack_constructed_variable(RhsScalar, blockW, sizeW, blocking.blockW());
00159
00160
00161
00162 for(Index k2=0; k2<depth; k2+=kc)
00163 {
00164 const Index actual_kc = (std::min)(k2+kc,depth)-k2;
00165
00166
00167
00168
00169
00170 pack_rhs(blockB, &rhs(k2,0), rhsStride, actual_kc, cols);
00171
00172
00173
00174
00175 for(Index i2=0; i2<rows; i2+=mc)
00176 {
00177 const Index actual_mc = (std::min)(i2+mc,rows)-i2;
00178
00179
00180
00181
00182 pack_lhs(blockA, &lhs(i2,k2), lhsStride, actual_kc, actual_mc);
00183
00184
00185 gebp(res+i2, resStride, blockA, blockB, actual_mc, actual_kc, cols, alpha, -1, -1, 0, 0, blockW);
00186
00187 }
00188 }
00189 }
00190 }
00191
00192 };
00193
00194
00195
00196
00197
00198
00199 template<typename Lhs, typename Rhs>
00200 struct traits<GeneralProduct<Lhs,Rhs,GemmProduct> >
00201 : traits<ProductBase<GeneralProduct<Lhs,Rhs,GemmProduct>, Lhs, Rhs> >
00202 {};
00203
00204 template<typename Scalar, typename Index, typename Gemm, typename Lhs, typename Rhs, typename Dest, typename BlockingType>
00205 struct gemm_functor
00206 {
00207 gemm_functor(const Lhs& lhs, const Rhs& rhs, Dest& dest, Scalar actualAlpha,
00208 BlockingType& blocking)
00209 : m_lhs(lhs), m_rhs(rhs), m_dest(dest), m_actualAlpha(actualAlpha), m_blocking(blocking)
00210 {}
00211
00212 void initParallelSession() const
00213 {
00214 m_blocking.allocateB();
00215 }
00216
00217 void operator() (Index row, Index rows, Index col=0, Index cols=-1, GemmParallelInfo<Index>* info=0) const
00218 {
00219 if(cols==-1)
00220 cols = m_rhs.cols();
00221
00222 Gemm::run(rows, cols, m_lhs.cols(),
00223 &m_lhs.coeffRef(row,0), m_lhs.outerStride(),
00224 &m_rhs.coeffRef(0,col), m_rhs.outerStride(),
00225 (Scalar*)&(m_dest.coeffRef(row,col)), m_dest.outerStride(),
00226 m_actualAlpha, m_blocking, info);
00227 }
00228
00229 protected:
00230 const Lhs& m_lhs;
00231 const Rhs& m_rhs;
00232 Dest& m_dest;
00233 Scalar m_actualAlpha;
00234 BlockingType& m_blocking;
00235 };
00236
00237 template<int StorageOrder, typename LhsScalar, typename RhsScalar, int MaxRows, int MaxCols, int MaxDepth, int KcFactor=1,
00238 bool FiniteAtCompileTime = MaxRows!=Dynamic && MaxCols!=Dynamic && MaxDepth != Dynamic> class gemm_blocking_space;
00239
00240 template<typename _LhsScalar, typename _RhsScalar>
00241 class level3_blocking
00242 {
00243 typedef _LhsScalar LhsScalar;
00244 typedef _RhsScalar RhsScalar;
00245
00246 protected:
00247 LhsScalar* m_blockA;
00248 RhsScalar* m_blockB;
00249 RhsScalar* m_blockW;
00250
00251 DenseIndex m_mc;
00252 DenseIndex m_nc;
00253 DenseIndex m_kc;
00254
00255 public:
00256
00257 level3_blocking()
00258 : m_blockA(0), m_blockB(0), m_blockW(0), m_mc(0), m_nc(0), m_kc(0)
00259 {}
00260
00261 inline DenseIndex mc() const { return m_mc; }
00262 inline DenseIndex nc() const { return m_nc; }
00263 inline DenseIndex kc() const { return m_kc; }
00264
00265 inline LhsScalar* blockA() { return m_blockA; }
00266 inline RhsScalar* blockB() { return m_blockB; }
00267 inline RhsScalar* blockW() { return m_blockW; }
00268 };
00269
00270 template<int StorageOrder, typename _LhsScalar, typename _RhsScalar, int MaxRows, int MaxCols, int MaxDepth, int KcFactor>
00271 class gemm_blocking_space<StorageOrder,_LhsScalar,_RhsScalar,MaxRows, MaxCols, MaxDepth, KcFactor, true>
00272 : public level3_blocking<
00273 typename conditional<StorageOrder==RowMajor,_RhsScalar,_LhsScalar>::type,
00274 typename conditional<StorageOrder==RowMajor,_LhsScalar,_RhsScalar>::type>
00275 {
00276 enum {
00277 Transpose = StorageOrder==RowMajor,
00278 ActualRows = Transpose ? MaxCols : MaxRows,
00279 ActualCols = Transpose ? MaxRows : MaxCols
00280 };
00281 typedef typename conditional<Transpose,_RhsScalar,_LhsScalar>::type LhsScalar;
00282 typedef typename conditional<Transpose,_LhsScalar,_RhsScalar>::type RhsScalar;
00283 typedef gebp_traits<LhsScalar,RhsScalar> Traits;
00284 enum {
00285 SizeA = ActualRows * MaxDepth,
00286 SizeB = ActualCols * MaxDepth,
00287 SizeW = MaxDepth * Traits::WorkSpaceFactor
00288 };
00289
00290 EIGEN_ALIGN16 LhsScalar m_staticA[SizeA];
00291 EIGEN_ALIGN16 RhsScalar m_staticB[SizeB];
00292 EIGEN_ALIGN16 RhsScalar m_staticW[SizeW];
00293
00294 public:
00295
00296 gemm_blocking_space(DenseIndex , DenseIndex , DenseIndex )
00297 {
00298 this->m_mc = ActualRows;
00299 this->m_nc = ActualCols;
00300 this->m_kc = MaxDepth;
00301 this->m_blockA = m_staticA;
00302 this->m_blockB = m_staticB;
00303 this->m_blockW = m_staticW;
00304 }
00305
00306 inline void allocateA() {}
00307 inline void allocateB() {}
00308 inline void allocateW() {}
00309 inline void allocateAll() {}
00310 };
00311
00312 template<int StorageOrder, typename _LhsScalar, typename _RhsScalar, int MaxRows, int MaxCols, int MaxDepth, int KcFactor>
00313 class gemm_blocking_space<StorageOrder,_LhsScalar,_RhsScalar,MaxRows, MaxCols, MaxDepth, KcFactor, false>
00314 : public level3_blocking<
00315 typename conditional<StorageOrder==RowMajor,_RhsScalar,_LhsScalar>::type,
00316 typename conditional<StorageOrder==RowMajor,_LhsScalar,_RhsScalar>::type>
00317 {
00318 enum {
00319 Transpose = StorageOrder==RowMajor
00320 };
00321 typedef typename conditional<Transpose,_RhsScalar,_LhsScalar>::type LhsScalar;
00322 typedef typename conditional<Transpose,_LhsScalar,_RhsScalar>::type RhsScalar;
00323 typedef gebp_traits<LhsScalar,RhsScalar> Traits;
00324
00325 DenseIndex m_sizeA;
00326 DenseIndex m_sizeB;
00327 DenseIndex m_sizeW;
00328
00329 public:
00330
00331 gemm_blocking_space(DenseIndex rows, DenseIndex cols, DenseIndex depth)
00332 {
00333 this->m_mc = Transpose ? cols : rows;
00334 this->m_nc = Transpose ? rows : cols;
00335 this->m_kc = depth;
00336
00337 computeProductBlockingSizes<LhsScalar,RhsScalar,KcFactor>(this->m_kc, this->m_mc, this->m_nc);
00338 m_sizeA = this->m_mc * this->m_kc;
00339 m_sizeB = this->m_kc * this->m_nc;
00340 m_sizeW = this->m_kc*Traits::WorkSpaceFactor;
00341 }
00342
00343 void allocateA()
00344 {
00345 if(this->m_blockA==0)
00346 this->m_blockA = aligned_new<LhsScalar>(m_sizeA);
00347 }
00348
00349 void allocateB()
00350 {
00351 if(this->m_blockB==0)
00352 this->m_blockB = aligned_new<RhsScalar>(m_sizeB);
00353 }
00354
00355 void allocateW()
00356 {
00357 if(this->m_blockW==0)
00358 this->m_blockW = aligned_new<RhsScalar>(m_sizeW);
00359 }
00360
00361 void allocateAll()
00362 {
00363 allocateA();
00364 allocateB();
00365 allocateW();
00366 }
00367
00368 ~gemm_blocking_space()
00369 {
00370 aligned_delete(this->m_blockA, m_sizeA);
00371 aligned_delete(this->m_blockB, m_sizeB);
00372 aligned_delete(this->m_blockW, m_sizeW);
00373 }
00374 };
00375
00376 }
00377
00378 template<typename Lhs, typename Rhs>
00379 class GeneralProduct<Lhs, Rhs, GemmProduct>
00380 : public ProductBase<GeneralProduct<Lhs,Rhs,GemmProduct>, Lhs, Rhs>
00381 {
00382 enum {
00383 MaxDepthAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(Lhs::MaxColsAtCompileTime,Rhs::MaxRowsAtCompileTime)
00384 };
00385 public:
00386 EIGEN_PRODUCT_PUBLIC_INTERFACE(GeneralProduct)
00387
00388 typedef typename Lhs::Scalar LhsScalar;
00389 typedef typename Rhs::Scalar RhsScalar;
00390 typedef Scalar ResScalar;
00391
00392 GeneralProduct(const Lhs& lhs, const Rhs& rhs) : Base(lhs,rhs)
00393 {
00394 typedef internal::scalar_product_op<LhsScalar,RhsScalar> BinOp;
00395 EIGEN_CHECK_BINARY_COMPATIBILIY(BinOp,LhsScalar,RhsScalar);
00396 }
00397
00398 template<typename Dest> void scaleAndAddTo(Dest& dst, Scalar alpha) const
00399 {
00400 eigen_assert(dst.rows()==m_lhs.rows() && dst.cols()==m_rhs.cols());
00401
00402 typename internal::add_const_on_value_type<ActualLhsType>::type lhs = LhsBlasTraits::extract(m_lhs);
00403 typename internal::add_const_on_value_type<ActualRhsType>::type rhs = RhsBlasTraits::extract(m_rhs);
00404
00405 Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(m_lhs)
00406 * RhsBlasTraits::extractScalarFactor(m_rhs);
00407
00408 typedef internal::gemm_blocking_space<(Dest::Flags&RowMajorBit) ? RowMajor : ColMajor,LhsScalar,RhsScalar,
00409 Dest::MaxRowsAtCompileTime,Dest::MaxColsAtCompileTime,MaxDepthAtCompileTime> BlockingType;
00410
00411 typedef internal::gemm_functor<
00412 Scalar, Index,
00413 internal::general_matrix_matrix_product<
00414 Index,
00415 LhsScalar, (_ActualLhsType::Flags&RowMajorBit) ? RowMajor : ColMajor, bool(LhsBlasTraits::NeedToConjugate),
00416 RhsScalar, (_ActualRhsType::Flags&RowMajorBit) ? RowMajor : ColMajor, bool(RhsBlasTraits::NeedToConjugate),
00417 (Dest::Flags&RowMajorBit) ? RowMajor : ColMajor>,
00418 _ActualLhsType, _ActualRhsType, Dest, BlockingType> GemmFunctor;
00419
00420 BlockingType blocking(dst.rows(), dst.cols(), lhs.cols());
00421
00422 internal::parallelize_gemm<(Dest::MaxRowsAtCompileTime>32 || Dest::MaxRowsAtCompileTime==Dynamic)>(GemmFunctor(lhs, rhs, dst, actualAlpha, blocking), this->rows(), this->cols(), Dest::Flags&RowMajorBit);
00423 }
00424 };
00425
00426 }
00427
00428 #endif // EIGEN_GENERAL_MATRIX_MATRIX_H