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 // 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 #include "../InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
17 template<typename Scalar, typename Index, int StorageOrder, int UpLo, bool ConjLhs, bool ConjRhs>
19 
20 namespace internal {
21 
22 
29 // forward declarations (defined at the end of this file)
30 template<typename LhsScalar, typename RhsScalar, typename Index, int mr, int nr, bool ConjLhs, bool ConjRhs, int ResInnerStride, int UpLo>
31 struct tribb_kernel;
32 
33 /* Optimized matrix-matrix product evaluating only one triangular half */
34 template <typename Index,
35  typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
36  typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs,
37  int ResStorageOrder, int ResInnerStride, int UpLo, int Version = Specialized>
38 struct general_matrix_matrix_triangular_product;
39 
40 // as usual if the result is row major => we transpose the product
41 template <typename Index, typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
42  typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs,
43  int ResInnerStride, int UpLo, int Version>
44 struct general_matrix_matrix_triangular_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,RowMajor,ResInnerStride,UpLo,Version>
45 {
47  static EIGEN_STRONG_INLINE void run(Index size, Index depth,const LhsScalar* lhs, Index lhsStride,
48  const RhsScalar* rhs, Index rhsStride, ResScalar* res, Index resIncr, Index resStride,
49  const ResScalar& alpha, level3_blocking<RhsScalar,LhsScalar>& blocking)
50  {
51  general_matrix_matrix_triangular_product<Index,
52  RhsScalar, RhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateRhs,
53  LhsScalar, LhsStorageOrder==RowMajor ? ColMajor : RowMajor, ConjugateLhs,
54  ColMajor, ResInnerStride, UpLo==Lower?Upper:Lower>
55  ::run(size,depth,rhs,rhsStride,lhs,lhsStride,res,resIncr,resStride,alpha,blocking);
56  }
57 };
58 
59 template <typename Index, typename LhsScalar, int LhsStorageOrder, bool ConjugateLhs,
60  typename RhsScalar, int RhsStorageOrder, bool ConjugateRhs,
61  int ResInnerStride, int UpLo, int Version>
62 struct general_matrix_matrix_triangular_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,RhsStorageOrder,ConjugateRhs,ColMajor,ResInnerStride,UpLo,Version>
63 {
64  typedef typename ScalarBinaryOpTraits<LhsScalar, RhsScalar>::ReturnType ResScalar;
65  static EIGEN_STRONG_INLINE void run(Index size, Index depth,const LhsScalar* _lhs, Index lhsStride,
66  const RhsScalar* _rhs, Index rhsStride,
67  ResScalar* _res, Index resIncr, Index resStride,
68  const ResScalar& alpha, level3_blocking<LhsScalar,RhsScalar>& blocking)
69  {
70  typedef gebp_traits<LhsScalar,RhsScalar> Traits;
71 
72  typedef const_blas_data_mapper<LhsScalar, Index, LhsStorageOrder> LhsMapper;
73  typedef const_blas_data_mapper<RhsScalar, Index, RhsStorageOrder> RhsMapper;
74  typedef blas_data_mapper<typename Traits::ResScalar, Index, ColMajor, Unaligned, ResInnerStride> ResMapper;
75  LhsMapper lhs(_lhs,lhsStride);
76  RhsMapper rhs(_rhs,rhsStride);
77  ResMapper res(_res, resStride, resIncr);
78 
79  Index kc = blocking.kc();
80  Index mc = (std::min)(size,blocking.mc());
81 
82  // !!! mc must be a multiple of nr:
83  if(mc > Traits::nr)
84  mc = (mc/Traits::nr)*Traits::nr;
85 
86  std::size_t sizeA = kc*mc;
87  std::size_t sizeB = kc*size;
88 
89  ei_declare_aligned_stack_constructed_variable(LhsScalar, blockA, sizeA, blocking.blockA());
90  ei_declare_aligned_stack_constructed_variable(RhsScalar, blockB, sizeB, blocking.blockB());
91 
92  gemm_pack_lhs<LhsScalar, Index, LhsMapper, Traits::mr, Traits::LhsProgress, typename Traits::LhsPacket4Packing, LhsStorageOrder> pack_lhs;
93  gemm_pack_rhs<RhsScalar, Index, RhsMapper, Traits::nr, RhsStorageOrder> pack_rhs;
94  gebp_kernel<LhsScalar, RhsScalar, Index, ResMapper, Traits::mr, Traits::nr, ConjugateLhs, ConjugateRhs> gebp;
95  tribb_kernel<LhsScalar, RhsScalar, Index, Traits::mr, Traits::nr, ConjugateLhs, ConjugateRhs, ResInnerStride, UpLo> sybb;
96 
97  for(Index k2=0; k2<depth; k2+=kc)
98  {
99  const Index actual_kc = (std::min)(k2+kc,depth)-k2;
100 
101  // note that the actual rhs is the transpose/adjoint of mat
102  pack_rhs(blockB, rhs.getSubMapper(k2,0), actual_kc, size);
103 
104  for(Index i2=0; i2<size; i2+=mc)
105  {
106  const Index actual_mc = (std::min)(i2+mc,size)-i2;
107 
108  pack_lhs(blockA, lhs.getSubMapper(i2, k2), actual_kc, actual_mc);
109 
110  // the selected actual_mc * size panel of res is split into three different part:
111  // 1 - before the diagonal => processed with gebp or skipped
112  // 2 - the actual_mc x actual_mc symmetric block => processed with a special kernel
113  // 3 - after the diagonal => processed with gebp or skipped
114  if (UpLo==Lower)
115  gebp(res.getSubMapper(i2, 0), blockA, blockB, actual_mc, actual_kc,
116  (std::min)(size,i2), alpha, -1, -1, 0, 0);
117 
118  sybb(_res+resStride*i2 + resIncr*i2, resIncr, resStride, blockA, blockB + actual_kc*i2, actual_mc, actual_kc, alpha);
119 
120  if (UpLo==Upper)
121  {
122  Index j2 = i2+actual_mc;
123  gebp(res.getSubMapper(i2, j2), blockA, blockB+actual_kc*j2, actual_mc,
124  actual_kc, (std::max)(Index(0), size-j2), alpha, -1, -1, 0, 0);
125  }
126  }
127  }
128  }
129 };
130 
131 // Optimized packed Block * packed Block product kernel evaluating only one given triangular part
132 // This kernel is built on top of the gebp kernel:
133 // - the current destination block is processed per panel of actual_mc x BlockSize
134 // where BlockSize is set to the minimal value allowing gebp to be as fast as possible
135 // - then, as usual, each panel is split into three parts along the diagonal,
136 // the sub blocks above and below the diagonal are processed as usual,
137 // while the triangular block overlapping the diagonal is evaluated into a
138 // small temporary buffer which is then accumulated into the result using a
139 // triangular traversal.
140 template<typename LhsScalar, typename RhsScalar, typename Index, int mr, int nr, bool ConjLhs, bool ConjRhs, int ResInnerStride, int UpLo>
141 struct tribb_kernel
142 {
143  typedef gebp_traits<LhsScalar,RhsScalar,ConjLhs,ConjRhs> Traits;
144  typedef typename Traits::ResScalar ResScalar;
145 
146  enum {
147  BlockSize = meta_least_common_multiple<plain_enum_max(mr, nr), plain_enum_min(mr,nr)>::ret
148  };
149  void operator()(ResScalar* _res, Index resIncr, Index resStride, const LhsScalar* blockA, const RhsScalar* blockB, Index size, Index depth, const ResScalar& alpha)
150  {
151  typedef blas_data_mapper<ResScalar, Index, ColMajor, Unaligned, ResInnerStride> ResMapper;
152  typedef blas_data_mapper<ResScalar, Index, ColMajor, Unaligned> BufferMapper;
153  ResMapper res(_res, resStride, resIncr);
154  gebp_kernel<LhsScalar, RhsScalar, Index, ResMapper, mr, nr, ConjLhs, ConjRhs> gebp_kernel1;
155  gebp_kernel<LhsScalar, RhsScalar, Index, BufferMapper, mr, nr, ConjLhs, ConjRhs> gebp_kernel2;
156 
157  Matrix<ResScalar,BlockSize,BlockSize,ColMajor> buffer((internal::constructor_without_unaligned_array_assert()));
158 
159  // let's process the block per panel of actual_mc x BlockSize,
160  // again, each is split into three parts, etc.
161  for (Index j=0; j<size; j+=BlockSize)
162  {
163  Index actualBlockSize = std::min<Index>(BlockSize,size - j);
164  const RhsScalar* actual_b = blockB+j*depth;
165 
166  if(UpLo==Upper)
167  gebp_kernel1(res.getSubMapper(0, j), blockA, actual_b, j, depth, actualBlockSize, alpha,
168  -1, -1, 0, 0);
169 
170  // selfadjoint micro block
171  {
172  Index i = j;
173  buffer.setZero();
174  // 1 - apply the kernel on the temporary buffer
175  gebp_kernel2(BufferMapper(buffer.data(), BlockSize), blockA+depth*i, actual_b, actualBlockSize, depth, actualBlockSize, alpha,
176  -1, -1, 0, 0);
177 
178  // 2 - triangular accumulation
179  for(Index j1=0; j1<actualBlockSize; ++j1)
180  {
181  typename ResMapper::LinearMapper r = res.getLinearMapper(i,j+j1);
182  for(Index i1=UpLo==Lower ? j1 : 0;
183  UpLo==Lower ? i1<actualBlockSize : i1<=j1; ++i1)
184  r(i1) += buffer(i1,j1);
185  }
186  }
187 
188  if(UpLo==Lower)
189  {
190  Index i = j+actualBlockSize;
191  gebp_kernel1(res.getSubMapper(i, j), blockA+depth*i, actual_b, size-i,
192  depth, actualBlockSize, alpha, -1, -1, 0, 0);
193  }
194  }
195  }
196 };
197 
198 } // end namespace internal
199 
200 // high level API
201 
202 template<typename MatrixType, typename ProductType, int UpLo, bool IsOuterProduct>
204 
205 
206 template<typename MatrixType, typename ProductType, int UpLo>
208 {
209  static void run(MatrixType& mat, const ProductType& prod, const typename MatrixType::Scalar& alpha, bool beta)
210  {
211  typedef typename MatrixType::Scalar Scalar;
212 
214  typedef internal::blas_traits<Lhs> LhsBlasTraits;
215  typedef typename LhsBlasTraits::DirectLinearAccessType ActualLhs;
216  typedef internal::remove_all_t<ActualLhs> ActualLhs_;
217  internal::add_const_on_value_type_t<ActualLhs> actualLhs = LhsBlasTraits::extract(prod.lhs());
218 
220  typedef internal::blas_traits<Rhs> RhsBlasTraits;
221  typedef typename RhsBlasTraits::DirectLinearAccessType ActualRhs;
222  typedef internal::remove_all_t<ActualRhs> ActualRhs_;
223  internal::add_const_on_value_type_t<ActualRhs> actualRhs = RhsBlasTraits::extract(prod.rhs());
224 
225  Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(prod.lhs().derived()) * RhsBlasTraits::extractScalarFactor(prod.rhs().derived());
226 
227  if(!beta)
228  mat.template triangularView<UpLo>().setZero();
229 
230  enum {
231  StorageOrder = (internal::traits<MatrixType>::Flags&RowMajorBit) ? RowMajor : ColMajor,
232  UseLhsDirectly = ActualLhs_::InnerStrideAtCompileTime==1,
233  UseRhsDirectly = ActualRhs_::InnerStrideAtCompileTime==1
234  };
235 
236  internal::gemv_static_vector_if<Scalar,Lhs::SizeAtCompileTime,Lhs::MaxSizeAtCompileTime,!UseLhsDirectly> static_lhs;
237  ei_declare_aligned_stack_constructed_variable(Scalar, actualLhsPtr, actualLhs.size(),
238  (UseLhsDirectly ? const_cast<Scalar*>(actualLhs.data()) : static_lhs.data()));
239  if(!UseLhsDirectly) Map<typename ActualLhs_::PlainObject>(actualLhsPtr, actualLhs.size()) = actualLhs;
240 
241  internal::gemv_static_vector_if<Scalar,Rhs::SizeAtCompileTime,Rhs::MaxSizeAtCompileTime,!UseRhsDirectly> static_rhs;
242  ei_declare_aligned_stack_constructed_variable(Scalar, actualRhsPtr, actualRhs.size(),
243  (UseRhsDirectly ? const_cast<Scalar*>(actualRhs.data()) : static_rhs.data()));
244  if(!UseRhsDirectly) Map<typename ActualRhs_::PlainObject>(actualRhsPtr, actualRhs.size()) = actualRhs;
245 
246 
247  selfadjoint_rank1_update<Scalar,Index,StorageOrder,UpLo,
248  LhsBlasTraits::NeedToConjugate && NumTraits<Scalar>::IsComplex,
249  RhsBlasTraits::NeedToConjugate && NumTraits<Scalar>::IsComplex>
250  ::run(actualLhs.size(), mat.data(), mat.outerStride(), actualLhsPtr, actualRhsPtr, actualAlpha);
251  }
252 };
253 
254 template<typename MatrixType, typename ProductType, int UpLo>
256 {
257  static void run(MatrixType& mat, const ProductType& prod, const typename MatrixType::Scalar& alpha, bool beta)
258  {
260  typedef internal::blas_traits<Lhs> LhsBlasTraits;
261  typedef typename LhsBlasTraits::DirectLinearAccessType ActualLhs;
262  typedef internal::remove_all_t<ActualLhs> ActualLhs_;
263  internal::add_const_on_value_type_t<ActualLhs> actualLhs = LhsBlasTraits::extract(prod.lhs());
264 
266  typedef internal::blas_traits<Rhs> RhsBlasTraits;
267  typedef typename RhsBlasTraits::DirectLinearAccessType ActualRhs;
268  typedef internal::remove_all_t<ActualRhs> ActualRhs_;
269  internal::add_const_on_value_type_t<ActualRhs> actualRhs = RhsBlasTraits::extract(prod.rhs());
270 
271  typename ProductType::Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(prod.lhs().derived()) * RhsBlasTraits::extractScalarFactor(prod.rhs().derived());
272 
273  if(!beta)
274  mat.template triangularView<UpLo>().setZero();
275 
276  enum {
277  IsRowMajor = (internal::traits<MatrixType>::Flags&RowMajorBit) ? 1 : 0,
278  LhsIsRowMajor = ActualLhs_::Flags&RowMajorBit ? 1 : 0,
279  RhsIsRowMajor = ActualRhs_::Flags&RowMajorBit ? 1 : 0,
280  SkipDiag = (UpLo&(UnitDiag|ZeroDiag))!=0
281  };
282 
283  Index size = mat.cols();
284  if(SkipDiag)
285  size--;
286  Index depth = actualLhs.cols();
287 
288  typedef internal::gemm_blocking_space<IsRowMajor ? RowMajor : ColMajor,typename Lhs::Scalar,typename Rhs::Scalar,
289  MatrixType::MaxColsAtCompileTime, MatrixType::MaxColsAtCompileTime, ActualRhs_::MaxColsAtCompileTime> BlockingType;
290 
291  BlockingType blocking(size, size, depth, 1, false);
292 
293  internal::general_matrix_matrix_triangular_product<Index,
294  typename Lhs::Scalar, LhsIsRowMajor ? RowMajor : ColMajor, LhsBlasTraits::NeedToConjugate,
295  typename Rhs::Scalar, RhsIsRowMajor ? RowMajor : ColMajor, RhsBlasTraits::NeedToConjugate,
296  IsRowMajor ? RowMajor : ColMajor, MatrixType::InnerStrideAtCompileTime, UpLo&(Lower|Upper)>
297  ::run(size, depth,
298  &actualLhs.coeffRef(SkipDiag&&(UpLo&Lower)==Lower ? 1 : 0,0), actualLhs.outerStride(),
299  &actualRhs.coeffRef(0,SkipDiag&&(UpLo&Upper)==Upper ? 1 : 0), actualRhs.outerStride(),
300  mat.data() + (SkipDiag ? (bool(IsRowMajor) != ((UpLo&Lower)==Lower) ? mat.innerStride() : mat.outerStride() ) : 0),
301  mat.innerStride(), mat.outerStride(), actualAlpha, blocking);
302  }
303 };
304 
305 template<typename MatrixType, unsigned int UpLo>
306 template<typename ProductType>
307 EIGEN_DEVICE_FUNC TriangularView<MatrixType,UpLo>& TriangularViewImpl<MatrixType,UpLo,Dense>::_assignProduct(const ProductType& prod, const Scalar& alpha, bool beta)
308 {
309  EIGEN_STATIC_ASSERT((UpLo&UnitDiag)==0, WRITING_TO_TRIANGULAR_PART_WITH_UNIT_DIAGONAL_IS_NOT_SUPPORTED);
310  eigen_assert(derived().nestedExpression().rows() == prod.rows() && derived().cols() == prod.cols());
311 
312  general_product_to_triangular_selector<MatrixType, ProductType, UpLo, internal::traits<ProductType>::InnerSize==1>::run(derived().nestedExpression().const_cast_derived(), prod, alpha, beta);
313 
314  return derived();
315 }
316 
317 } // end namespace Eigen
318 
319 #endif // EIGEN_GENERAL_MATRIX_MATRIX_TRIANGULAR_H
IndexedView_or_Block operator()(const RowIndices &rowIndices, const ColIndices &colIndices)
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:883
#define eigen_assert(x)
Definition: Macros.h:902
#define ei_declare_aligned_stack_constructed_variable(TYPE, NAME, SIZE, BUFFER)
Definition: Memory.h:847
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:61
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:98
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182
@ UnitDiag
Definition: Constants.h:215
@ ZeroDiag
Definition: Constants.h:217
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ ColMajor
Definition: Constants.h:321
@ RowMajor
Definition: Constants.h:323
const unsigned int RowMajorBit
Definition: Constants.h:68
bfloat16() max(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:690
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
constexpr int plain_enum_min(A a, B b)
Definition: Meta.h:516
constexpr int plain_enum_max(A a, B b)
Definition: Meta.h:524
typename remove_all< T >::type remove_all_t
Definition: Meta.h:119
typename add_const_on_value_type< T >::type add_const_on_value_type_t
Definition: Meta.h:176
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
@ Specialized
Definition: Constants.h:312
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231
Determines whether the given binary operation of two numeric types is allowed and what the scalar ret...
Definition: XprHelper.h:856
static void run(MatrixType &mat, const ProductType &prod, const typename MatrixType::Scalar &alpha, bool beta)
static void run(MatrixType &mat, const ProductType &prod, const typename MatrixType::Scalar &alpha, bool beta)
std::ptrdiff_t j