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