TriangularSolverVector.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) 2008-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_TRIANGULAR_SOLVER_VECTOR_H
11 #define EIGEN_TRIANGULAR_SOLVER_VECTOR_H
12 
13 #include "../InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
17 namespace internal {
18 
19 template<typename LhsScalar, typename RhsScalar, typename Index, int Mode, bool Conjugate, int StorageOrder>
20 struct triangular_solve_vector<LhsScalar, RhsScalar, Index, OnTheRight, Mode, Conjugate, StorageOrder>
21 {
22  static void run(Index size, const LhsScalar* _lhs, Index lhsStride, RhsScalar* rhs)
23  {
24  triangular_solve_vector<LhsScalar,RhsScalar,Index,OnTheLeft,
25  ((Mode&Upper)==Upper ? Lower : Upper) | (Mode&UnitDiag),
26  Conjugate,StorageOrder==RowMajor?ColMajor:RowMajor
27  >::run(size, _lhs, lhsStride, rhs);
28  }
29 };
30 
31 // forward and backward substitution, row-major, rhs is a vector
32 template<typename LhsScalar, typename RhsScalar, typename Index, int Mode, bool Conjugate>
33 struct triangular_solve_vector<LhsScalar, RhsScalar, Index, OnTheLeft, Mode, Conjugate, RowMajor>
34 {
35  enum {
36  IsLower = ((Mode&Lower)==Lower)
37  };
38  static void run(Index size, const LhsScalar* _lhs, Index lhsStride, RhsScalar* rhs)
39  {
40  typedef Map<const Matrix<LhsScalar,Dynamic,Dynamic,RowMajor>, 0, OuterStride<> > LhsMap;
41  const LhsMap lhs(_lhs,size,size,OuterStride<>(lhsStride));
42 
43  typedef const_blas_data_mapper<LhsScalar,Index,RowMajor> LhsMapper;
44  typedef const_blas_data_mapper<RhsScalar,Index,ColMajor> RhsMapper;
45 
46  std::conditional_t<
47  Conjugate,
48  const CwiseUnaryOp<typename internal::scalar_conjugate_op<LhsScalar>,LhsMap>,
49  const LhsMap&> cjLhs(lhs);
50  static const Index PanelWidth = EIGEN_TUNE_TRIANGULAR_PANEL_WIDTH;
51  for(Index pi=IsLower ? 0 : size;
52  IsLower ? pi<size : pi>0;
53  IsLower ? pi+=PanelWidth : pi-=PanelWidth)
54  {
55  Index actualPanelWidth = (std::min)(IsLower ? size - pi : pi, PanelWidth);
56 
57  Index r = IsLower ? pi : size - pi; // remaining size
58  if (r > 0)
59  {
60  // let's directly call the low level product function because:
61  // 1 - it is faster to compile
62  // 2 - it is slightly faster at runtime
63  Index startRow = IsLower ? pi : pi-actualPanelWidth;
64  Index startCol = IsLower ? 0 : pi;
65 
66  general_matrix_vector_product<Index,LhsScalar,LhsMapper,RowMajor,Conjugate,RhsScalar,RhsMapper,false>::run(
67  actualPanelWidth, r,
68  LhsMapper(&lhs.coeffRef(startRow,startCol), lhsStride),
69  RhsMapper(rhs + startCol, 1),
70  rhs + startRow, 1,
71  RhsScalar(-1));
72  }
73 
74  for(Index k=0; k<actualPanelWidth; ++k)
75  {
76  Index i = IsLower ? pi+k : pi-k-1;
77  Index s = IsLower ? pi : i+1;
78  if (k>0)
79  rhs[i] -= (cjLhs.row(i).segment(s,k).transpose().cwiseProduct(Map<const Matrix<RhsScalar,Dynamic,1> >(rhs+s,k))).sum();
80 
81  if((!(Mode & UnitDiag)) && !is_identically_zero(rhs[i]))
82  rhs[i] /= cjLhs(i,i);
83  }
84  }
85  }
86 };
87 
88 // forward and backward substitution, column-major, rhs is a vector
89 template<typename LhsScalar, typename RhsScalar, typename Index, int Mode, bool Conjugate>
90 struct triangular_solve_vector<LhsScalar, RhsScalar, Index, OnTheLeft, Mode, Conjugate, ColMajor>
91 {
92  enum {
93  IsLower = ((Mode&Lower)==Lower)
94  };
95  static void run(Index size, const LhsScalar* _lhs, Index lhsStride, RhsScalar* rhs)
96  {
97  typedef Map<const Matrix<LhsScalar,Dynamic,Dynamic,ColMajor>, 0, OuterStride<> > LhsMap;
98  const LhsMap lhs(_lhs,size,size,OuterStride<>(lhsStride));
99  typedef const_blas_data_mapper<LhsScalar,Index,ColMajor> LhsMapper;
100  typedef const_blas_data_mapper<RhsScalar,Index,ColMajor> RhsMapper;
101  std::conditional_t<Conjugate,
102  const CwiseUnaryOp<typename internal::scalar_conjugate_op<LhsScalar>,LhsMap>,
103  const LhsMap&
104  > cjLhs(lhs);
105  static const Index PanelWidth = EIGEN_TUNE_TRIANGULAR_PANEL_WIDTH;
106 
107  for(Index pi=IsLower ? 0 : size;
108  IsLower ? pi<size : pi>0;
109  IsLower ? pi+=PanelWidth : pi-=PanelWidth)
110  {
111  Index actualPanelWidth = (std::min)(IsLower ? size - pi : pi, PanelWidth);
112  Index startBlock = IsLower ? pi : pi-actualPanelWidth;
113  Index endBlock = IsLower ? pi + actualPanelWidth : 0;
114 
115  for(Index k=0; k<actualPanelWidth; ++k)
116  {
117  Index i = IsLower ? pi+k : pi-k-1;
118  if(!is_identically_zero(rhs[i]))
119  {
120  if(!(Mode & UnitDiag))
121  rhs[i] /= cjLhs.coeff(i,i);
122 
123  Index r = actualPanelWidth - k - 1; // remaining size
124  Index s = IsLower ? i+1 : i-r;
125  if (r>0)
126  Map<Matrix<RhsScalar,Dynamic,1> >(rhs+s,r) -= rhs[i] * cjLhs.col(i).segment(s,r);
127  }
128  }
129  Index r = IsLower ? size - endBlock : startBlock; // remaining size
130  if (r > 0)
131  {
132  // let's directly call the low level product function because:
133  // 1 - it is faster to compile
134  // 2 - it is slightly faster at runtime
135  general_matrix_vector_product<Index,LhsScalar,LhsMapper,ColMajor,Conjugate,RhsScalar,RhsMapper,false>::run(
136  r, actualPanelWidth,
137  LhsMapper(&lhs.coeffRef(endBlock,startBlock), lhsStride),
138  RhsMapper(rhs+startBlock, 1),
139  rhs+endBlock, 1, RhsScalar(-1));
140  }
141  }
142  }
143 };
144 
145 } // end namespace internal
146 
147 } // end namespace Eigen
148 
149 #endif // EIGEN_TRIANGULAR_SOLVER_VECTOR_H
#define EIGEN_TUNE_TRIANGULAR_PANEL_WIDTH
Definition: Settings.h:38
@ UnitDiag
Definition: Constants.h:215
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ ColMajor
Definition: Constants.h:321
@ RowMajor
Definition: Constants.h:323
@ OnTheLeft
Definition: Constants.h:334
@ OnTheRight
Definition: Constants.h:336
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
bool is_identically_zero(const Scalar &s)
Definition: Meta.h:507
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82