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