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;
21template <
bool Conjugate,
class SparseLUType>
28 typedef typename SparseLUType::Scalar
Scalar;
40 this->m_sparseLU =
view.m_sparseLU;
43 void setSparseLU(SparseLUType* sparseLU) {m_sparseLU = sparseLU;}
45 template<
typename Rhs,
typename Dest>
48 Dest& X(X_base.derived());
51 THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
55 for(
Index j = 0; j < B.cols(); ++j){
56 X.col(j) = m_sparseLU->colsPermutation() * B.const_cast_derived().col(j);
59 m_sparseLU->matrixU().template solveTransposedInPlace<Conjugate>(X);
62 m_sparseLU->matrixL().template solveTransposedInPlace<Conjugate>(X);
65 for (
Index j = 0; j < B.cols(); ++j)
66 X.col(j) = m_sparseLU->rowsPermutation().transpose() * X.col(j);
69 inline Index rows()
const {
return m_sparseLU->rows(); }
70 inline Index cols()
const {
return m_sparseLU->cols(); }
73 SparseLUType *m_sparseLU;
130template <
typename _MatrixType,
typename _OrderingType>
141 typedef typename MatrixType::Scalar
Scalar;
205 return transposeView;
280#ifdef EIGEN_PARSED_BY_DOXYGEN
287 template<
typename Rhs>
313 template<
typename Rhs,
typename Dest>
316 Dest& X(X_base.derived());
319 THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
323 X.resize(B.rows(),B.cols());
326 for(
Index j = 0; j < B.cols(); ++j)
330 this->
matrixL().solveInPlace(X);
331 this->
matrixU().solveInPlace(X);
334 for (
Index j = 0; j < B.cols(); ++j)
360 for (
typename SCMatrix::InnerIterator it(
m_Lstore, j); it; ++it)
364 det *=
abs(it.value());
389 for (
typename SCMatrix::InnerIterator it(
m_Lstore, j); it; ++it)
391 if(it.row() < j)
continue;
394 det +=
log(
abs(it.value()));
415 for (
typename SCMatrix::InnerIterator it(
m_Lstore, j); it; ++it)
421 else if(it.value()==0)
443 for (
typename SCMatrix::InnerIterator it(
m_Lstore, j); it; ++it)
509template <
typename MatrixType,
typename OrderingType>
530 if(!mat.isCompressed())
531 IndexVector::Map(outerIndexPtr, mat.cols()+1) = IndexVector::Map(m_mat.outerIndexPtr(),mat.cols()+1);
534 for (
Index i = 0; i < mat.cols(); i++)
536 m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];
537 m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];
546 if (!m_symmetricmode) {
553 Index m = m_mat.cols();
555 for (
Index i = 0; i < m; ++i) iwork(post(i)) = post(m_etree(i));
560 for (
Index i = 0; i < m; i++)
561 post_perm.
indices()(i) = post(i);
564 if(m_perm_c.size()) {
565 m_perm_c = post_perm * m_perm_c;
570 m_analysisIsOk =
true;
594template <
typename MatrixType,
typename OrderingType>
598 eigen_assert(m_analysisIsOk &&
"analyzePattern() should be called first");
599 eigen_assert((matrix.rows() == matrix.cols()) &&
"Only for squared matrices");
601 m_isInitialized =
true;
611 if (matrix.isCompressed()) outerIndexPtr = matrix.outerIndexPtr();
615 for(
Index i = 0; i <= matrix.cols(); i++) outerIndexPtr_t[i] = m_mat.outerIndexPtr()[i];
616 outerIndexPtr = outerIndexPtr_t;
618 for (
Index i = 0; i < matrix.cols(); i++)
620 m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];
621 m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];
623 if(!matrix.isCompressed())
delete[] outerIndexPtr;
627 m_perm_c.resize(matrix.cols());
628 for(
StorageIndex i = 0; i < matrix.cols(); ++i) m_perm_c.indices()(i) = i;
631 Index m = m_mat.rows();
632 Index n = m_mat.cols();
633 Index nnz = m_mat.nonZeros();
634 Index maxpanel = m_perfv.panel_size * m;
637 Index info = Base::memInit(m, n, nnz, lwork, m_perfv.fillfactor, m_perfv.panel_size, m_glu);
640 m_lastError =
"UNABLE TO ALLOCATE WORKING MEMORY\n\n" ;
641 m_factorizationIsOk =
false;
668 if ( m_symmetricmode ==
true )
669 Base::heap_relax_snode(n, m_etree, m_perfv.relax, marker, relax_end);
671 Base::relax_snode(n, m_etree, m_perfv.relax, marker, relax_end);
675 m_perm_r.indices().setConstant(-1);
679 m_glu.supno(0) =
emptyIdxLU; m_glu.xsup.setConstant(0);
680 m_glu.xsup(0) = m_glu.xlsub(0) = m_glu.xusub(0) = m_glu.xlusup(0) =
Index(0);
691 for (jcol = 0; jcol < n; )
694 Index panel_size = m_perfv.panel_size;
695 for (k = jcol + 1; k < (
std::min)(jcol+panel_size, n); k++)
699 panel_size = k - jcol;
704 panel_size = n - jcol;
707 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);
710 Base::panel_bmod(m, panel_size, jcol, nseg1, dense, tempv, segrep, repfnz, m_glu);
713 for ( jj = jcol; jj< jcol + panel_size; jj++)
721 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);
724 m_lastError =
"UNABLE TO EXPAND MEMORY IN COLUMN_DFS() ";
726 m_factorizationIsOk =
false;
732 info = Base::column_bmod(jj, (nseg - nseg1), dense_k, tempv, segrep_k, repfnz_k, jcol, m_glu);
735 m_lastError =
"UNABLE TO EXPAND MEMORY IN COLUMN_BMOD() ";
737 m_factorizationIsOk =
false;
742 info = Base::copy_to_ucol(jj, nseg, segrep, repfnz_k ,m_perm_r.indices(), dense_k, m_glu);
745 m_lastError =
"UNABLE TO EXPAND MEMORY IN COPY_TO_UCOL() ";
747 m_factorizationIsOk =
false;
752 info = Base::pivotL(jj, m_diagpivotthresh, m_perm_r.indices(), iperm_c.
indices(), pivrow, m_glu);
755 m_lastError =
"THE MATRIX IS STRUCTURALLY SINGULAR ... ZERO COLUMN AT ";
756 std::ostringstream returnInfo;
758 m_lastError += returnInfo.str();
760 m_factorizationIsOk =
false;
766 if (pivrow != jj) m_detPermR = -m_detPermR;
769 Base::pruneL(jj, m_perm_r.indices(), pivrow, nseg, segrep, repfnz_k, xprune, m_glu);
772 for (i = 0; i < nseg; i++)
781 m_detPermR = m_perm_r.determinant();
782 m_detPermC = m_perm_c.determinant();
785 Base::countnz(n, m_nnzL, m_nnzU, m_glu);
787 Base::fixupL(n, m_perm_r.indices(), m_glu);
790 m_Lstore.setInfos(m, n, m_glu.lusup, m_glu.xlusup, m_glu.lsub, m_glu.xlsub, m_glu.supno, m_glu.xsup);
795 m_factorizationIsOk =
true;
798template<
typename MappedSupernodalType>
801 typedef typename MappedSupernodalType::Scalar
Scalar;
806 template<
typename Dest>
811 template<
bool Conjugate,
typename Dest>
814 m_mapL.template solveTransposedInPlace<Conjugate>(X);
820template<
typename MatrixLType,
typename MatrixUType>
823 typedef typename MatrixLType::Scalar
Scalar;
832 Index nrhs = X.cols();
844 for (
Index j = 0; j < nrhs; j++)
846 X(fsupc, j) /=
m_mapL.valuePtr()[luptr];
854 U = A.template triangularView<Upper>().solve(U);
857 for (
Index j = 0; j < nrhs; ++j)
859 for (
Index jcol = fsupc; jcol < fsupc + nsupc; jcol++)
861 typename MatrixUType::InnerIterator it(
m_mapU, jcol);
864 Index irow = it.index();
865 X(irow, j) -= X(jcol, j) * it.
value();
875 Index nrhs = X.cols();
885 for (
Index j = 0; j < nrhs; ++j)
887 for (
Index jcol = fsupc; jcol < fsupc + nsupc; jcol++)
889 typename MatrixUType::InnerIterator it(
m_mapU, jcol);
892 Index irow = it.index();
893 X(jcol, j) -= X(irow, j) * (
Conjugate? conj(it.value()): it.value());
899 for (
Index j = 0; j < nrhs; j++)
909 U = A.adjoint().template triangularView<Lower>().solve(U);
911 U = A.transpose().template triangularView<Lower>().solve(U);
#define eigen_assert(x)
Definition: Macros.h:1047
#define ei_declare_aligned_stack_constructed_variable(TYPE, NAME, SIZE, BUFFER)
Definition: Memory.h:768
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Definition: StaticAssert.h:127
Definition: ForwardDeclarations.h:87
EIGEN_DEVICE_FUNC CoeffReturnType value() const
Definition: DenseBase.h:526
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:50
Convenience specialization of Stride to specify only an outer stride See class Map for some examples.
Definition: Stride.h:107
InverseReturnType inverse() const
Definition: PermutationMatrix.h:185
const IndicesType & indices() const
const version of indices().
Definition: PermutationMatrix.h:360
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Resizes *this to a rows x cols matrix.
Definition: PlainObjectBase.h:271
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Resizes to the given size, and sets all coefficients in this expression to zero.
Definition: CwiseNullaryOp.h:562
EIGEN_DEVICE_FUNC Derived & setConstant(Index size, const Scalar &val)
Resizes to the given size, and sets all coefficients in this expression to the given value val.
Definition: CwiseNullaryOp.h:361
Pseudo expression representing a solving operation.
Definition: Solve.h:63
Sparse supernodal LU factorization for general matrices.
Definition: SparseLU.h:132
Scalar determinant()
Definition: SparseLU.h:434
NCMatrix m_mat
Definition: SparseLU.h:475
@ ColsAtCompileTime
Definition: SparseLU.h:152
@ MaxColsAtCompileTime
Definition: SparseLU.h:153
Scalar absDeterminant()
Definition: SparseLU.h:350
const PermutationType & colsPermutation() const
Definition: SparseLU.h:270
bool m_factorizationIsOk
Definition: SparseLU.h:472
MatrixType::Scalar Scalar
Definition: SparseLU.h:141
Index m_detPermC
Definition: SparseLU.h:490
internal::MappedSuperNodalMatrix< Scalar, StorageIndex > SCMatrix
Definition: SparseLU.h:145
PermutationType m_perm_c
Definition: SparseLU.h:478
PermutationType m_perm_r
Definition: SparseLU.h:479
const SparseLUTransposeView< false, SparseLU< _MatrixType, _OrderingType > > transpose()
Definition: SparseLU.h:200
Matrix< StorageIndex, Dynamic, 1 > IndexVector
Definition: SparseLU.h:147
_OrderingType OrderingType
Definition: SparseLU.h:140
SCMatrix m_Lstore
Definition: SparseLU.h:476
void initperfvalues()
Definition: SparseLU.h:460
void factorize(const MatrixType &matrix)
Definition: SparseLU.h:595
ComputationInfo m_info
Definition: SparseLU.h:471
internal::perfvalues m_perfv
Definition: SparseLU.h:487
void simplicialfactorize(const MatrixType &matrix)
SparseLUMatrixLReturnType< SCMatrix > matrixL() const
Definition: SparseLU.h:243
Index m_nnzL
Definition: SparseLU.h:489
std::string lastErrorMessage() const
Definition: SparseLU.h:308
PermutationMatrix< Dynamic, Dynamic, StorageIndex > PermutationType
Definition: SparseLU.h:148
MappedSparseMatrix< Scalar, ColMajor, StorageIndex > m_Ustore
Definition: SparseLU.h:477
Index rows() const
Definition: SparseLU.h:229
Index m_nnzU
Definition: SparseLU.h:489
Index cols() const
Definition: SparseLU.h:230
Scalar signDeterminant()
Definition: SparseLU.h:406
Base::GlobalLU_t m_glu
Definition: SparseLU.h:482
MatrixType::RealScalar RealScalar
Definition: SparseLU.h:142
Index m_detPermR
Definition: SparseLU.h:490
MatrixType::StorageIndex StorageIndex
Definition: SparseLU.h:143
Scalar logAbsDeterminant() const
Definition: SparseLU.h:380
SparseLU()
Definition: SparseLU.h:158
IndexVector m_etree
Definition: SparseLU.h:480
void setPivotThreshold(const RealScalar &thresh)
Set the threshold used for a diagonal entry to be an acceptable pivot.
Definition: SparseLU.h:275
void compute(const MatrixType &matrix)
Compute the symbolic and numeric factorization of the input sparse matrix.
Definition: SparseLU.h:182
bool m_analysisIsOk
Definition: SparseLU.h:473
SparseSolverBase< SparseLU< _MatrixType, _OrderingType > > APIBase
Definition: SparseLU.h:134
~SparseLU()
Definition: SparseLU.h:169
void analyzePattern(const MatrixType &matrix)
Compute the column permutation to minimize the fill-in.
Definition: SparseLU.h:510
const SparseLUTransposeView< true, SparseLU< _MatrixType, _OrderingType > > adjoint()
Definition: SparseLU.h:221
Index nnzU() const
Definition: SparseLU.h:456
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SparseLU.h:299
const PermutationType & rowsPermutation() const
Definition: SparseLU.h:262
bool m_symmetricmode
Definition: SparseLU.h:485
std::string m_lastError
Definition: SparseLU.h:474
SparseLU(const MatrixType &matrix)
Definition: SparseLU.h:162
Matrix< Scalar, Dynamic, 1 > ScalarVector
Definition: SparseLU.h:146
bool _solve_impl(const MatrixBase< Rhs > &B, MatrixBase< Dest > &X_base) const
Definition: SparseLU.h:314
Index nnzL() const
Definition: SparseLU.h:455
internal::SparseLUImpl< Scalar, StorageIndex > Base
Definition: SparseLU.h:149
RealScalar m_diagpivotthresh
Definition: SparseLU.h:488
SparseLUMatrixUReturnType< SCMatrix, MappedSparseMatrix< Scalar, ColMajor, StorageIndex > > matrixU() const
Definition: SparseLU.h:253
SparseMatrix< Scalar, ColMajor, StorageIndex > NCMatrix
Definition: SparseLU.h:144
_MatrixType MatrixType
Definition: SparseLU.h:139
void isSymmetric(bool sym)
Indicate that the pattern of the input matrix is symmetric.
Definition: SparseLU.h:232
Definition: SparseLU.h:23
SparseLUType::MatrixType MatrixType
Definition: SparseLU.h:30
Index cols() const
Definition: SparseLU.h:70
void setSparseLU(SparseLUType *sparseLU)
Definition: SparseLU.h:43
SparseLUType::Scalar Scalar
Definition: SparseLU.h:28
bool _solve_impl(const MatrixBase< Rhs > &B, MatrixBase< Dest > &X_base) const
Definition: SparseLU.h:46
SparseSolverBase< SparseLUTransposeView< Conjugate, SparseLUType > > APIBase
Definition: SparseLU.h:25
void setIsInitialized(const bool isInitialized)
Definition: SparseLU.h:42
@ MaxColsAtCompileTime
Definition: SparseLU.h:35
@ ColsAtCompileTime
Definition: SparseLU.h:34
SparseLUTransposeView()
Definition: SparseLU.h:38
Index rows() const
Definition: SparseLU.h:69
SparseLUType::OrderingType OrderingType
Definition: SparseLU.h:31
SparseLUTransposeView(const SparseLUTransposeView &view)
Definition: SparseLU.h:39
SparseLUType::StorageIndex StorageIndex
Definition: SparseLU.h:29
Index rows() const
Definition: SparseMatrix.h:138
Index cols() const
Definition: SparseMatrix.h:140
A base class for sparse solvers.
Definition: SparseSolverBase.h:68
const Solve< SparseLU< _MatrixType, _OrderingType >, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SparseSolverBase.h:88
void _solve_impl(const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
Definition: SparseSolverBase.h:111
bool m_isInitialized
Definition: SparseSolverBase.h:119
Expression of a fixed-size or dynamic-size sub-vector.
Definition: VectorBlock.h:60
Base class for sparseLU.
Definition: SparseLUImpl.h:21
Definition: XprHelper.h:110
UnitType abs(const UnitType x) noexcept
Compute absolute value.
Definition: math.h:721
dimensionless::scalar_t log(const ScalarUnit x) noexcept
Compute natural logarithm.
Definition: math.h:349
ComputationInfo
Enum for reporting the status of a computation.
Definition: Constants.h:440
@ NumericalIssue
The provided data did not satisfy the prerequisites.
Definition: Constants.h:444
@ Success
Computation was successful.
Definition: Constants.h:442
const unsigned int RowMajorBit
for a matrix, this means that the storage order is row-major.
Definition: Constants.h:66
constexpr common_t< T1, T2 > min(const T1 x, const T2 y) noexcept
Compile-time pairwise minimum function.
Definition: min.hpp:35
Index LUnumTempV(Index &m, Index &w, Index &t, Index &b)
Definition: SparseLU_Memory.h:39
@ emptyIdxLU
Definition: SparseLU_Memory.h:38
int coletree(const MatrixType &mat, IndexVector &parent, IndexVector &firstRowElt, typename MatrixType::StorageIndex *perm=0)
Compute the column elimination tree of a sparse matrix.
Definition: SparseColEtree.h:61
@ LUNoMarker
Definition: SparseLU_Memory.h:37
void treePostorder(typename IndexVector::Scalar n, IndexVector &parent, IndexVector &post)
Post order a tree.
Definition: SparseColEtree.h:178
Namespace containing all symbols from the Eigen library.
Definition: MatrixExponential.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
Definition: SparseLU.h:800
void solveTransposedInPlace(MatrixBase< Dest > &X) const
Definition: SparseLU.h:812
const MappedSupernodalType & m_mapL
Definition: SparseLU.h:817
Index rows() const
Definition: SparseLU.h:804
Index cols() const
Definition: SparseLU.h:805
SparseLUMatrixLReturnType(const MappedSupernodalType &mapL)
Definition: SparseLU.h:802
void solveInPlace(MatrixBase< Dest > &X) const
Definition: SparseLU.h:807
MappedSupernodalType::Scalar Scalar
Definition: SparseLU.h:801
Definition: SparseLU.h:822
const MatrixUType & m_mapU
Definition: SparseLU.h:918
Index cols() const
Definition: SparseLU.h:828
SparseLUMatrixUReturnType(const MatrixLType &mapL, const MatrixUType &mapU)
Definition: SparseLU.h:824
Index rows() const
Definition: SparseLU.h:827
void solveInPlace(MatrixBase< Dest > &X) const
Definition: SparseLU.h:830
const MatrixLType & m_mapL
Definition: SparseLU.h:917
void solveTransposedInPlace(MatrixBase< Dest > &X) const
Definition: SparseLU.h:872
MatrixLType::Scalar Scalar
Definition: SparseLU.h:823
Definition: SparseLU_Structs.h:77
Definition: SparseLU_Structs.h:96
Index relax
Definition: SparseLU_Structs.h:98
Index panel_size
Definition: SparseLU_Structs.h:97
Index maxsuper
Definition: SparseLU_Structs.h:101
Index rowblk
Definition: SparseLU_Structs.h:102
Index colblk
Definition: SparseLU_Structs.h:103
Index fillfactor
Definition: SparseLU_Structs.h:104