BandTriangularSolver.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) 2011 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_BAND_TRIANGULARSOLVER_H
11 #define EIGEN_BAND_TRIANGULARSOLVER_H
12 
13 namespace Eigen {
14 namespace internal {
15 
16 /* \internal
17  * Solve Ax=b with A a band triangular matrix
18  * TODO: extend it to matrices for x abd b */
19 template <typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar, int StorageOrder>
21 
22 template <typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar>
23 struct band_solve_triangular_selector<Index, Mode, LhsScalar, ConjLhs, RhsScalar, RowMajor> {
26  enum { IsLower = (Mode & Lower) ? 1 : 0 };
27  static void run(Index size, Index k, const LhsScalar* _lhs, Index lhsStride, RhsScalar* _other) {
28  const LhsMap lhs(_lhs, size, k + 1, OuterStride<>(lhsStride));
29  RhsMap other(_other, size, 1);
30  std::conditional_t<ConjLhs, const CwiseUnaryOp<typename internal::scalar_conjugate_op<LhsScalar>, LhsMap>,
31  const LhsMap&>
32  cjLhs(lhs);
33 
34  for (int col = 0; col < other.cols(); ++col) {
35  for (int ii = 0; ii < size; ++ii) {
36  int i = IsLower ? ii : size - ii - 1;
37  int actual_k = (std::min)(k, ii);
38  int actual_start = IsLower ? k - actual_k : 1;
39 
40  if (actual_k > 0)
41  other.coeffRef(i, col) -= cjLhs.row(i)
42  .segment(actual_start, actual_k)
43  .transpose()
44  .cwiseProduct(other.col(col).segment(IsLower ? i - actual_k : i + 1, actual_k))
45  .sum();
46 
47  if ((Mode & UnitDiag) == 0) other.coeffRef(i, col) /= cjLhs(i, IsLower ? k : 0);
48  }
49  }
50  }
51 };
52 
53 template <typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar>
54 struct band_solve_triangular_selector<Index, Mode, LhsScalar, ConjLhs, RhsScalar, ColMajor> {
57  enum { IsLower = (Mode & Lower) ? 1 : 0 };
58  static void run(Index size, Index k, const LhsScalar* _lhs, Index lhsStride, RhsScalar* _other) {
59  const LhsMap lhs(_lhs, k + 1, size, OuterStride<>(lhsStride));
60  RhsMap other(_other, size, 1);
61  std::conditional_t<ConjLhs, const CwiseUnaryOp<typename internal::scalar_conjugate_op<LhsScalar>, LhsMap>,
62  const LhsMap&>
63  cjLhs(lhs);
64 
65  for (int col = 0; col < other.cols(); ++col) {
66  for (int ii = 0; ii < size; ++ii) {
67  int i = IsLower ? ii : size - ii - 1;
68  int actual_k = (std::min)(k, size - ii - 1);
69  int actual_start = IsLower ? 1 : k - actual_k;
70 
71  if ((Mode & UnitDiag) == 0) other.coeffRef(i, col) /= cjLhs(IsLower ? 0 : k, i);
72 
73  if (actual_k > 0)
74  other.col(col).segment(IsLower ? i + 1 : i - actual_k, actual_k) -=
75  other.coeff(i, col) * cjLhs.col(i).segment(actual_start, actual_k);
76  }
77  }
78  }
79 };
80 
81 } // namespace internal
82 } // namespace Eigen
83 
84 #endif // EIGEN_BAND_TRIANGULARSOLVER_H
int i
Definition: BiCGSTAB_step_by_step.cpp:9
m col(1)
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
Convenience specialization of Stride to specify only an outer stride See class Map for some examples.
Definition: Stride.h:104
#define min(a, b)
Definition: datatypes.h:22
@ UnitDiag
Definition: Constants.h:215
@ Lower
Definition: Constants.h:211
@ ColMajor
Definition: Constants.h:318
@ RowMajor
Definition: Constants.h:320
char char char int int * k
Definition: level2_impl.h:374
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
Definition: Eigen_Colamd.h:49
static void run(Index size, Index k, const LhsScalar *_lhs, Index lhsStride, RhsScalar *_other)
Definition: BandTriangularSolver.h:58
Map< Matrix< RhsScalar, Dynamic, 1 > > RhsMap
Definition: BandTriangularSolver.h:56
Map< const Matrix< LhsScalar, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > LhsMap
Definition: BandTriangularSolver.h:55
Map< const Matrix< LhsScalar, Dynamic, Dynamic, RowMajor >, 0, OuterStride<> > LhsMap
Definition: BandTriangularSolver.h:24
Map< Matrix< RhsScalar, Dynamic, 1 > > RhsMap
Definition: BandTriangularSolver.h:25
static void run(Index size, Index k, const LhsScalar *_lhs, Index lhsStride, RhsScalar *_other)
Definition: BandTriangularSolver.h:27
Definition: BandTriangularSolver.h:20