12 #ifndef EIGEN_SPARSE_LU_H
13 #define EIGEN_SPARSE_LU_H
15 #include "./InternalHeaderCheck.h"
19 template <
typename MatrixType_,
typename OrderingType_ = COLAMDOrdering<
typename MatrixType_::StorageIndex> >
class SparseLU;
20 template <
typename MappedSparseMatrixType>
struct SparseLUMatrixLReturnType;
21 template <
typename MatrixLType,
typename MatrixUType>
struct SparseLUMatrixUReturnType;
23 template <
bool Conjugate,
class SparseLUType>
24 class SparseLUTransposeView :
public SparseSolverBase<SparseLUTransposeView<Conjugate,SparseLUType> >
27 typedef SparseSolverBase<SparseLUTransposeView<Conjugate,SparseLUType> > APIBase;
28 using APIBase::m_isInitialized;
30 typedef typename SparseLUType::Scalar Scalar;
31 typedef typename SparseLUType::StorageIndex StorageIndex;
32 typedef typename SparseLUType::MatrixType MatrixType;
33 typedef typename SparseLUType::OrderingType OrderingType;
36 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
37 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
40 SparseLUTransposeView() : m_sparseLU(NULL) {}
41 SparseLUTransposeView(
const SparseLUTransposeView& view) {
42 this->m_sparseLU = view.m_sparseLU;
44 void setIsInitialized(
const bool isInitialized) {this->m_isInitialized = isInitialized;}
45 void setSparseLU(SparseLUType* sparseLU) {m_sparseLU = sparseLU;}
46 using APIBase::_solve_impl;
47 template<
typename Rhs,
typename Dest>
48 bool _solve_impl(
const MatrixBase<Rhs> &B, MatrixBase<Dest> &X_base)
const
50 Dest& X(X_base.derived());
51 eigen_assert(m_sparseLU->info() ==
Success &&
"The matrix should be factorized first");
53 THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
57 for(
Index j = 0; j < B.cols(); ++j){
58 X.col(j) = m_sparseLU->colsPermutation() * B.const_cast_derived().col(j);
61 m_sparseLU->matrixU().template solveTransposedInPlace<Conjugate>(X);
64 m_sparseLU->matrixL().template solveTransposedInPlace<Conjugate>(X);
67 for (
Index j = 0; j < B.cols(); ++j)
68 X.col(j) = m_sparseLU->rowsPermutation().transpose() * X.col(j);
71 inline Index rows()
const {
return m_sparseLU->rows(); }
72 inline Index cols()
const {
return m_sparseLU->cols(); }
75 SparseLUType *m_sparseLU;
76 SparseLUTransposeView& operator=(
const SparseLUTransposeView&);
132 template <
typename MatrixType_,
typename OrderingType_>
133 class SparseLU :
public SparseSolverBase<SparseLU<MatrixType_,OrderingType_> >,
public internal::SparseLUImpl<typename MatrixType_::Scalar, typename MatrixType_::StorageIndex>
137 using APIBase::m_isInitialized;
139 using APIBase::_solve_impl;
141 typedef MatrixType_ MatrixType;
142 typedef OrderingType_ OrderingType;
143 typedef typename MatrixType::Scalar Scalar;
144 typedef typename MatrixType::RealScalar RealScalar;
145 typedef typename MatrixType::StorageIndex StorageIndex;
147 typedef internal::MappedSuperNodalMatrix<Scalar, StorageIndex> SCMatrix;
151 typedef internal::SparseLUImpl<Scalar, StorageIndex> Base;
154 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
155 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
160 SparseLU():m_lastError(
""),m_Ustore(0,0,0,0,0,0),m_symmetricmode(
false),m_diagpivotthresh(1.0),m_detPermR(1)
164 explicit SparseLU(
const MatrixType& matrix)
165 : m_lastError(
""),m_Ustore(0,0,0,0,0,0),m_symmetricmode(
false),m_diagpivotthresh(1.0),m_detPermR(1)
177 void factorize (
const MatrixType& matrix);
178 void simplicialfactorize(
const MatrixType& matrix);
202 const SparseLUTransposeView<false,SparseLU<MatrixType_,OrderingType_> >
transpose()
204 SparseLUTransposeView<false, SparseLU<MatrixType_,OrderingType_> > transposeView;
205 transposeView.setSparseLU(
this);
206 transposeView.setIsInitialized(this->m_isInitialized);
207 return transposeView;
223 const SparseLUTransposeView<true, SparseLU<MatrixType_,OrderingType_> >
adjoint()
225 SparseLUTransposeView<true, SparseLU<MatrixType_,OrderingType_> > adjointView;
226 adjointView.setSparseLU(
this);
227 adjointView.setIsInitialized(this->m_isInitialized);
231 inline Index rows()
const {
return m_mat.
rows(); }
232 inline Index cols()
const {
return m_mat.
cols(); }
236 m_symmetricmode = sym;
245 SparseLUMatrixLReturnType<SCMatrix>
matrixL()
const
247 return SparseLUMatrixLReturnType<SCMatrix>(m_Lstore);
255 SparseLUMatrixUReturnType<SCMatrix,Map<SparseMatrix<Scalar,ColMajor,StorageIndex> > >
matrixU()
const
257 return SparseLUMatrixUReturnType<SCMatrix, Map<SparseMatrix<Scalar,ColMajor,StorageIndex> > >(m_Lstore, m_Ustore);
279 m_diagpivotthresh = thresh;
282 #ifdef EIGEN_PARSED_BY_DOXYGEN
289 template<
typename Rhs>
303 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
315 template<
typename Rhs,
typename Dest>
319 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first");
321 THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
332 this->
matrixL().solveInPlace(X);
333 this->
matrixU().solveInPlace(X);
355 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
357 Scalar det = Scalar(1.);
360 for (
Index j = 0; j < this->cols(); ++j)
362 for (
typename SCMatrix::InnerIterator it(m_Lstore, j); it; ++it)
366 det *=
abs(it.value());
387 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
388 Scalar det = Scalar(0.);
389 for (
Index j = 0; j < this->cols(); ++j)
391 for (
typename SCMatrix::InnerIterator it(m_Lstore, j); it; ++it)
393 if(it.row() < j)
continue;
396 det +=
log(
abs(it.value()));
410 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
415 for (
Index j = 0; j < this->cols(); ++j)
417 for (
typename SCMatrix::InnerIterator it(m_Lstore, j); it; ++it)
423 else if(it.value()==0)
429 return det * m_detPermR * m_detPermC;
438 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
440 Scalar det = Scalar(1.);
443 for (
Index j = 0; j < this->cols(); ++j)
445 for (
typename SCMatrix::InnerIterator it(m_Lstore, j); it; ++it)
454 return (m_detPermR * m_detPermC) > 0 ? det : -det;
457 Index nnzL()
const {
return m_nnzL; }
458 Index nnzU()
const {
return m_nnzU; }
462 void initperfvalues()
464 m_perfv.panel_size = 16;
466 m_perfv.maxsuper = 128;
469 m_perfv.fillfactor = 20;
474 bool m_factorizationIsOk;
476 std::string m_lastError;
479 Map<SparseMatrix<Scalar,ColMajor,StorageIndex>> m_Ustore;
480 PermutationType m_perm_c;
481 PermutationType m_perm_r ;
484 typename Base::GlobalLU_t m_glu;
487 bool m_symmetricmode;
489 internal::perfvalues m_perfv;
490 RealScalar m_diagpivotthresh;
491 Index m_nnzL, m_nnzU;
492 Index m_detPermR, m_detPermC;
495 SparseLU (
const SparseLU& );
511 template <
typename MatrixType,
typename OrderingType>
529 ei_declare_aligned_stack_constructed_variable(StorageIndex,outerIndexPtr,mat.cols()+1,mat.isCompressed()?
const_cast<StorageIndex*
>(mat.outerIndexPtr()):0);
532 if(!mat.isCompressed())
533 IndexVector::Map(outerIndexPtr, mat.cols()+1) = IndexVector::Map(m_mat.outerIndexPtr(),mat.cols()+1);
536 for (
Index i = 0; i < mat.cols(); i++)
538 m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];
539 m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];
545 internal::coletree(m_mat, m_etree,firstRowElt);
548 if (!m_symmetricmode) {
551 internal::treePostorder(StorageIndex(m_mat.cols()), m_etree, post);
555 Index m = m_mat.cols();
557 for (
Index i = 0; i < m; ++i) iwork(post(i)) = post(m_etree(i));
562 for (
Index i = 0; i < m; i++)
563 post_perm.
indices()(i) = post(i);
566 if(m_perm_c.size()) {
567 m_perm_c = post_perm * m_perm_c;
572 m_analysisIsOk =
true;
596 template <
typename MatrixType,
typename OrderingType>
599 using internal::emptyIdxLU;
600 eigen_assert(m_analysisIsOk &&
"analyzePattern() should be called first");
601 eigen_assert((matrix.rows() == matrix.cols()) &&
"Only for squared matrices");
603 m_isInitialized =
true;
612 const StorageIndex * outerIndexPtr;
613 if (matrix.isCompressed()) outerIndexPtr = matrix.outerIndexPtr();
616 StorageIndex* outerIndexPtr_t =
new StorageIndex[matrix.cols()+1];
617 for(
Index i = 0; i <= matrix.cols(); i++) outerIndexPtr_t[i] = m_mat.outerIndexPtr()[i];
618 outerIndexPtr = outerIndexPtr_t;
620 for (
Index i = 0; i < matrix.cols(); i++)
622 m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];
623 m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];
625 if(!matrix.isCompressed())
delete[] outerIndexPtr;
629 m_perm_c.resize(matrix.cols());
630 for(StorageIndex i = 0; i < matrix.cols(); ++i) m_perm_c.indices()(i) = i;
633 Index m = m_mat.rows();
634 Index n = m_mat.cols();
635 Index nnz = m_mat.nonZeros();
636 Index maxpanel = m_perfv.panel_size * m;
639 Index info = Base::memInit(m, n, nnz, lwork, m_perfv.fillfactor, m_perfv.panel_size, m_glu);
642 m_lastError =
"UNABLE TO ALLOCATE WORKING MEMORY\n\n" ;
643 m_factorizationIsOk =
false;
663 tempv.
setZero(internal::LUnumTempV(m, m_perfv.panel_size, m_perfv.maxsuper, m) );
670 if ( m_symmetricmode ==
true )
671 Base::heap_relax_snode(n, m_etree, m_perfv.relax, marker, relax_end);
673 Base::relax_snode(n, m_etree, m_perfv.relax, marker, relax_end);
677 m_perm_r.indices().setConstant(-1);
681 m_glu.supno(0) = emptyIdxLU; m_glu.xsup.setConstant(0);
682 m_glu.xsup(0) = m_glu.xlsub(0) = m_glu.xusub(0) = m_glu.xlusup(0) =
Index(0);
693 for (jcol = 0; jcol < n; )
696 Index panel_size = m_perfv.panel_size;
697 for (k = jcol + 1; k < (std::min)(jcol+panel_size, n); k++)
699 if (relax_end(k) != emptyIdxLU)
701 panel_size = k - jcol;
706 panel_size = n - jcol;
709 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);
712 Base::panel_bmod(m, panel_size, jcol, nseg1, dense, tempv, segrep, repfnz, m_glu);
715 for ( jj = jcol; jj< jcol + panel_size; jj++)
723 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);
726 m_lastError =
"UNABLE TO EXPAND MEMORY IN COLUMN_DFS() ";
728 m_factorizationIsOk =
false;
734 info = Base::column_bmod(jj, (nseg - nseg1), dense_k, tempv, segrep_k, repfnz_k, jcol, m_glu);
737 m_lastError =
"UNABLE TO EXPAND MEMORY IN COLUMN_BMOD() ";
739 m_factorizationIsOk =
false;
744 info = Base::copy_to_ucol(jj, nseg, segrep, repfnz_k ,m_perm_r.indices(), dense_k, m_glu);
747 m_lastError =
"UNABLE TO EXPAND MEMORY IN COPY_TO_UCOL() ";
749 m_factorizationIsOk =
false;
754 info = Base::pivotL(jj, m_diagpivotthresh, m_perm_r.indices(), iperm_c.
indices(), pivrow, m_glu);
757 m_lastError =
"THE MATRIX IS STRUCTURALLY SINGULAR ... ZERO COLUMN AT ";
758 std::ostringstream returnInfo;
760 m_lastError += returnInfo.str();
762 m_factorizationIsOk =
false;
768 if (pivrow != jj) m_detPermR = -m_detPermR;
771 Base::pruneL(jj, m_perm_r.indices(), pivrow, nseg, segrep, repfnz_k, xprune, m_glu);
774 for (i = 0; i < nseg; i++)
777 repfnz_k(irep) = emptyIdxLU;
783 m_detPermR = m_perm_r.determinant();
784 m_detPermC = m_perm_c.determinant();
787 Base::countnz(n, m_nnzL, m_nnzU, m_glu);
789 Base::fixupL(n, m_perm_r.indices(), m_glu);
792 m_Lstore.setInfos(m, n, m_glu.lusup, m_glu.xlusup, m_glu.lsub, m_glu.xlsub, m_glu.supno, m_glu.xsup);
797 m_factorizationIsOk =
true;
800 template<
typename MappedSupernodalType>
801 struct SparseLUMatrixLReturnType : internal::no_assignment_operator
803 typedef typename MappedSupernodalType::Scalar Scalar;
804 explicit SparseLUMatrixLReturnType(
const MappedSupernodalType& mapL) : m_mapL(mapL)
806 Index rows()
const {
return m_mapL.rows(); }
807 Index cols()
const {
return m_mapL.cols(); }
808 template<
typename Dest>
809 void solveInPlace( MatrixBase<Dest> &X)
const
811 m_mapL.solveInPlace(X);
813 template<
bool Conjugate,
typename Dest>
814 void solveTransposedInPlace( MatrixBase<Dest> &X)
const
816 m_mapL.template solveTransposedInPlace<Conjugate>(X);
819 const MappedSupernodalType& m_mapL;
822 template<
typename MatrixLType,
typename MatrixUType>
823 struct SparseLUMatrixUReturnType : internal::no_assignment_operator
825 typedef typename MatrixLType::Scalar Scalar;
826 SparseLUMatrixUReturnType(
const MatrixLType& mapL,
const MatrixUType& mapU)
827 : m_mapL(mapL),m_mapU(mapU)
829 Index rows()
const {
return m_mapL.rows(); }
830 Index cols()
const {
return m_mapL.cols(); }
832 template<
typename Dest>
void solveInPlace(MatrixBase<Dest> &X)
const
834 Index nrhs = X.cols();
837 for (
Index k = m_mapL.nsuper(); k >= 0; k--)
839 Index fsupc = m_mapL.supToCol()[k];
840 Index lda = m_mapL.colIndexPtr()[fsupc+1] - m_mapL.colIndexPtr()[fsupc];
841 Index nsupc = m_mapL.supToCol()[k+1] - fsupc;
842 Index luptr = m_mapL.colIndexPtr()[fsupc];
846 for (
Index j = 0; j < nrhs; j++)
848 X(fsupc, j) /= m_mapL.valuePtr()[luptr];
854 Map<const Matrix<Scalar,Dynamic,Dynamic, ColMajor>, 0, OuterStride<> > A( &(m_mapL.valuePtr()[luptr]), nsupc, nsupc, OuterStride<>(lda) );
855 Map< Matrix<Scalar,Dynamic,Dest::ColsAtCompileTime, ColMajor>, 0, OuterStride<> > U (&(X.coeffRef(fsupc,0)), nsupc, nrhs, OuterStride<>(n) );
856 U = A.template triangularView<Upper>().solve(U);
859 for (
Index j = 0; j < nrhs; ++j)
861 for (
Index jcol = fsupc; jcol < fsupc + nsupc; jcol++)
863 typename MatrixUType::InnerIterator it(m_mapU, jcol);
866 Index irow = it.index();
867 X(irow, j) -= X(jcol, j) * it.value();
874 template<
bool Conjugate,
typename Dest>
void solveTransposedInPlace(MatrixBase<Dest> &X)
const
877 Index nrhs = X.cols();
880 for (
Index k = 0; k <= m_mapL.nsuper(); k++)
882 Index fsupc = m_mapL.supToCol()[k];
883 Index lda = m_mapL.colIndexPtr()[fsupc+1] - m_mapL.colIndexPtr()[fsupc];
884 Index nsupc = m_mapL.supToCol()[k+1] - fsupc;
885 Index luptr = m_mapL.colIndexPtr()[fsupc];
887 for (
Index j = 0; j < nrhs; ++j)
889 for (
Index jcol = fsupc; jcol < fsupc + nsupc; jcol++)
891 typename MatrixUType::InnerIterator it(m_mapU, jcol);
894 Index irow = it.index();
895 X(jcol, j) -= X(irow, j) * (Conjugate?
conj(it.value()): it.value());
901 for (
Index j = 0; j < nrhs; j++)
903 X(fsupc, j) /= (Conjugate?
conj(m_mapL.valuePtr()[luptr]) : m_mapL.valuePtr()[luptr]);
908 Map<const Matrix<Scalar,Dynamic,Dynamic, ColMajor>, 0, OuterStride<> > A( &(m_mapL.valuePtr()[luptr]), nsupc, nsupc, OuterStride<>(lda) );
909 Map< Matrix<Scalar,Dynamic,Dest::ColsAtCompileTime, ColMajor>, 0, OuterStride<> > U (&(X(fsupc,0)), nsupc, nrhs, OuterStride<>(n) );
911 U = A.adjoint().template triangularView<Lower>().solve(U);
913 U = A.transpose().template triangularView<Lower>().solve(U);
919 const MatrixLType& m_mapL;
920 const MatrixUType& m_mapU;
Derived & derived()
Definition: EigenBase.h:48
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: EigenBase.h:65
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:62
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:98
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
InverseReturnType inverse() const
Definition: PermutationMatrix.h:187
const IndicesType & indices() const
Definition: PermutationMatrix.h:362
void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:283
Derived & setConstant(Index size, const Scalar &val)
Definition: CwiseNullaryOp.h:363
Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:564
Pseudo expression representing a solving operation.
Definition: Solve.h:65
Sparse supernodal LU factorization for general matrices.
Definition: SparseLU.h:134
void setPivotThreshold(const RealScalar &thresh)
Definition: SparseLU.h:277
Scalar logAbsDeterminant() const
Definition: SparseLU.h:382
SparseLUMatrixUReturnType< SCMatrix, Map< SparseMatrix< Scalar, ColMajor, StorageIndex > > > matrixU() const
Definition: SparseLU.h:255
void factorize(const MatrixType &matrix)
Definition: SparseLU.h:597
const PermutationType & rowsPermutation() const
Definition: SparseLU.h:264
std::string lastErrorMessage() const
Definition: SparseLU.h:310
const SparseLUTransposeView< true, SparseLU< MatrixType_, OrderingType_ > > adjoint()
Definition: SparseLU.h:223
SparseLUMatrixLReturnType< SCMatrix > matrixL() const
Definition: SparseLU.h:245
void compute(const MatrixType &matrix)
Definition: SparseLU.h:184
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SparseLU.h:301
Scalar signDeterminant()
Definition: SparseLU.h:408
Scalar absDeterminant()
Definition: SparseLU.h:352
void analyzePattern(const MatrixType &matrix)
Definition: SparseLU.h:512
const Solve< SparseLU, Rhs > solve(const MatrixBase< Rhs > &B) const
const SparseLUTransposeView< false, SparseLU< MatrixType_, OrderingType_ > > transpose()
Definition: SparseLU.h:202
void isSymmetric(bool sym)
Definition: SparseLU.h:234
const PermutationType & colsPermutation() const
Definition: SparseLU.h:272
Scalar determinant()
Definition: SparseLU.h:436
Index cols() const
Definition: SparseMatrix.h:142
Index rows() const
Definition: SparseMatrix.h:140
A base class for sparse solvers.
Definition: SparseSolverBase.h:70
Expression of a fixed-size or dynamic-size sub-vector.
Definition: VectorBlock.h:62
ComputationInfo
Definition: Constants.h:442
@ NumericalIssue
Definition: Constants.h:446
@ Success
Definition: Constants.h:444
const unsigned int RowMajorBit
Definition: Constants.h:68
Namespace containing all symbols from the Eigen library.
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_log_op< typename Derived::Scalar >, const Derived > log(const Eigen::ArrayBase< Derived > &x)