12#ifndef EIGEN_SPARSE_LU_H
13#define EIGEN_SPARSE_LU_H
17template <
typename _MatrixType,
typename _OrderingType = COLAMDOrdering<
typename _MatrixType::StorageIndex> >
class SparseLU;
18template <
typename MappedSparseMatrixType>
struct SparseLUMatrixLReturnType;
19template <
typename MatrixLType,
typename MatrixUType>
struct SparseLUMatrixUReturnType;
73template <
typename _MatrixType,
typename _OrderingType>
78 using APIBase::m_isInitialized;
80 using APIBase::_solve_impl;
82 typedef _MatrixType MatrixType;
83 typedef _OrderingType OrderingType;
84 typedef typename MatrixType::Scalar Scalar;
85 typedef typename MatrixType::RealScalar RealScalar;
86 typedef typename MatrixType::StorageIndex StorageIndex;
95 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
96 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
100 SparseLU():m_lastError(
""),m_Ustore(0,0,0,0,0,0),m_symmetricmode(
false),m_diagpivotthresh(1.0),m_detPermR(1)
104 explicit SparseLU(
const MatrixType& matrix)
105 : m_lastError(
""),m_Ustore(0,0,0,0,0,0),m_symmetricmode(
false),m_diagpivotthresh(1.0),m_detPermR(1)
117 void factorize (
const MatrixType& matrix);
118 void simplicialfactorize(
const MatrixType& matrix);
132 inline Index rows()
const {
return m_mat.
rows(); }
133 inline Index cols()
const {
return m_mat.
cols(); }
137 m_symmetricmode = sym;
146 SparseLUMatrixLReturnType<SCMatrix>
matrixL()
const
148 return SparseLUMatrixLReturnType<SCMatrix>(m_Lstore);
156 SparseLUMatrixUReturnType<SCMatrix,MappedSparseMatrix<Scalar,ColMajor,StorageIndex> >
matrixU()
const
158 return SparseLUMatrixUReturnType<SCMatrix, MappedSparseMatrix<Scalar,ColMajor,StorageIndex> >(m_Lstore, m_Ustore);
180 m_diagpivotthresh = thresh;
183#ifdef EIGEN_PARSED_BY_DOXYGEN
190 template<
typename Rhs>
204 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
216 template<
typename Rhs,
typename Dest>
219 Dest& X(X_base.derived());
220 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first");
222 THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
226 X.resize(B.rows(),B.cols());
229 for(
Index j = 0; j < B.cols(); ++j)
233 this->
matrixL().solveInPlace(X);
234 this->
matrixU().solveInPlace(X);
237 for (
Index j = 0; j < B.cols(); ++j)
256 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
258 Scalar det = Scalar(1.);
261 for (
Index j = 0; j < this->cols(); ++j)
267 det *=
abs(it.value());
288 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
289 Scalar det = Scalar(0.);
290 for (
Index j = 0; j < this->cols(); ++j)
294 if(it.row() < j)
continue;
297 det +=
log(
abs(it.value()));
311 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
316 for (
Index j = 0; j < this->cols(); ++j)
324 else if(it.value()==0)
330 return det * m_detPermR * m_detPermC;
339 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
341 Scalar det = Scalar(1.);
344 for (
Index j = 0; j < this->cols(); ++j)
355 return (m_detPermR * m_detPermC) > 0 ? det : -det;
360 void initperfvalues()
362 m_perfv.panel_size = 16;
364 m_perfv.maxsuper = 128;
367 m_perfv.fillfactor = 20;
371 mutable ComputationInfo m_info;
372 bool m_factorizationIsOk;
374 std::string m_lastError;
377 MappedSparseMatrix<Scalar,ColMajor,StorageIndex> m_Ustore;
378 PermutationType m_perm_c;
379 PermutationType m_perm_r ;
382 typename Base::GlobalLU_t m_glu;
385 bool m_symmetricmode;
387 internal::perfvalues m_perfv;
388 RealScalar m_diagpivotthresh;
389 Index m_nnzL, m_nnzU;
390 Index m_detPermR, m_detPermC;
393 SparseLU (
const SparseLU& );
410template <
typename MatrixType,
typename OrderingType>
428 ei_declare_aligned_stack_constructed_variable(StorageIndex,outerIndexPtr,mat.cols()+1,mat.isCompressed()?
const_cast<StorageIndex*
>(mat.outerIndexPtr()):0);
431 if(!mat.isCompressed())
432 IndexVector::Map(outerIndexPtr, mat.cols()+1) = IndexVector::Map(m_mat.outerIndexPtr(),mat.cols()+1);
435 for (
Index i = 0; i < mat.cols(); i++)
437 m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];
438 m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];
443 IndexVector firstRowElt;
444 internal::coletree(m_mat, m_etree,firstRowElt);
447 if (!m_symmetricmode) {
448 IndexVector post, iwork;
450 internal::treePostorder(StorageIndex(m_mat.cols()), m_etree, post);
454 Index m = m_mat.cols();
456 for (
Index i = 0; i < m; ++i) iwork(post(i)) = post(m_etree(i));
460 PermutationType post_perm(m);
461 for (
Index i = 0; i < m; i++)
462 post_perm.
indices()(i) = post(i);
465 if(m_perm_c.size()) {
466 m_perm_c = post_perm * m_perm_c;
471 m_analysisIsOk =
true;
495template <
typename MatrixType,
typename OrderingType>
498 using internal::emptyIdxLU;
499 eigen_assert(m_analysisIsOk &&
"analyzePattern() should be called first");
500 eigen_assert((matrix.rows() == matrix.cols()) &&
"Only for squared matrices");
502 m_isInitialized =
true;
512 const StorageIndex * outerIndexPtr;
513 if (matrix.isCompressed()) outerIndexPtr = matrix.outerIndexPtr();
516 StorageIndex* outerIndexPtr_t =
new StorageIndex[matrix.cols()+1];
517 for(
Index i = 0; i <= matrix.cols(); i++) outerIndexPtr_t[i] = m_mat.outerIndexPtr()[i];
518 outerIndexPtr = outerIndexPtr_t;
520 for (
Index i = 0; i < matrix.cols(); i++)
522 m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];
523 m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];
525 if(!matrix.isCompressed())
delete[] outerIndexPtr;
529 m_perm_c.resize(matrix.cols());
530 for(StorageIndex i = 0; i < matrix.cols(); ++i) m_perm_c.indices()(i) = i;
533 Index m = m_mat.rows();
534 Index n = m_mat.cols();
535 Index nnz = m_mat.nonZeros();
536 Index maxpanel = m_perfv.panel_size * m;
542 m_lastError =
"UNABLE TO ALLOCATE WORKING MEMORY\n\n" ;
543 m_factorizationIsOk =
false;
548 IndexVector segrep(m); segrep.
setZero();
549 IndexVector parent(m); parent.
setZero();
550 IndexVector xplore(m); xplore.
setZero();
551 IndexVector repfnz(maxpanel);
552 IndexVector panel_lsub(maxpanel);
553 IndexVector xprune(n); xprune.
setZero();
554 IndexVector marker(m*internal::LUNoMarker); marker.
setZero();
563 tempv.
setZero(internal::LUnumTempV(m, m_perfv.panel_size, m_perfv.maxsuper, m) );
566 PermutationType iperm_c(m_perm_c.inverse());
569 IndexVector relax_end(n);
570 if ( m_symmetricmode ==
true )
577 m_perm_r.indices().setConstant(-1);
581 m_glu.supno(0) = emptyIdxLU; m_glu.xsup.setConstant(0);
582 m_glu.xsup(0) = m_glu.xlsub(0) = m_glu.xusub(0) = m_glu.xlusup(0) =
Index(0);
588 IndexVector panel_histo(n);
594 for (jcol = 0; jcol < n; )
597 Index panel_size = m_perfv.panel_size;
598 for (k = jcol + 1; k < (std::min)(jcol+panel_size, n); k++)
600 if (relax_end(k) != emptyIdxLU)
602 panel_size = k - jcol;
607 panel_size = n - jcol;
610 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);
613 Base::panel_bmod(m, panel_size, jcol, nseg1, dense, tempv, segrep, repfnz, m_glu);
616 for ( jj = jcol; jj< jcol + panel_size; jj++)
624 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);
627 m_lastError =
"UNABLE TO EXPAND MEMORY IN COLUMN_DFS() ";
629 m_factorizationIsOk =
false;
638 m_lastError =
"UNABLE TO EXPAND MEMORY IN COLUMN_BMOD() ";
640 m_factorizationIsOk =
false;
648 m_lastError =
"UNABLE TO EXPAND MEMORY IN COPY_TO_UCOL() ";
650 m_factorizationIsOk =
false;
658 m_lastError =
"THE MATRIX IS STRUCTURALLY SINGULAR ... ZERO COLUMN AT ";
659 std::ostringstream returnInfo;
661 m_lastError += returnInfo.str();
663 m_factorizationIsOk =
false;
669 if (pivrow != jj) m_detPermR = -m_detPermR;
672 Base::pruneL(jj, m_perm_r.indices(), pivrow, nseg, segrep, repfnz_k, xprune, m_glu);
675 for (i = 0; i < nseg; i++)
678 repfnz_k(irep) = emptyIdxLU;
684 m_detPermR = m_perm_r.determinant();
685 m_detPermC = m_perm_c.determinant();
693 m_Lstore.setInfos(m, n, m_glu.lusup, m_glu.xlusup, m_glu.lsub, m_glu.xlsub, m_glu.supno, m_glu.xsup);
698 m_factorizationIsOk =
true;
701template<
typename MappedSupernodalType>
702struct SparseLUMatrixLReturnType : internal::no_assignment_operator
704 typedef typename MappedSupernodalType::Scalar
Scalar;
705 explicit SparseLUMatrixLReturnType(
const MappedSupernodalType& mapL) : m_mapL(mapL)
707 Index rows() {
return m_mapL.rows(); }
708 Index cols() {
return m_mapL.cols(); }
709 template<
typename Dest>
710 void solveInPlace( MatrixBase<Dest> &X)
const
712 m_mapL.solveInPlace(X);
714 const MappedSupernodalType& m_mapL;
717template<
typename MatrixLType,
typename MatrixUType>
718struct SparseLUMatrixUReturnType : internal::no_assignment_operator
720 typedef typename MatrixLType::Scalar Scalar;
721 SparseLUMatrixUReturnType(
const MatrixLType& mapL,
const MatrixUType& mapU)
722 : m_mapL(mapL),m_mapU(mapU)
724 Index rows() {
return m_mapL.rows(); }
725 Index cols() {
return m_mapL.cols(); }
727 template<
typename Dest>
void solveInPlace(MatrixBase<Dest> &X)
const
729 Index nrhs = X.cols();
732 for (
Index k = m_mapL.nsuper(); k >= 0; k--)
734 Index fsupc = m_mapL.supToCol()[k];
735 Index lda = m_mapL.colIndexPtr()[fsupc+1] - m_mapL.colIndexPtr()[fsupc];
736 Index nsupc = m_mapL.supToCol()[k+1] - fsupc;
737 Index luptr = m_mapL.colIndexPtr()[fsupc];
741 for (
Index j = 0; j < nrhs; j++)
743 X(fsupc, j) /= m_mapL.valuePtr()[luptr];
748 Map<const Matrix<Scalar,Dynamic,Dynamic, ColMajor>, 0, OuterStride<> > A( &(m_mapL.valuePtr()[luptr]), nsupc, nsupc, OuterStride<>(lda) );
749 Map< Matrix<Scalar,Dynamic,Dest::ColsAtCompileTime, ColMajor>, 0, OuterStride<> > U (&(X(fsupc,0)), nsupc, nrhs, OuterStride<>(n) );
750 U = A.template triangularView<Upper>().solve(U);
753 for (
Index j = 0; j < nrhs; ++j)
755 for (
Index jcol = fsupc; jcol < fsupc + nsupc; jcol++)
757 typename MatrixUType::InnerIterator it(m_mapU, jcol);
760 Index irow = it.index();
761 X(irow, j) -= X(jcol, j) * it.value();
767 const MatrixLType& m_mapL;
768 const MatrixUType& m_mapU;
ColXpr col(Index i)
Definition DenseBase.h:839
Sparse matrix.
Definition MappedSparseMatrix.h:34
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
InverseReturnType inverse() const
Definition PermutationMatrix.h:185
Permutation matrix.
Definition PermutationMatrix.h:298
const IndicesType & indices() const
Definition PermutationMatrix.h:360
void resize(Index rows, Index cols)
Definition PlainObjectBase.h:279
Derived & setConstant(Index size, const Scalar &val)
Definition CwiseNullaryOp.h:341
Derived & setZero(Index size)
Definition CwiseNullaryOp.h:515
Pseudo expression representing a solving operation.
Definition Solve.h:63
Sparse supernodal LU factorization for general matrices.
Definition SparseLU.h:75
Scalar determinant()
Definition SparseLU.h:337
Scalar absDeterminant()
Definition SparseLU.h:253
const PermutationType & colsPermutation() const
Definition SparseLU.h:173
void factorize(const MatrixType &matrix)
Definition SparseLU.h:496
SparseLUMatrixLReturnType< SCMatrix > matrixL() const
Definition SparseLU.h:146
std::string lastErrorMessage() const
Definition SparseLU.h:211
Scalar signDeterminant()
Definition SparseLU.h:309
const Solve< SparseLU, Rhs > solve(const MatrixBase< Rhs > &B) const
Scalar logAbsDeterminant() const
Definition SparseLU.h:283
void setPivotThreshold(const RealScalar &thresh)
Definition SparseLU.h:178
void compute(const MatrixType &matrix)
Definition SparseLU.h:124
void analyzePattern(const MatrixType &matrix)
Definition SparseLU.h:411
ComputationInfo info() const
Reports whether previous computation was successful.
Definition SparseLU.h:202
const PermutationType & rowsPermutation() const
Definition SparseLU.h:165
SparseLUMatrixUReturnType< SCMatrix, MappedSparseMatrix< Scalar, ColMajor, StorageIndex > > matrixU() const
Definition SparseLU.h:156
void isSymmetric(bool sym)
Definition SparseLU.h:135
A versatible sparse matrix representation.
Definition SparseMatrix.h:98
Index rows() const
Definition SparseMatrix.h:136
Index cols() const
Definition SparseMatrix.h:138
SparseSolverBase()
Definition SparseSolverBase.h:72
Expression of a fixed-size or dynamic-size sub-vector.
Definition VectorBlock.h:60
a class to manipulate the L supernodal factor from the SparseLU factorization
Definition SparseLU_SupernodalMatrix.h:34
Definition SparseLUImpl.h:21
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 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
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
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
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
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
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:151
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
void countnz(const Index n, Index &nnzL, Index &nnzU, GlobalLU_t &glu)
Count Nonzero elements in the factors.
Definition SparseLU_Utils.h:21
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
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
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
ComputationInfo
Definition Constants.h:430
@ NumericalIssue
Definition Constants.h:434
@ Success
Definition Constants.h:432
const unsigned int RowMajorBit
Definition Constants.h:61
Namespace containing all symbols from the Eigen library.
Definition A05_PortingFrom2To3.dox:1
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:65
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_log_op< typename Derived::Scalar >, const Derived > log(const Eigen::ArrayBase< Derived > &x)