$darkmode
Eigen  5.0.1-dev
TriangularSolverVector.h
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 // IWYU pragma: private
14 #include "../InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 namespace internal {
19 
20 template <typename LhsScalar, typename RhsScalar, typename Index, int Mode, bool Conjugate, int StorageOrder>
21 struct triangular_solve_vector<LhsScalar, RhsScalar, Index, OnTheRight, Mode, Conjugate, StorageOrder> {
22  static void run(Index size, const LhsScalar* _lhs, Index lhsStride, RhsScalar* rhs) {
23  triangular_solve_vector<LhsScalar, RhsScalar, Index, OnTheLeft,
24  ((Mode & Upper) == Upper ? Lower : Upper) | (Mode & UnitDiag), Conjugate,
25  StorageOrder == RowMajor ? ColMajor : RowMajor>::run(size, _lhs, lhsStride, rhs);
26  }
27 };
28 
29 // forward and backward substitution, row-major, rhs is a vector
30 template <typename LhsScalar, typename RhsScalar, typename Index, int Mode, bool Conjugate>
31 struct triangular_solve_vector<LhsScalar, RhsScalar, Index, OnTheLeft, Mode, Conjugate, RowMajor> {
32  enum { IsLower = ((Mode & Lower) == Lower) };
33  static void run(Index size, const LhsScalar* _lhs, Index lhsStride, RhsScalar* rhs) {
34  typedef Map<const Matrix<LhsScalar, Dynamic, Dynamic, RowMajor>, 0, OuterStride<> > LhsMap;
35  const LhsMap lhs(_lhs, size, size, OuterStride<>(lhsStride));
36 
37  typedef const_blas_data_mapper<LhsScalar, Index, RowMajor> LhsMapper;
38  typedef const_blas_data_mapper<RhsScalar, Index, ColMajor> RhsMapper;
39 
40  std::conditional_t<Conjugate, const CwiseUnaryOp<typename internal::scalar_conjugate_op<LhsScalar>, LhsMap>,
41  const LhsMap&>
42  cjLhs(lhs);
43  static const Index PanelWidth = EIGEN_TUNE_TRIANGULAR_PANEL_WIDTH;
44  for (Index pi = IsLower ? 0 : size; IsLower ? pi < size : pi > 0; IsLower ? pi += PanelWidth : pi -= PanelWidth) {
45  Index actualPanelWidth = (std::min)(IsLower ? size - pi : pi, PanelWidth);
46 
47  Index r = IsLower ? pi : size - pi; // remaining size
48  if (r > 0) {
49  // let's directly call the low level product function because:
50  // 1 - it is faster to compile
51  // 2 - it is slightly faster at runtime
52  Index startRow = IsLower ? pi : pi - actualPanelWidth;
53  Index startCol = IsLower ? 0 : pi;
54 
55  general_matrix_vector_product<Index, LhsScalar, LhsMapper, RowMajor, Conjugate, RhsScalar, RhsMapper,
56  false>::run(actualPanelWidth, r,
57  LhsMapper(&lhs.coeffRef(startRow, startCol), lhsStride),
58  RhsMapper(rhs + startCol, 1), rhs + startRow, 1, RhsScalar(-1));
59  }
60 
61  for (Index k = 0; k < actualPanelWidth; ++k) {
62  Index i = IsLower ? pi + k : pi - k - 1;
63  Index s = IsLower ? pi : i + 1;
64  if (k > 0)
65  rhs[i] -= (cjLhs.row(i).segment(s, k).transpose().cwiseProduct(
66  Map<const Matrix<RhsScalar, Dynamic, 1> >(rhs + s, k)))
67  .sum();
68 
69  if ((!(Mode & UnitDiag)) && !is_identically_zero(rhs[i])) rhs[i] /= cjLhs(i, i);
70  }
71  }
72  }
73 };
74 
75 // forward and backward substitution, column-major, rhs is a vector
76 template <typename LhsScalar, typename RhsScalar, typename Index, int Mode, bool Conjugate>
77 struct triangular_solve_vector<LhsScalar, RhsScalar, Index, OnTheLeft, Mode, Conjugate, ColMajor> {
78  enum { IsLower = ((Mode & Lower) == Lower) };
79  static void run(Index size, const LhsScalar* _lhs, Index lhsStride, RhsScalar* rhs) {
80  typedef Map<const Matrix<LhsScalar, Dynamic, Dynamic, ColMajor>, 0, OuterStride<> > LhsMap;
81  const LhsMap lhs(_lhs, size, size, OuterStride<>(lhsStride));
82  typedef const_blas_data_mapper<LhsScalar, Index, ColMajor> LhsMapper;
83  typedef const_blas_data_mapper<RhsScalar, Index, ColMajor> RhsMapper;
84  std::conditional_t<Conjugate, const CwiseUnaryOp<typename internal::scalar_conjugate_op<LhsScalar>, LhsMap>,
85  const LhsMap&>
86  cjLhs(lhs);
87  static const Index PanelWidth = EIGEN_TUNE_TRIANGULAR_PANEL_WIDTH;
88 
89  for (Index pi = IsLower ? 0 : size; IsLower ? pi < size : pi > 0; IsLower ? pi += PanelWidth : pi -= PanelWidth) {
90  Index actualPanelWidth = (std::min)(IsLower ? size - pi : pi, PanelWidth);
91  Index startBlock = IsLower ? pi : pi - actualPanelWidth;
92  Index endBlock = IsLower ? pi + actualPanelWidth : 0;
93 
94  for (Index k = 0; k < actualPanelWidth; ++k) {
95  Index i = IsLower ? pi + k : pi - k - 1;
96  if (!is_identically_zero(rhs[i])) {
97  if (!(Mode & UnitDiag)) rhs[i] /= cjLhs.coeff(i, i);
98 
99  Index r = actualPanelWidth - k - 1; // remaining size
100  Index s = IsLower ? i + 1 : i - r;
101  if (r > 0) Map<Matrix<RhsScalar, Dynamic, 1> >(rhs + s, r) -= rhs[i] * cjLhs.col(i).segment(s, r);
102  }
103  }
104  Index r = IsLower ? size - endBlock : startBlock; // remaining size
105  if (r > 0) {
106  // let's directly call the low level product function because:
107  // 1 - it is faster to compile
108  // 2 - it is slightly faster at runtime
109  general_matrix_vector_product<Index, LhsScalar, LhsMapper, ColMajor, Conjugate, RhsScalar, RhsMapper,
110  false>::run(r, actualPanelWidth,
111  LhsMapper(&lhs.coeffRef(endBlock, startBlock), lhsStride),
112  RhsMapper(rhs + startBlock, 1), rhs + endBlock, 1, RhsScalar(-1));
113  }
114  }
115  }
116 };
117 
118 } // end namespace internal
119 
120 } // end namespace Eigen
121 
122 #endif // EIGEN_TRIANGULAR_SOLVER_VECTOR_H
Definition: Constants.h:318
Definition: Constants.h:333
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
Definition: Constants.h:211
Definition: Constants.h:215
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Definition: Constants.h:331
Definition: Constants.h:213
Definition: Constants.h:320