10 #ifndef EIGEN_SPARSEINVERSE_H
11 #define EIGEN_SPARSEINVERSE_H
15 #include "../../../../Eigen/Sparse"
16 #include "../../../../Eigen/SparseLU"
33 template <
typename Scalar>
43 const Scalar correctedIncrement = increment +
_correction;
44 const Scalar previousSum =
_sum;
45 _sum += correctedIncrement;
49 template <
typename Scalar, Index W
idth = 16>
79 template <
typename Derived,
typename OtherDerived>
81 typedef typename Derived::Scalar Scalar;
85 static_assert(internal::is_same<Scalar, typename OtherDerived::Scalar>::value,
"mismatched types");
87 internal::evaluator<Derived> thisEval(
A.
derived());
88 typename Derived::ReverseInnerIterator
i(thisEval, 0);
90 internal::evaluator<OtherDerived> otherEval(other.
derived());
91 typename OtherDerived::ReverseInnerIterator
j(otherEval, 0);
95 if (
i.index() ==
j.index()) {
96 res += numext::conj(
i.value()) *
j.value();
99 }
else if (
i.index() >
j.index())
126 template <
typename Scalar>
172 invD = DU.
diagonal().cwiseInverse();
173 Upper = (invD.asDiagonal() * DU).
template triangularView<StrictlyUpper>();
188 MatrixType colInv =
Lower.transpose().template triangularView<UnitUpper>();
189 colInv +=
Upper.transpose();
193 rowInv +=
Lower.transpose();
197 for (
Index recurseLevel =
Upper.cols() - 1; recurseLevel >= 0; recurseLevel--) {
198 const auto&
col =
Lower.col(recurseLevel);
199 const auto&
row =
Upper.row(recurseLevel);
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;
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;
218 const Scalar diag = inverseDiagonal(recurseLevel) -
accurateDot(
row, colInv.col(recurseLevel));
219 rowIter.valueRef() = diag;
220 colIter.valueRef() = diag;
SparseMatrix< double > A(n, n)
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)
void operator+=(Scalar increment)
KahanSum< Scalar > _totalSum
Matrix< Scalar, Width, 1 > _block
Kahan algorithm based accumulator.
void operator+=(Scalar increment)
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
SparseMatrix< double, Options_, StorageIndex_ > & derived()