10 #ifndef EIGEN_SPARSEINVERSE_H
11 #define EIGEN_SPARSEINVERSE_H
16 #include "../../../../Eigen/Sparse"
17 #include "../../../../Eigen/SparseLU"
34 template <
typename Scalar>
46 _sum += correctedIncrement;
50 template <
typename Scalar, Index W
idth = 16>
80 template <
typename Derived,
typename OtherDerived>
89 typename Derived::ReverseInnerIterator
i(thisEval, 0);
92 typename OtherDerived::ReverseInnerIterator
j(otherEval, 0);
96 if (
i.index() ==
j.index()) {
100 }
else if (
i.index() >
j.index())
127 template <
typename Scalar>
173 invD = DU.
diagonal().cwiseInverse();
174 Upper = (invD.asDiagonal() * DU).
template triangularView<StrictlyUpper>();
189 MatrixType colInv =
Lower.transpose().template triangularView<UnitUpper>();
190 colInv +=
Upper.transpose();
194 rowInv +=
Lower.transpose();
198 for (
Index recurseLevel =
Upper.cols() - 1; recurseLevel >= 0; recurseLevel--) {
199 const auto&
col =
Lower.col(recurseLevel);
200 const auto&
row =
Upper.row(recurseLevel);
204 for (; recurseLevel < colIter.index(); --colIter) {
206 colIter.valueRef() = element;
207 rowInv.
coeffRef(colIter.index(), recurseLevel) = element;
212 for (; recurseLevel < rowIter.index(); --rowIter) {
214 rowIter.valueRef() = element;
215 colInv.
coeffRef(recurseLevel, rowIter.index()) = element;
220 rowIter.valueRef() =
diag;
221 colIter.valueRef() =
diag;
AnnoyingScalar conj(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:133
int i
Definition: BiCGSTAB_step_by_step.cpp:9
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
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2