GeneralMatrixMatrixTriangular.h revision 7faaa9f3f0df9d23790277834d426c3d992ac3ba
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@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
13namespace Eigen {
14
15template<typename Scalar, typename Index, int StorageOrder, int UpLo, bool ConjLhs, bool ConjRhs>
16struct selfadjoint_rank1_update;
17
18namespace internal {
19
20/**********************************************************************
21* This file implements a general A * B product while
22* evaluating only one triangular part of the product.
23* This is more general version of self adjoint product (C += A A^T)
24* as the level 3 SYRK Blas routine.
25**********************************************************************/
26
27// forward declarations (defined at the end of this file)
28template<typename LhsScalar, typename RhsScalar, typename Index, int mr, int nr, bool ConjLhs, bool ConjRhs, int UpLo>
29struct tribb_kernel;
30
31/* Optimized matrix-matrix product evaluating only one triangular half */
32template <typename Index,
33          typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
34          typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs,
35                              int ResStorageOrder, int  UpLo, int Version = Specialized>
36struct general_matrix_matrix_triangular_product;
37
38// as usual if the result is row major => we transpose the product
39template <typename Index, typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
40                          typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs, int  UpLo, int Version>
41struct general_matrix_matrix_triangular_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,RowMajor,UpLo,Version>
42{
43  typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
44  static EIGEN_STRONG_INLINE void run(Index size, Index depth,const LhsScalar* lhs, Index lhsStride,
45                                      const RhsScalar* rhs, Index rhsStride, ResScalar* res, Index resStride, const ResScalar& alpha)
46  {
47    general_matrix_matrix_triangular_product<Index,
48        RhsScalar, RhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateRhs,
49        LhsScalar, LhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateLhs,
50        ColMajor, UpLo==Lower?Upper:Lower>
51      ::run(size,depth,rhs,rhsStride,lhs,lhsStride,res,resStride,alpha);
52  }
53};
54
55template <typename Index, typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
56                          typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs, int  UpLo, int Version>
57struct general_matrix_matrix_triangular_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,ColMajor,UpLo,Version>
58{
59  typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
60  static EIGEN_STRONG_INLINE void run(Index size, Index depth,const LhsScalar* _lhs, Index lhsStride,
61                                      const RhsScalar* _rhs, Index rhsStride, ResScalar* res, Index resStride, const ResScalar& alpha)
62  {
63    const_blas_data_mapper<LhsScalar, Index, LhsStorageOrder> lhs(_lhs,lhsStride);
64    const_blas_data_mapper<RhsScalar, Index, RhsStorageOrder> rhs(_rhs,rhsStride);
65
66    typedef gebp_traits<LhsScalar,RhsScalar> Traits;
67
68    Index kc = depth; // cache block size along the K direction
69    Index mc = size;  // cache block size along the M direction
70    Index nc = size;  // cache block size along the N direction
71    computeProductBlockingSizes<LhsScalar,RhsScalar>(kc, mc, nc);
72    // !!! mc must be a multiple of nr:
73    if(mc > Traits::nr)
74      mc = (mc/Traits::nr)*Traits::nr;
75
76    std::size_t sizeW = kc*Traits::WorkSpaceFactor;
77    std::size_t sizeB = sizeW + kc*size;
78    ei_declare_aligned_stack_constructed_variable(LhsScalar, blockA, kc*mc, 0);
79    ei_declare_aligned_stack_constructed_variable(RhsScalar, allocatedBlockB, sizeB, 0);
80    RhsScalar* blockB = allocatedBlockB + sizeW;
81
82    gemm_pack_lhs<LhsScalar, Index, Traits::mr, Traits::LhsProgress, LhsStorageOrder> pack_lhs;
83    gemm_pack_rhs<RhsScalar, Index, Traits::nr, RhsStorageOrder> pack_rhs;
84    gebp_kernel <LhsScalar, RhsScalar, Index, Traits::mr, Traits::nr, ConjugateLhs, ConjugateRhs> gebp;
85    tribb_kernel<LhsScalar, RhsScalar, Index, Traits::mr, Traits::nr, ConjugateLhs, ConjugateRhs, UpLo> sybb;
86
87    for(Index k2=0; k2<depth; k2+=kc)
88    {
89      const Index actual_kc = (std::min)(k2+kc,depth)-k2;
90
91      // note that the actual rhs is the transpose/adjoint of mat
92      pack_rhs(blockB, &rhs(k2,0), rhsStride, actual_kc, size);
93
94      for(Index i2=0; i2<size; i2+=mc)
95      {
96        const Index actual_mc = (std::min)(i2+mc,size)-i2;
97
98        pack_lhs(blockA, &lhs(i2, k2), lhsStride, actual_kc, actual_mc);
99
100        // the selected actual_mc * size panel of res is split into three different part:
101        //  1 - before the diagonal => processed with gebp or skipped
102        //  2 - the actual_mc x actual_mc symmetric block => processed with a special kernel
103        //  3 - after the diagonal => processed with gebp or skipped
104        if (UpLo==Lower)
105          gebp(res+i2, resStride, blockA, blockB, actual_mc, actual_kc, (std::min)(size,i2), alpha,
106               -1, -1, 0, 0, allocatedBlockB);
107
108        sybb(res+resStride*i2 + i2, resStride, blockA, blockB + actual_kc*i2, actual_mc, actual_kc, alpha, allocatedBlockB);
109
110        if (UpLo==Upper)
111        {
112          Index j2 = i2+actual_mc;
113          gebp(res+resStride*j2+i2, resStride, blockA, blockB+actual_kc*j2, actual_mc, actual_kc, (std::max)(Index(0), size-j2), alpha,
114               -1, -1, 0, 0, allocatedBlockB);
115        }
116      }
117    }
118  }
119};
120
121// Optimized packed Block * packed Block product kernel evaluating only one given triangular part
122// This kernel is built on top of the gebp kernel:
123// - the current destination block is processed per panel of actual_mc x BlockSize
124//   where BlockSize is set to the minimal value allowing gebp to be as fast as possible
125// - then, as usual, each panel is split into three parts along the diagonal,
126//   the sub blocks above and below the diagonal are processed as usual,
127//   while the triangular block overlapping the diagonal is evaluated into a
128//   small temporary buffer which is then accumulated into the result using a
129//   triangular traversal.
130template<typename LhsScalar, typename RhsScalar, typename Index, int mr, int nr, bool ConjLhs, bool ConjRhs, int UpLo>
131struct tribb_kernel
132{
133  typedef gebp_traits<LhsScalar,RhsScalar,ConjLhs,ConjRhs> Traits;
134  typedef typename Traits::ResScalar ResScalar;
135
136  enum {
137    BlockSize  = EIGEN_PLAIN_ENUM_MAX(mr,nr)
138  };
139  void operator()(ResScalar* res, Index resStride, const LhsScalar* blockA, const RhsScalar* blockB, Index size, Index depth, const ResScalar& alpha, RhsScalar* workspace)
140  {
141    gebp_kernel<LhsScalar, RhsScalar, Index, mr, nr, ConjLhs, ConjRhs> gebp_kernel;
142    Matrix<ResScalar,BlockSize,BlockSize,ColMajor> buffer;
143
144    // let's process the block per panel of actual_mc x BlockSize,
145    // again, each is split into three parts, etc.
146    for (Index j=0; j<size; j+=BlockSize)
147    {
148      Index actualBlockSize = std::min<Index>(BlockSize,size - j);
149      const RhsScalar* actual_b = blockB+j*depth;
150
151      if(UpLo==Upper)
152        gebp_kernel(res+j*resStride, resStride, blockA, actual_b, j, depth, actualBlockSize, alpha,
153                    -1, -1, 0, 0, workspace);
154
155      // selfadjoint micro block
156      {
157        Index i = j;
158        buffer.setZero();
159        // 1 - apply the kernel on the temporary buffer
160        gebp_kernel(buffer.data(), BlockSize, blockA+depth*i, actual_b, actualBlockSize, depth, actualBlockSize, alpha,
161                    -1, -1, 0, 0, workspace);
162        // 2 - triangular accumulation
163        for(Index j1=0; j1<actualBlockSize; ++j1)
164        {
165          ResScalar* r = res + (j+j1)*resStride + i;
166          for(Index i1=UpLo==Lower ? j1 : 0;
167              UpLo==Lower ? i1<actualBlockSize : i1<=j1; ++i1)
168            r[i1] += buffer(i1,j1);
169        }
170      }
171
172      if(UpLo==Lower)
173      {
174        Index i = j+actualBlockSize;
175        gebp_kernel(res+j*resStride+i, resStride, blockA+depth*i, actual_b, size-i, depth, actualBlockSize, alpha,
176                    -1, -1, 0, 0, workspace);
177      }
178    }
179  }
180};
181
182} // end namespace internal
183
184// high level API
185
186template<typename MatrixType, typename ProductType, int UpLo, bool IsOuterProduct>
187struct general_product_to_triangular_selector;
188
189
190template<typename MatrixType, typename ProductType, int UpLo>
191struct general_product_to_triangular_selector<MatrixType,ProductType,UpLo,true>
192{
193  static void run(MatrixType& mat, const ProductType& prod, const typename MatrixType::Scalar& alpha)
194  {
195    typedef typename MatrixType::Scalar Scalar;
196    typedef typename MatrixType::Index Index;
197
198    typedef typename internal::remove_all<typename ProductType::LhsNested>::type Lhs;
199    typedef internal::blas_traits<Lhs> LhsBlasTraits;
200    typedef typename LhsBlasTraits::DirectLinearAccessType ActualLhs;
201    typedef typename internal::remove_all<ActualLhs>::type _ActualLhs;
202    typename internal::add_const_on_value_type<ActualLhs>::type actualLhs = LhsBlasTraits::extract(prod.lhs());
203
204    typedef typename internal::remove_all<typename ProductType::RhsNested>::type Rhs;
205    typedef internal::blas_traits<Rhs> RhsBlasTraits;
206    typedef typename RhsBlasTraits::DirectLinearAccessType ActualRhs;
207    typedef typename internal::remove_all<ActualRhs>::type _ActualRhs;
208    typename internal::add_const_on_value_type<ActualRhs>::type actualRhs = RhsBlasTraits::extract(prod.rhs());
209
210    Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(prod.lhs().derived()) * RhsBlasTraits::extractScalarFactor(prod.rhs().derived());
211
212    enum {
213      StorageOrder = (internal::traits<MatrixType>::Flags&RowMajorBit) ? RowMajor : ColMajor,
214      UseLhsDirectly = _ActualLhs::InnerStrideAtCompileTime==1,
215      UseRhsDirectly = _ActualRhs::InnerStrideAtCompileTime==1
216    };
217
218    internal::gemv_static_vector_if<Scalar,Lhs::SizeAtCompileTime,Lhs::MaxSizeAtCompileTime,!UseLhsDirectly> static_lhs;
219    ei_declare_aligned_stack_constructed_variable(Scalar, actualLhsPtr, actualLhs.size(),
220      (UseLhsDirectly ? const_cast<Scalar*>(actualLhs.data()) : static_lhs.data()));
221    if(!UseLhsDirectly) Map<typename _ActualLhs::PlainObject>(actualLhsPtr, actualLhs.size()) = actualLhs;
222
223    internal::gemv_static_vector_if<Scalar,Rhs::SizeAtCompileTime,Rhs::MaxSizeAtCompileTime,!UseRhsDirectly> static_rhs;
224    ei_declare_aligned_stack_constructed_variable(Scalar, actualRhsPtr, actualRhs.size(),
225      (UseRhsDirectly ? const_cast<Scalar*>(actualRhs.data()) : static_rhs.data()));
226    if(!UseRhsDirectly) Map<typename _ActualRhs::PlainObject>(actualRhsPtr, actualRhs.size()) = actualRhs;
227
228
229    selfadjoint_rank1_update<Scalar,Index,StorageOrder,UpLo,
230                              LhsBlasTraits::NeedToConjugate && NumTraits<Scalar>::IsComplex,
231                              RhsBlasTraits::NeedToConjugate && NumTraits<Scalar>::IsComplex>
232          ::run(actualLhs.size(), mat.data(), mat.outerStride(), actualLhsPtr, actualRhsPtr, actualAlpha);
233  }
234};
235
236template<typename MatrixType, typename ProductType, int UpLo>
237struct general_product_to_triangular_selector<MatrixType,ProductType,UpLo,false>
238{
239  static void run(MatrixType& mat, const ProductType& prod, const typename MatrixType::Scalar& alpha)
240  {
241    typedef typename MatrixType::Index Index;
242
243    typedef typename internal::remove_all<typename ProductType::LhsNested>::type Lhs;
244    typedef internal::blas_traits<Lhs> LhsBlasTraits;
245    typedef typename LhsBlasTraits::DirectLinearAccessType ActualLhs;
246    typedef typename internal::remove_all<ActualLhs>::type _ActualLhs;
247    typename internal::add_const_on_value_type<ActualLhs>::type actualLhs = LhsBlasTraits::extract(prod.lhs());
248
249    typedef typename internal::remove_all<typename ProductType::RhsNested>::type Rhs;
250    typedef internal::blas_traits<Rhs> RhsBlasTraits;
251    typedef typename RhsBlasTraits::DirectLinearAccessType ActualRhs;
252    typedef typename internal::remove_all<ActualRhs>::type _ActualRhs;
253    typename internal::add_const_on_value_type<ActualRhs>::type actualRhs = RhsBlasTraits::extract(prod.rhs());
254
255    typename ProductType::Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(prod.lhs().derived()) * RhsBlasTraits::extractScalarFactor(prod.rhs().derived());
256
257    internal::general_matrix_matrix_triangular_product<Index,
258      typename Lhs::Scalar, _ActualLhs::Flags&RowMajorBit ? RowMajor : ColMajor, LhsBlasTraits::NeedToConjugate,
259      typename Rhs::Scalar, _ActualRhs::Flags&RowMajorBit ? RowMajor : ColMajor, RhsBlasTraits::NeedToConjugate,
260      MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor, UpLo>
261      ::run(mat.cols(), actualLhs.cols(),
262            &actualLhs.coeffRef(0,0), actualLhs.outerStride(), &actualRhs.coeffRef(0,0), actualRhs.outerStride(),
263            mat.data(), mat.outerStride(), actualAlpha);
264  }
265};
266
267template<typename MatrixType, unsigned int UpLo>
268template<typename ProductDerived, typename _Lhs, typename _Rhs>
269TriangularView<MatrixType,UpLo>& TriangularView<MatrixType,UpLo>::assignProduct(const ProductBase<ProductDerived, _Lhs,_Rhs>& prod, const Scalar& alpha)
270{
271  general_product_to_triangular_selector<MatrixType, ProductDerived, UpLo, (_Lhs::ColsAtCompileTime==1) || (_Rhs::RowsAtCompileTime==1)>::run(m_matrix.const_cast_derived(), prod.derived(), alpha);
272
273  return *this;
274}
275
276} // end namespace Eigen
277
278#endif // EIGEN_GENERAL_MATRIX_MATRIX_TRIANGULAR_H
279