Eigen-unsupported  5.0.1-dev+284dcc12
 
Loading...
Searching...
No Matches
SparseInverse.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2022 Julian Kent <jkflying@gmail.com>
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_SPARSEINVERSE_H
11#define EIGEN_SPARSEINVERSE_H
12
13// IWYU pragma: private
14#include "./InternalHeaderCheck.h"
15
16#include "../../../../Eigen/Sparse"
17#include "../../../../Eigen/SparseLU"
18
19namespace Eigen {
20
34template <typename Scalar>
35class KahanSum {
36 // Straightforward Kahan summation for accurate accumulation of a sum of numbers
37 Scalar _sum{};
38 Scalar _correction{};
39
40 public:
41 Scalar value() { return _sum; }
42
43 void operator+=(Scalar increment) {
44 const Scalar correctedIncrement = increment + _correction;
45 const Scalar previousSum = _sum;
46 _sum += correctedIncrement;
47 _correction = correctedIncrement - (_sum - previousSum);
48 }
49};
50template <typename Scalar, Index Width = 16>
51class FABSum {
52 // https://epubs.siam.org/doi/pdf/10.1137/19M1257780
53 // Fast and Accurate Blocked Summation
54 // Uses naive summation for the fast sum, and Kahan summation for the accurate sum
55 // Theoretically SIMD sum could be changed to a tree sum which would improve accuracy
56 // over naive summation
57 KahanSum<Scalar> _totalSum;
59 Index _blockUsed{};
60
61 public:
62 Scalar value() { return _block.topRows(_blockUsed).sum() + _totalSum.value(); }
63
64 void operator+=(Scalar increment) {
65 _block(_blockUsed++, 0) = increment;
66 if (_blockUsed == Width) {
67 _totalSum += _block.sum();
68 _blockUsed = 0;
69 }
70 }
71};
72
80template <typename Derived, typename OtherDerived>
81typename Derived::Scalar accurateDot(const SparseMatrixBase<Derived>& A, const SparseMatrixBase<OtherDerived>& other) {
82 typedef typename Derived::Scalar Scalar;
83 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
84 EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
85 EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Derived, OtherDerived)
86 static_assert(internal::is_same<Scalar, typename OtherDerived::Scalar>::value, "mismatched types");
87
88 internal::evaluator<Derived> thisEval(A.derived());
89 typename Derived::ReverseInnerIterator i(thisEval, 0);
90
91 internal::evaluator<OtherDerived> otherEval(other.derived());
92 typename OtherDerived::ReverseInnerIterator j(otherEval, 0);
93
94 FABSum<Scalar> res;
95 while (i && j) {
96 if (i.index() == j.index()) {
97 res += numext::conj(i.value()) * j.value();
98 --i;
99 --j;
100 } else if (i.index() > j.index())
101 --i;
102 else
103 --j;
104 }
105 return res.value();
106}
107
127template <typename Scalar>
128class SparseInverse {
129 public:
130 typedef SparseMatrix<Scalar, ColMajor> MatrixType;
131 typedef SparseMatrix<Scalar, RowMajor> RowMatrixType;
132
133 SparseInverse() {}
134
142 SparseInverse(const SparseLU<MatrixType>& slu) { _result = computeInverse(slu); }
143
147 SparseInverse& compute(const SparseMatrix<Scalar>& A) {
149 slu.compute(A);
150 _result = computeInverse(slu);
151 return *this;
152 }
153
157 const MatrixType& inverse() const { return _result; }
158
163 static MatrixType computeInverse(const SparseLU<MatrixType>& slu) {
164 if (slu.info() != Success) {
165 return MatrixType(0, 0);
166 }
167
168 // Extract from SparseLU and decompose into L, inverse D and U terms
170 RowMatrixType Upper;
171 {
172 RowMatrixType DU = slu.matrixU().toSparse();
173 invD = DU.diagonal().cwiseInverse();
174 Upper = (invD.asDiagonal() * DU).template triangularView<StrictlyUpper>();
175 }
176 MatrixType Lower = slu.matrixL().toSparse().template triangularView<StrictlyLower>();
177
178 // Compute the inverse and reapply the permutation matrix from the LU decomposition
179 return slu.colsPermutation().transpose() * computeInverse(Upper, invD, Lower) * slu.rowsPermutation();
180 }
181
185 static MatrixType computeInverse(const RowMatrixType& Upper, const Matrix<Scalar, Dynamic, 1>& inverseDiagonal,
186 const MatrixType& Lower) {
187 // Calculate the 'minimal set', which is the nonzeros of (L+U).transpose()
188 // It could be zeroed, but we will overwrite all non-zeros anyways.
189 MatrixType colInv = Lower.transpose().template triangularView<UnitUpper>();
190 colInv += Upper.transpose();
191
192 // We also need rowmajor representation in order to do efficient row-wise dot products
193 RowMatrixType rowInv = Upper.transpose().template triangularView<UnitLower>();
194 rowInv += Lower.transpose();
195
196 // Use the Takahashi algorithm to build the supporting elements of the inverse
197 // upwards and to the left, from the bottom right element, 1 col/row at a time
198 for (Index recurseLevel = Upper.cols() - 1; recurseLevel >= 0; recurseLevel--) {
199 const auto& col = Lower.col(recurseLevel);
200 const auto& row = Upper.row(recurseLevel);
201
202 // Calculate the inverse values for the nonzeros in this column
203 typename MatrixType::ReverseInnerIterator colIter(colInv, recurseLevel);
204 for (; recurseLevel < colIter.index(); --colIter) {
205 const Scalar element = -accurateDot(col, rowInv.row(colIter.index()));
206 colIter.valueRef() = element;
207 rowInv.coeffRef(colIter.index(), recurseLevel) = element;
208 }
209
210 // Calculate the inverse values for the nonzeros in this row
211 typename RowMatrixType::ReverseInnerIterator rowIter(rowInv, recurseLevel);
212 for (; recurseLevel < rowIter.index(); --rowIter) {
213 const Scalar element = -accurateDot(row, colInv.col(rowIter.index()));
214 rowIter.valueRef() = element;
215 colInv.coeffRef(recurseLevel, rowIter.index()) = element;
216 }
217
218 // And finally the diagonal, which corresponds to both row and col iterator now
219 const Scalar diag = inverseDiagonal(recurseLevel) - accurateDot(row, colInv.col(recurseLevel));
220 rowIter.valueRef() = diag;
221 colIter.valueRef() = diag;
222 }
223
224 return colInv;
225 }
226
227 private:
228 MatrixType _result;
229};
230
231} // namespace Eigen
232#endif
Kahan algorithm based accumulator.
Definition SparseInverse.h:35
static MatrixType computeInverse(const RowMatrixType &Upper, const Matrix< Scalar, Dynamic, 1 > &inverseDiagonal, const MatrixType &Lower)
Internal function to calculate the inverse from strictly upper, diagonal and strictly lower component...
Definition SparseInverse.h:185
SparseInverse(const SparseLU< MatrixType > &slu)
This Constructor is for if you already have a factored SparseLU and would like to use it to calculate...
Definition SparseInverse.h:142
static MatrixType computeInverse(const SparseLU< MatrixType > &slu)
Internal function to calculate the sparse inverse in a functional way.
Definition SparseInverse.h:163
const MatrixType & inverse() const
return the already-calculated sparse inverse, or a 0x0 matrix if it could not be computed
Definition SparseInverse.h:157
SparseInverse & compute(const SparseMatrix< Scalar > &A)
Calculate the sparse inverse from a given sparse input.
Definition SparseInverse.h:147
SparseLUMatrixUReturnType< SCMatrix, Map< SparseMatrix< Scalar, ColMajor, StorageIndex > > > matrixU() const
SparseLUMatrixLReturnType< SCMatrix > matrixL() const
void compute(const MatrixType &matrix)
ComputationInfo info() const
const PermutationType & colsPermutation() const
const PermutationType & rowsPermutation() const
RowXpr row(Index i)
ColXpr col(Index i)
Scalar & coeffRef(Index row, Index col)
DiagonalReturnType diagonal()
Namespace containing all symbols from the Eigen library.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
Derived::Scalar accurateDot(const SparseMatrixBase< Derived > &A, const SparseMatrixBase< OtherDerived > &other)
computes an accurate dot product on two sparse vectors
Definition SparseInverse.h:81
constexpr Derived & derived()