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 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 #include "../../../../Eigen/Sparse"
17 #include "../../../../Eigen/SparseLU"
18 
19 namespace Eigen {
20 
34 template <typename Scalar>
35 class KahanSum {
36  // Straightforward Kahan summation for accurate accumulation of a sum of numbers
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 };
50 template <typename Scalar, Index Width = 16>
51 class 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
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 
80 template <typename Derived, typename OtherDerived>
82  typedef typename Derived::Scalar Scalar;
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 
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 
127 template <typename Scalar>
129  public:
132 
134 
143 
149  slu.compute(A);
150  _result = computeInverse(slu);
151  return *this;
152  }
153 
157  const MatrixType& inverse() const { return _result; }
158 
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
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 
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:
229 };
230 
231 } // namespace Eigen
232 #endif
AnnoyingScalar conj(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:133
int i
Definition: BiCGSTAB_step_by_step.cpp:9
m col(1)
m row(1)
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition: PartialRedux_count.cpp:3
#define EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(TYPE0, TYPE1)
Definition: StaticAssert.h:60
#define EIGEN_STATIC_ASSERT_VECTOR_ONLY(TYPE)
Definition: StaticAssert.h:36
SCALAR Scalar
Definition: bench_gemm.cpp:45
Definition: SparseInverse.h:51
Index _blockUsed
Definition: SparseInverse.h:59
Scalar value()
Definition: SparseInverse.h:62
void operator+=(Scalar increment)
Definition: SparseInverse.h:64
KahanSum< Scalar > _totalSum
Definition: SparseInverse.h:57
Matrix< Scalar, Width, 1 > _block
Definition: SparseInverse.h:58
Kahan algorithm based accumulator.
Definition: SparseInverse.h:35
Scalar value()
Definition: SparseInverse.h:41
void operator+=(Scalar increment)
Definition: SparseInverse.h:43
Scalar _correction
Definition: SparseInverse.h:38
Scalar _sum
Definition: SparseInverse.h:37
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
InverseReturnType transpose() const
Definition: PermutationMatrix.h:177
calculate sparse subset of inverse of sparse matrix
Definition: SparseInverse.h:128
SparseInverse()
Definition: SparseInverse.h:133
const MatrixType & inverse() const
return the already-calculated sparse inverse, or a 0x0 matrix if it could not be computed
Definition: SparseInverse.h:157
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 & compute(const SparseMatrix< Scalar > &A)
Calculate the sparse inverse from a given sparse input.
Definition: SparseInverse.h:147
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
MatrixType _result
Definition: SparseInverse.h:228
SparseMatrix< Scalar, RowMajor > RowMatrixType
Definition: SparseInverse.h:131
SparseMatrix< Scalar, ColMajor > MatrixType
Definition: SparseInverse.h:130
Sparse supernodal LU factorization for general matrices.
Definition: SparseLU.h:151
SparseLUMatrixUReturnType< SCMatrix, Map< SparseMatrix< Scalar, ColMajor, StorageIndex > > > matrixU() const
Give the MatrixU.
Definition: SparseLU.h:284
const PermutationType & rowsPermutation() const
Give the row matrix permutation.
Definition: SparseLU.h:293
SparseLUMatrixLReturnType< SCMatrix > matrixL() const
Give the matrixL.
Definition: SparseLU.h:275
void compute(const MatrixType &matrix)
Analyze and factorize the matrix so the solver is ready to solve.
Definition: SparseLU.h:210
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SparseLU.h:326
const PermutationType & colsPermutation() const
Give the column matrix permutation.
Definition: SparseLU.h:299
Base class of any sparse matrices or sparse expressions.
Definition: SparseMatrixBase.h:30
const Derived & derived() const
Definition: SparseMatrixBase.h:144
const ConstDiagonalReturnType diagonal() const
Definition: SparseMatrix.h:757
Scalar & coeffRef(Index row, Index col)
Definition: SparseMatrix.h:275
Base::ReverseInnerIterator ReverseInnerIterator
Definition: SparseMatrix.h:139
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ Success
Definition: Constants.h:440
const char const char const char * diag
Definition: level2_impl.h:86
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
Derived::Scalar accurateDot(const SparseMatrixBase< Derived > &A, const SparseMatrixBase< OtherDerived > &other)
computes an accurate dot product on two sparse vectors
Definition: SparseInverse.h:81
Definition: CoreEvaluators.h:104
Definition: Meta.h:205
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2