blob: 014af24e4af1a9a6f570a5fc89b5be1a8026b6a9 [file] [log] [blame]
Googler45874d82019-08-21 12:06:47 -07001// 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
13namespace internal {
14
Rasmus Munk Larsen2434cfd2023-12-06 12:02:41 -080015/* \internal
16 * Solve Ax=b with A a band triangular matrix
17 * TODO: extend it to matrices for x abd b */
18template <typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar, int StorageOrder>
Googler45874d82019-08-21 12:06:47 -070019struct band_solve_triangular_selector;
20
Rasmus Munk Larsen2434cfd2023-12-06 12:02:41 -080021template <typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar>
22struct band_solve_triangular_selector<Index, Mode, LhsScalar, ConjLhs, RhsScalar, RowMajor> {
23 typedef Map<const Matrix<LhsScalar, Dynamic, Dynamic, RowMajor>, 0, OuterStride<> > LhsMap;
24 typedef Map<Matrix<RhsScalar, Dynamic, 1> > RhsMap;
25 enum { IsLower = (Mode & Lower) ? 1 : 0 };
26 static void run(Index size, Index k, const LhsScalar* _lhs, Index lhsStride, RhsScalar* _other) {
27 const LhsMap lhs(_lhs, size, k + 1, OuterStride<>(lhsStride));
28 RhsMap other(_other, size, 1);
29 std::conditional_t<ConjLhs, const CwiseUnaryOp<typename internal::scalar_conjugate_op<LhsScalar>, LhsMap>,
30 const LhsMap&>
31 cjLhs(lhs);
Googler45874d82019-08-21 12:06:47 -070032
Rasmus Munk Larsen2434cfd2023-12-06 12:02:41 -080033 for (int col = 0; col < other.cols(); ++col) {
34 for (int ii = 0; ii < size; ++ii) {
35 int i = IsLower ? ii : size - ii - 1;
36 int actual_k = (std::min)(k, ii);
37 int actual_start = IsLower ? k - actual_k : 1;
Googler45874d82019-08-21 12:06:47 -070038
Rasmus Munk Larsen2434cfd2023-12-06 12:02:41 -080039 if (actual_k > 0)
40 other.coeffRef(i, col) -= cjLhs.row(i)
41 .segment(actual_start, actual_k)
42 .transpose()
43 .cwiseProduct(other.col(col).segment(IsLower ? i - actual_k : i + 1, actual_k))
44 .sum();
Googler45874d82019-08-21 12:06:47 -070045
Rasmus Munk Larsen2434cfd2023-12-06 12:02:41 -080046 if ((Mode & UnitDiag) == 0) other.coeffRef(i, col) /= cjLhs(i, IsLower ? k : 0);
Googler45874d82019-08-21 12:06:47 -070047 }
48 }
49 }
50};
51
Rasmus Munk Larsen2434cfd2023-12-06 12:02:41 -080052template <typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar>
53struct band_solve_triangular_selector<Index, Mode, LhsScalar, ConjLhs, RhsScalar, ColMajor> {
54 typedef Map<const Matrix<LhsScalar, Dynamic, Dynamic, ColMajor>, 0, OuterStride<> > LhsMap;
55 typedef Map<Matrix<RhsScalar, Dynamic, 1> > RhsMap;
56 enum { IsLower = (Mode & Lower) ? 1 : 0 };
57 static void run(Index size, Index k, const LhsScalar* _lhs, Index lhsStride, RhsScalar* _other) {
58 const LhsMap lhs(_lhs, k + 1, size, OuterStride<>(lhsStride));
59 RhsMap other(_other, size, 1);
60 std::conditional_t<ConjLhs, const CwiseUnaryOp<typename internal::scalar_conjugate_op<LhsScalar>, LhsMap>,
61 const LhsMap&>
62 cjLhs(lhs);
Googler45874d82019-08-21 12:06:47 -070063
Rasmus Munk Larsen2434cfd2023-12-06 12:02:41 -080064 for (int col = 0; col < other.cols(); ++col) {
65 for (int ii = 0; ii < size; ++ii) {
66 int i = IsLower ? ii : size - ii - 1;
67 int actual_k = (std::min)(k, size - ii - 1);
68 int actual_start = IsLower ? 1 : k - actual_k;
Googler45874d82019-08-21 12:06:47 -070069
Rasmus Munk Larsen2434cfd2023-12-06 12:02:41 -080070 if ((Mode & UnitDiag) == 0) other.coeffRef(i, col) /= cjLhs(IsLower ? 0 : k, i);
71
72 if (actual_k > 0)
73 other.col(col).segment(IsLower ? i + 1 : i - actual_k, actual_k) -=
74 other.coeff(i, col) * cjLhs.col(i).segment(actual_start, actual_k);
75 }
76 }
77 }
78};
79
Rasmus Munk Larsenbf93d0d2024-04-05 18:15:52 -070080} // end namespace internal
Rasmus Munk Larsen2434cfd2023-12-06 12:02:41 -080081
82#endif // EIGEN_BAND_TRIANGULARSOLVER_H