12 #ifndef EIGEN_SPARSE_LU_H 13 #define EIGEN_SPARSE_LU_H 17 template <
typename _MatrixType,
typename _OrderingType = COLAMDOrdering<
typename _MatrixType::Index> >
class SparseLU;
72 template <
typename _MatrixType,
typename _OrderingType>
77 typedef _OrderingType OrderingType;
79 typedef typename MatrixType::RealScalar RealScalar;
80 typedef typename MatrixType::Index Index;
89 SparseLU():m_isInitialized(true),m_lastError(
""),m_Ustore(0,0,0,0,0,0),m_symmetricmode(false),m_diagpivotthresh(1.0),m_detPermR(1)
93 SparseLU(
const MatrixType& matrix):m_isInitialized(true),m_lastError(
""),m_Ustore(0,0,0,0,0,0),m_symmetricmode(false),m_diagpivotthresh(1.0),m_detPermR(1)
105 void factorize (
const MatrixType& matrix);
106 void simplicialfactorize(
const MatrixType& matrix);
120 inline Index rows()
const {
return m_mat.
rows(); }
121 inline Index cols()
const {
return m_mat.
cols(); }
125 m_symmetricmode = sym;
168 m_diagpivotthresh = thresh;
177 template<
typename Rhs>
180 eigen_assert(m_factorizationIsOk &&
"SparseLU is not initialized.");
181 eigen_assert(rows()==B.rows()
182 &&
"SparseLU::solve(): invalid number of rows of the right hand side matrix B");
190 template<
typename Rhs>
193 eigen_assert(m_factorizationIsOk &&
"SparseLU is not initialized.");
194 eigen_assert(rows()==B.
rows()
195 &&
"SparseLU::solve(): invalid number of rows of the right hand side matrix B");
209 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
221 template<
typename Rhs,
typename Dest>
224 Dest& X(X_base.derived());
225 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first");
227 THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
231 X.resize(B.rows(),B.cols());
234 for(Index j = 0; j < B.cols(); ++j)
238 this->
matrixL().solveInPlace(X);
239 this->
matrixU().solveInPlace(X);
242 for (Index j = 0; j < B.cols(); ++j)
260 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
262 Scalar det = Scalar(1.);
265 for (Index j = 0; j < this->cols(); ++j)
267 for (
typename SCMatrix::InnerIterator it(m_Lstore, j); it; ++it)
272 det *= abs(it.value());
290 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
291 Scalar det = Scalar(0.);
292 for (Index j = 0; j < this->cols(); ++j)
294 for (
typename SCMatrix::InnerIterator it(m_Lstore, j); it; ++it)
296 if(it.row() < j)
continue;
299 using std::log;
using std::abs;
300 det += log(abs(it.value()));
314 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
319 for (Index j = 0; j < this->cols(); ++j)
321 for (
typename SCMatrix::InnerIterator it(m_Lstore, j); it; ++it)
327 else if(it.value()==0)
333 return det * m_detPermR * m_detPermC;
342 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
344 Scalar det = Scalar(1.);
347 for (Index j = 0; j < this->cols(); ++j)
349 for (
typename SCMatrix::InnerIterator it(m_Lstore, j); it; ++it)
358 return det * Scalar(m_detPermR * m_detPermC);
363 void initperfvalues()
365 m_perfv.panel_size = 16;
367 m_perfv.maxsuper = 128;
370 m_perfv.fillfactor = 20;
375 bool m_isInitialized;
376 bool m_factorizationIsOk;
378 std::string m_lastError;
382 PermutationType m_perm_c;
383 PermutationType m_perm_r ;
389 bool m_symmetricmode;
392 RealScalar m_diagpivotthresh;
393 Index m_nnzL, m_nnzU;
394 Index m_detPermR, m_detPermC;
414 template <
typename MatrixType,
typename OrderingType>
426 if (m_perm_c.
size()) {
429 const Index * outerIndexPtr;
430 if (mat.isCompressed()) outerIndexPtr = mat.outerIndexPtr();
433 Index *outerIndexPtr_t =
new Index[mat.cols()+1];
434 for(Index i = 0; i <= mat.cols(); i++) outerIndexPtr_t[i] = m_mat.
outerIndexPtr()[i];
435 outerIndexPtr = outerIndexPtr_t;
437 for (Index i = 0; i < mat.cols(); i++)
442 if(!mat.isCompressed())
delete[] outerIndexPtr;
449 if (!m_symmetricmode) {
456 Index m = m_mat.
cols();
458 for (Index i = 0; i < m; ++i) iwork(post(i)) = post(m_etree(i));
463 for (Index i = 0; i < m; i++)
464 post_perm.
indices()(i) = post(i);
467 if(m_perm_c.
size()) {
468 m_perm_c = post_perm * m_perm_c;
473 m_analysisIsOk =
true;
497 template <
typename MatrixType,
typename OrderingType>
500 using internal::emptyIdxLU;
501 eigen_assert(m_analysisIsOk &&
"analyzePattern() should be called first");
502 eigen_assert((matrix.rows() == matrix.cols()) &&
"Only for squared matrices");
504 typedef typename IndexVector::Scalar Index;
514 const Index * outerIndexPtr;
515 if (matrix.isCompressed()) outerIndexPtr = matrix.outerIndexPtr();
518 Index* outerIndexPtr_t =
new Index[matrix.cols()+1];
519 for(Index i = 0; i <= matrix.cols(); i++) outerIndexPtr_t[i] = m_mat.
outerIndexPtr()[i];
520 outerIndexPtr = outerIndexPtr_t;
522 for (Index i = 0; i < matrix.cols(); i++)
527 if(!matrix.isCompressed())
delete[] outerIndexPtr;
531 m_perm_c.
resize(matrix.cols());
532 for(Index i = 0; i < matrix.cols(); ++i) m_perm_c.
indices()(i) = i;
535 Index m = m_mat.
rows();
536 Index n = m_mat.
cols();
538 Index maxpanel = m_perfv.panel_size * m;
541 Index
info =
Base::memInit(m, n, nnz, lwork, m_perfv.fillfactor, m_perfv.panel_size, m_glu);
544 m_lastError =
"UNABLE TO ALLOCATE WORKING MEMORY\n\n" ;
545 m_factorizationIsOk =
false;
565 tempv.
setZero(internal::LUnumTempV(m, m_perfv.panel_size, m_perfv.maxsuper, m) );
572 if ( m_symmetricmode ==
true )
579 m_perm_r.
indices().setConstant(-1);
583 m_glu.supno(0) = emptyIdxLU; m_glu.xsup.
setConstant(0);
584 m_glu.xsup(0) = m_glu.xlsub(0) = m_glu.xusub(0) = m_glu.xlusup(0) = Index(0);
596 for (jcol = 0; jcol < n; )
599 Index panel_size = m_perfv.panel_size;
600 for (k = jcol + 1; k < (std::min)(jcol+panel_size, n); k++)
602 if (relax_end(k) != emptyIdxLU)
604 panel_size = k - jcol;
609 panel_size = n - jcol;
612 Base::panel_dfs(m, panel_size, jcol, m_mat, m_perm_r.
indices(), nseg1, dense, panel_lsub, segrep, repfnz, xprune, marker, parent, xplore, m_glu);
615 Base::panel_bmod(m, panel_size, jcol, nseg1, dense, tempv, segrep, repfnz, m_glu);
618 for ( jj = jcol; jj< jcol + panel_size; jj++)
626 info =
Base::column_dfs(m, jj, m_perm_r.
indices(), m_perfv.maxsuper, nseg, panel_lsubk, segrep, repfnz_k, xprune, marker, parent, xplore, m_glu);
629 m_lastError =
"UNABLE TO EXPAND MEMORY IN COLUMN_DFS() ";
631 m_factorizationIsOk =
false;
637 info =
Base::column_bmod(jj, (nseg - nseg1), dense_k, tempv, segrep_k, repfnz_k, jcol, m_glu);
640 m_lastError =
"UNABLE TO EXPAND MEMORY IN COLUMN_BMOD() ";
642 m_factorizationIsOk =
false;
650 m_lastError =
"UNABLE TO EXPAND MEMORY IN COPY_TO_UCOL() ";
652 m_factorizationIsOk =
false;
657 info =
Base::pivotL(jj, m_diagpivotthresh, m_perm_r.
indices(), iperm_c.indices(), pivrow, m_glu);
660 m_lastError =
"THE MATRIX IS STRUCTURALLY SINGULAR ... ZERO COLUMN AT ";
661 std::ostringstream returnInfo;
663 m_lastError += returnInfo.str();
665 m_factorizationIsOk =
false;
671 if (pivrow != jj) m_detPermR = -m_detPermR;
677 for (i = 0; i < nseg; i++)
680 repfnz_k(irep) = emptyIdxLU;
695 m_Lstore.
setInfos(m, n, m_glu.lusup, m_glu.xlusup, m_glu.lsub, m_glu.xlsub, m_glu.supno, m_glu.xsup);
700 m_factorizationIsOk =
true;
703 template<
typename MappedSupernodalType>
706 typedef typename MappedSupernodalType::Index Index;
710 Index rows() {
return m_mapL.rows(); }
711 Index cols() {
return m_mapL.cols(); }
712 template<
typename Dest>
715 m_mapL.solveInPlace(X);
717 const MappedSupernodalType& m_mapL;
720 template<
typename MatrixLType,
typename MatrixUType>
723 typedef typename MatrixLType::Index Index;
726 : m_mapL(mapL),m_mapU(mapU)
728 Index rows() {
return m_mapL.rows(); }
729 Index cols() {
return m_mapL.cols(); }
733 Index nrhs = X.cols();
736 for (Index k = m_mapL.nsuper(); k >= 0; k--)
738 Index fsupc = m_mapL.supToCol()[k];
739 Index lda = m_mapL.colIndexPtr()[fsupc+1] - m_mapL.colIndexPtr()[fsupc];
740 Index nsupc = m_mapL.supToCol()[k+1] - fsupc;
741 Index luptr = m_mapL.colIndexPtr()[fsupc];
745 for (Index j = 0; j < nrhs; j++)
747 X(fsupc, j) /= m_mapL.valuePtr()[luptr];
754 U = A.template triangularView<Upper>().
solve(U);
757 for (Index j = 0; j < nrhs; ++j)
759 for (Index jcol = fsupc; jcol < fsupc + nsupc; jcol++)
761 typename MatrixUType::InnerIterator it(m_mapU, jcol);
764 Index irow = it.index();
765 X(irow, j) -= X(jcol, j) * it.
value();
771 const MatrixLType& m_mapL;
772 const MatrixUType& m_mapU;
777 template<
typename _MatrixType,
typename Derived,
typename Rhs>
782 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
784 template<
typename Dest>
void evalTo(Dest& dst)
const 786 dec()._solve(rhs(),dst);
790 template<
typename _MatrixType,
typename Derived,
typename Rhs>
795 EIGEN_MAKE_SPARSE_SOLVE_HELPERS(Dec,Rhs)
797 template<
typename Dest>
void evalTo(Dest& dst)
const 799 this->defaultEvalTo(dst);
Definition: gtest_unittest.cc:5031
Definition: ForwardDeclarations.h:124
ColXpr col(Index i)
Definition: DenseBase.h:733
Index column_bmod(const Index jcol, const Index nseg, BlockScalarVector dense, ScalarVector &tempv, BlockIndexVector segrep, BlockIndexVector repfnz, Index fpanelc, GlobalLU_t &glu)
Performs numeric block updates (sup-col) in topological order.
Definition: SparseLU_column_bmod.h:53
void uncompress()
Turns the matrix into the uncompressed mode.
Definition: SparseMatrix.h:478
const PermutationType & rowsPermutation() const
Definition: SparseLU.h:153
EIGEN_STRONG_INLINE const Scalar * data() const
Definition: PlainObjectBase.h:212
const internal::sparse_solve_retval< SparseLU, Rhs > solve(const SparseMatrixBase< Rhs > &B) const
Definition: SparseLU.h:191
Index copy_to_ucol(const Index jcol, const Index nseg, IndexVector &segrep, BlockIndexVector repfnz, IndexVector &perm_r, BlockScalarVector dense, GlobalLU_t &glu)
Performs numeric block updates (sup-col) in topological order.
Definition: SparseLU_copy_to_ucol.h:50
Definition: SparseLU.h:19
iterative scaling algorithm to equilibrate rows and column norms in matrices
Definition: TestIMU_Common.h:87
Definition: SparseSolve.h:18
Index cols() const
Definition: SparseMatrix.h:121
void fixupL(const Index n, const IndexVector &perm_r, GlobalLU_t &glu)
Fix up the data storage lsub for L-subscripts.
Definition: SparseLU_Utils.h:52
void factorize(const MatrixType &matrix)
Definition: SparseLU.h:498
const unsigned int RowMajorBit
for a matrix, this means that the storage order is row-major.
Definition: Constants.h:53
Definition: XprHelper.h:32
Sparse supernodal LU factorization for general matrices.
Definition: SparseLU.h:17
const Index * outerIndexPtr() const
Definition: SparseMatrix.h:149
Index column_dfs(const Index m, const Index jcol, IndexVector &perm_r, Index maxsuper, Index &nseg, BlockIndexVector lsub_col, IndexVector &segrep, BlockIndexVector repfnz, IndexVector &xprune, IndexVector &marker, IndexVector &parent, IndexVector &xplore, GlobalLU_t &glu)
Performs a symbolic factorization on column jcol and decide the supernode boundary.
Definition: SparseLU_column_dfs.h:93
int coletree(const MatrixType &mat, IndexVector &parent, IndexVector &firstRowElt, typename MatrixType::Index *perm=0)
Compute the column elimination tree of a sparse matrix.
Definition: SparseColEtree.h:61
void setPivotThreshold(const RealScalar &thresh)
Set the threshold used for a diagonal entry to be an acceptable pivot.
Definition: SparseLU.h:166
Scalar absDeterminant()
Definition: SparseLU.h:258
void panel_bmod(const Index m, const Index w, const Index jcol, const Index nseg, ScalarVector &dense, ScalarVector &tempv, IndexVector &segrep, IndexVector &repfnz, GlobalLU_t &glu)
Performs numeric block updates (sup-panel) in topological order.
Definition: SparseLU_panel_bmod.h:56
void pruneL(const Index jcol, const IndexVector &perm_r, const Index pivrow, const Index nseg, const IndexVector &segrep, BlockIndexVector repfnz, IndexVector &xprune, GlobalLU_t &glu)
Prunes the L-structure.
Definition: SparseLU_pruneL.h:53
Definition: SparseLU_Structs.h:77
Index size() const
Definition: PermutationMatrix.h:114
SparseLUMatrixLReturnType< SCMatrix > matrixL() const
Definition: SparseLU.h:134
Expression of a fixed-size or dynamic-size sub-vector.
Definition: ForwardDeclarations.h:83
Index nonZeros() const
Definition: SparseMatrix.h:246
const internal::solve_retval< SparseLU, Rhs > solve(const MatrixBase< Rhs > &B) const
Definition: SparseLU.h:178
Index rows() const
Definition: SparseMatrixBase.h:159
Base class of any sparse matrices or sparse expressions.
Definition: ForwardDeclarations.h:239
void compute(const MatrixType &matrix)
Compute the symbolic and numeric factorization of the input sparse matrix.
Definition: SparseLU.h:112
void isSymmetric(bool sym)
Indicate that the pattern of the input matrix is symmetric.
Definition: SparseLU.h:123
The provided data did not satisfy the prerequisites.
Definition: Constants.h:378
EIGEN_STRONG_INLINE void resize(Index nbRows, Index nbCols)
Resizes *this to a rows x cols matrix.
Definition: PlainObjectBase.h:235
Index rows() const
Definition: SparseMatrix.h:119
const Index * innerNonZeroPtr() const
Definition: SparseMatrix.h:158
void setInfos(Index m, Index n, ScalarVector &nzval, IndexVector &nzval_colptr, IndexVector &rowind, IndexVector &rowind_colptr, IndexVector &col_to_sup, IndexVector &sup_to_col)
Set appropriate pointers for the lower triangular supernodal matrix These infos are available at the ...
Definition: SparseLU_SupernodalMatrix.h:61
Scalar determinant()
Definition: SparseLU.h:340
Derived & setZero(Index size)
Resizes to the given size, and sets all coefficients in this expression to zero.
Definition: CwiseNullaryOp.h:515
Derived & setConstant(Index size, const Scalar &value)
Resizes to the given size, and sets all coefficients in this expression to the given value...
Definition: CwiseNullaryOp.h:348
Computation was successful.
Definition: Constants.h:376
void relax_snode(const Index n, IndexVector &et, const Index relax_columns, IndexVector &descendants, IndexVector &relax_end)
Identify the initial relaxed supernodes.
Definition: SparseLU_relax_snode.h:47
Base class for sparseLU.
Definition: SparseLUImpl.h:20
Definition: BandTriangularSolver.h:13
Definition: SparseLU.h:18
Definition: SparseSolve.h:17
Transpose< PermutationBase > inverse() const
Definition: PermutationMatrix.h:201
const IndicesType & indices() const
const version of indices().
Definition: PermutationMatrix.h:387
std::string lastErrorMessage() const
Definition: SparseLU.h:216
Scalar logAbsDeterminant() const
Definition: SparseLU.h:288
void panel_dfs(const Index m, const Index w, const Index jcol, MatrixType &A, IndexVector &perm_r, Index &nseg, ScalarVector &dense, IndexVector &panel_lsub, IndexVector &segrep, IndexVector &repfnz, IndexVector &xprune, IndexVector &marker, IndexVector &parent, IndexVector &xplore, GlobalLU_t &glu)
Performs a symbolic factorization on a panel of columns [jcol, jcol+w)
Definition: SparseLU_panel_dfs.h:219
Definition: ForwardDeclarations.h:125
const PermutationType & colsPermutation() const
Definition: SparseLU.h:161
Index pivotL(const Index jcol, const RealScalar &diagpivotthresh, IndexVector &perm_r, IndexVector &iperm_c, Index &pivrow, GlobalLU_t &glu)
Performs the numerical pivotin on the current column of L, and the CDIV operation.
Definition: SparseLU_pivotL.h:60
CoeffReturnType value() const
Definition: DenseBase.h:422
void analyzePattern(const MatrixType &matrix)
Compute the column permutation to minimize the fill-in.
Definition: SparseLU.h:415
void countnz(const Index n, Index &nnzL, Index &nnzU, GlobalLU_t &glu)
Count Nonzero elements in the factors.
Definition: SparseLU_Utils.h:21
Definition: SparseLU_Structs.h:97
void resize(Index newSize)
Resizes to given size.
Definition: PermutationMatrix.h:142
ComputationInfo
Enum for reporting the status of a computation.
Definition: Constants.h:374
Index determinant() const
Definition: PermutationMatrix.h:258
Scalar signDeterminant()
Definition: SparseLU.h:312
SparseLUMatrixUReturnType< SCMatrix, MappedSparseMatrix< Scalar, ColMajor, Index > > matrixU() const
Definition: SparseLU.h:144
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Index memInit(Index m, Index n, Index annz, Index lwork, Index fillratio, Index panel_size, GlobalLU_t &glu)
Allocate various working space for the numerical factorization phase.
Definition: SparseLU_Memory.h:152
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SparseLU.h:207
double Scalar
Common scalar type.
Definition: FlexibleKalmanBase.h:48
void treePostorder(Index n, IndexVector &parent, IndexVector &post)
Post order a tree.
Definition: SparseColEtree.h:178
void heap_relax_snode(const Index n, IndexVector &et, const Index relax_columns, IndexVector &descendants, IndexVector &relax_end)
Identify the initial relaxed supernodes.
Definition: SparseLU_heap_relax_snode.h:46