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