Home | History | Annotate | Download | only in products
      1 // This file is part of Eigen, a lightweight C++ template library
      2 // for linear algebra.
      3 //
      4 // Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud (at) inria.fr>
      5 //
      6 // This Source Code Form is subject to the terms of the Mozilla
      7 // Public License v. 2.0. If a copy of the MPL was not distributed
      8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
      9 
     10 #ifndef EIGEN_GENERAL_MATRIX_MATRIX_TRIANGULAR_H
     11 #define EIGEN_GENERAL_MATRIX_MATRIX_TRIANGULAR_H
     12 
     13 namespace Eigen {
     14 
     15 namespace internal {
     16 
     17 /**********************************************************************
     18 * This file implements a general A * B product while
     19 * evaluating only one triangular part of the product.
     20 * This is more general version of self adjoint product (C += A A^T)
     21 * as the level 3 SYRK Blas routine.
     22 **********************************************************************/
     23 
     24 // forward declarations (defined at the end of this file)
     25 template<typename LhsScalar, typename RhsScalar, typename Index, int mr, int nr, bool ConjLhs, bool ConjRhs, int UpLo>
     26 struct tribb_kernel;
     27 
     28 /* Optimized matrix-matrix product evaluating only one triangular half */
     29 template <typename Index,
     30           typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
     31           typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs,
     32                               int ResStorageOrder, int  UpLo, int Version = Specialized>
     33 struct general_matrix_matrix_triangular_product;
     34 
     35 // as usual if the result is row major => we transpose the product
     36 template <typename Index, typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
     37                           typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs, int  UpLo, int Version>
     38 struct general_matrix_matrix_triangular_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,RowMajor,UpLo,Version>
     39 {
     40   typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
     41   static EIGEN_STRONG_INLINE void run(Index size, Index depth,const LhsScalar* lhs, Index lhsStride,
     42                                       const RhsScalar* rhs, Index rhsStride, ResScalar* res, Index resStride, ResScalar alpha)
     43   {
     44     general_matrix_matrix_triangular_product<Index,
     45         RhsScalar, RhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateRhs,
     46         LhsScalar, LhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateLhs,
     47         ColMajor, UpLo==Lower?Upper:Lower>
     48       ::run(size,depth,rhs,rhsStride,lhs,lhsStride,res,resStride,alpha);
     49   }
     50 };
     51 
     52 template <typename Index, typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
     53                           typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs, int  UpLo, int Version>
     54 struct general_matrix_matrix_triangular_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,ColMajor,UpLo,Version>
     55 {
     56   typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
     57   static EIGEN_STRONG_INLINE void run(Index size, Index depth,const LhsScalar* _lhs, Index lhsStride,
     58                                       const RhsScalar* _rhs, Index rhsStride, ResScalar* res, Index resStride, ResScalar alpha)
     59   {
     60     const_blas_data_mapper<LhsScalar, Index, LhsStorageOrder> lhs(_lhs,lhsStride);
     61     const_blas_data_mapper<RhsScalar, Index, RhsStorageOrder> rhs(_rhs,rhsStride);
     62 
     63     typedef gebp_traits<LhsScalar,RhsScalar> Traits;
     64 
     65     Index kc = depth; // cache block size along the K direction
     66     Index mc = size;  // cache block size along the M direction
     67     Index nc = size;  // cache block size along the N direction
     68     computeProductBlockingSizes<LhsScalar,RhsScalar>(kc, mc, nc);
     69     // !!! mc must be a multiple of nr:
     70     if(mc > Traits::nr)
     71       mc = (mc/Traits::nr)*Traits::nr;
     72 
     73     std::size_t sizeW = kc*Traits::WorkSpaceFactor;
     74     std::size_t sizeB = sizeW + kc*size;
     75     ei_declare_aligned_stack_constructed_variable(LhsScalar, blockA, kc*mc, 0);
     76     ei_declare_aligned_stack_constructed_variable(RhsScalar, allocatedBlockB, sizeB, 0);
     77     RhsScalar* blockB = allocatedBlockB + sizeW;
     78 
     79     gemm_pack_lhs<LhsScalar, Index, Traits::mr, Traits::LhsProgress, LhsStorageOrder> pack_lhs;
     80     gemm_pack_rhs<RhsScalar, Index, Traits::nr, RhsStorageOrder> pack_rhs;
     81     gebp_kernel <LhsScalar, RhsScalar, Index, Traits::mr, Traits::nr, ConjugateLhs, ConjugateRhs> gebp;
     82     tribb_kernel<LhsScalar, RhsScalar, Index, Traits::mr, Traits::nr, ConjugateLhs, ConjugateRhs, UpLo> sybb;
     83 
     84     for(Index k2=0; k2<depth; k2+=kc)
     85     {
     86       const Index actual_kc = (std::min)(k2+kc,depth)-k2;
     87 
     88       // note that the actual rhs is the transpose/adjoint of mat
     89       pack_rhs(blockB, &rhs(k2,0), rhsStride, actual_kc, size);
     90 
     91       for(Index i2=0; i2<size; i2+=mc)
     92       {
     93         const Index actual_mc = (std::min)(i2+mc,size)-i2;
     94 
     95         pack_lhs(blockA, &lhs(i2, k2), lhsStride, actual_kc, actual_mc);
     96 
     97         // the selected actual_mc * size panel of res is split into three different part:
     98         //  1 - before the diagonal => processed with gebp or skipped
     99         //  2 - the actual_mc x actual_mc symmetric block => processed with a special kernel
    100         //  3 - after the diagonal => processed with gebp or skipped
    101         if (UpLo==Lower)
    102           gebp(res+i2, resStride, blockA, blockB, actual_mc, actual_kc, (std::min)(size,i2), alpha,
    103                -1, -1, 0, 0, allocatedBlockB);
    104 
    105         sybb(res+resStride*i2 + i2, resStride, blockA, blockB + actual_kc*i2, actual_mc, actual_kc, alpha, allocatedBlockB);
    106 
    107         if (UpLo==Upper)
    108         {
    109           Index j2 = i2+actual_mc;
    110           gebp(res+resStride*j2+i2, resStride, blockA, blockB+actual_kc*j2, actual_mc, actual_kc, (std::max)(Index(0), size-j2), alpha,
    111                -1, -1, 0, 0, allocatedBlockB);
    112         }
    113       }
    114     }
    115   }
    116 };
    117 
    118 // Optimized packed Block * packed Block product kernel evaluating only one given triangular part
    119 // This kernel is built on top of the gebp kernel:
    120 // - the current destination block is processed per panel of actual_mc x BlockSize
    121 //   where BlockSize is set to the minimal value allowing gebp to be as fast as possible
    122 // - then, as usual, each panel is split into three parts along the diagonal,
    123 //   the sub blocks above and below the diagonal are processed as usual,
    124 //   while the triangular block overlapping the diagonal is evaluated into a
    125 //   small temporary buffer which is then accumulated into the result using a
    126 //   triangular traversal.
    127 template<typename LhsScalar, typename RhsScalar, typename Index, int mr, int nr, bool ConjLhs, bool ConjRhs, int UpLo>
    128 struct tribb_kernel
    129 {
    130   typedef gebp_traits<LhsScalar,RhsScalar,ConjLhs,ConjRhs> Traits;
    131   typedef typename Traits::ResScalar ResScalar;
    132 
    133   enum {
    134     BlockSize  = EIGEN_PLAIN_ENUM_MAX(mr,nr)
    135   };
    136   void operator()(ResScalar* res, Index resStride, const LhsScalar* blockA, const RhsScalar* blockB, Index size, Index depth, ResScalar alpha, RhsScalar* workspace)
    137   {
    138     gebp_kernel<LhsScalar, RhsScalar, Index, mr, nr, ConjLhs, ConjRhs> gebp_kernel;
    139     Matrix<ResScalar,BlockSize,BlockSize,ColMajor> buffer;
    140 
    141     // let's process the block per panel of actual_mc x BlockSize,
    142     // again, each is split into three parts, etc.
    143     for (Index j=0; j<size; j+=BlockSize)
    144     {
    145       Index actualBlockSize = std::min<Index>(BlockSize,size - j);
    146       const RhsScalar* actual_b = blockB+j*depth;
    147 
    148       if(UpLo==Upper)
    149         gebp_kernel(res+j*resStride, resStride, blockA, actual_b, j, depth, actualBlockSize, alpha,
    150                     -1, -1, 0, 0, workspace);
    151 
    152       // selfadjoint micro block
    153       {
    154         Index i = j;
    155         buffer.setZero();
    156         // 1 - apply the kernel on the temporary buffer
    157         gebp_kernel(buffer.data(), BlockSize, blockA+depth*i, actual_b, actualBlockSize, depth, actualBlockSize, alpha,
    158                     -1, -1, 0, 0, workspace);
    159         // 2 - triangular accumulation
    160         for(Index j1=0; j1<actualBlockSize; ++j1)
    161         {
    162           ResScalar* r = res + (j+j1)*resStride + i;
    163           for(Index i1=UpLo==Lower ? j1 : 0;
    164               UpLo==Lower ? i1<actualBlockSize : i1<=j1; ++i1)
    165             r[i1] += buffer(i1,j1);
    166         }
    167       }
    168 
    169       if(UpLo==Lower)
    170       {
    171         Index i = j+actualBlockSize;
    172         gebp_kernel(res+j*resStride+i, resStride, blockA+depth*i, actual_b, size-i, depth, actualBlockSize, alpha,
    173                     -1, -1, 0, 0, workspace);
    174       }
    175     }
    176   }
    177 };
    178 
    179 } // end namespace internal
    180 
    181 // high level API
    182 
    183 template<typename MatrixType, unsigned int UpLo>
    184 template<typename ProductDerived, typename _Lhs, typename _Rhs>
    185 TriangularView<MatrixType,UpLo>& TriangularView<MatrixType,UpLo>::assignProduct(const ProductBase<ProductDerived, _Lhs,_Rhs>& prod, const Scalar& alpha)
    186 {
    187   typedef typename internal::remove_all<typename ProductDerived::LhsNested>::type Lhs;
    188   typedef internal::blas_traits<Lhs> LhsBlasTraits;
    189   typedef typename LhsBlasTraits::DirectLinearAccessType ActualLhs;
    190   typedef typename internal::remove_all<ActualLhs>::type _ActualLhs;
    191   typename internal::add_const_on_value_type<ActualLhs>::type actualLhs = LhsBlasTraits::extract(prod.lhs());
    192 
    193   typedef typename internal::remove_all<typename ProductDerived::RhsNested>::type Rhs;
    194   typedef internal::blas_traits<Rhs> RhsBlasTraits;
    195   typedef typename RhsBlasTraits::DirectLinearAccessType ActualRhs;
    196   typedef typename internal::remove_all<ActualRhs>::type _ActualRhs;
    197   typename internal::add_const_on_value_type<ActualRhs>::type actualRhs = RhsBlasTraits::extract(prod.rhs());
    198 
    199   typename ProductDerived::Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(prod.lhs().derived()) * RhsBlasTraits::extractScalarFactor(prod.rhs().derived());
    200 
    201   internal::general_matrix_matrix_triangular_product<Index,
    202     typename Lhs::Scalar, _ActualLhs::Flags&RowMajorBit ? RowMajor : ColMajor, LhsBlasTraits::NeedToConjugate,
    203     typename Rhs::Scalar, _ActualRhs::Flags&RowMajorBit ? RowMajor : ColMajor, RhsBlasTraits::NeedToConjugate,
    204     MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor, UpLo>
    205     ::run(m_matrix.cols(), actualLhs.cols(),
    206           &actualLhs.coeffRef(0,0), actualLhs.outerStride(), &actualRhs.coeffRef(0,0), actualRhs.outerStride(),
    207           const_cast<Scalar*>(m_matrix.data()), m_matrix.outerStride(), actualAlpha);
    208 
    209   return *this;
    210 }
    211 
    212 } // end namespace Eigen
    213 
    214 #endif // EIGEN_GENERAL_MATRIX_MATRIX_TRIANGULAR_H
    215