mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-09-12 17:33:15 +08:00
SparseLU: make COLAMDOrdering the default ordering method.
This commit is contained in:
parent
bd689ccc28
commit
bbaef8ebba
@ -14,9 +14,10 @@
|
||||
|
||||
namespace Eigen {
|
||||
|
||||
template <typename _MatrixType, typename _OrderingType> class SparseLU;
|
||||
template <typename _MatrixType, typename _OrderingType = COLAMDOrdering<typename _MatrixType::Index> > class SparseLU;
|
||||
template <typename MappedSparseMatrixType> struct SparseLUMatrixLReturnType;
|
||||
template <typename MatrixLType, typename MatrixUType> struct SparseLUMatrixUReturnType;
|
||||
|
||||
/** \ingroup SparseLU_Module
|
||||
* \class SparseLU
|
||||
*
|
||||
@ -62,7 +63,7 @@ template <typename MatrixLType, typename MatrixUType> struct SparseLUMatrixURetu
|
||||
* "unsupported/Eigen/src/IterativeSolvers/Scaling.h"
|
||||
*
|
||||
* \tparam _MatrixType The type of the sparse matrix. It must be a column-major SparseMatrix<>
|
||||
* \tparam _OrderingType The ordering method to use, either AMD, COLAMD or METIS
|
||||
* \tparam _OrderingType The ordering method to use, either AMD, COLAMD or METIS. Default is COLMAD
|
||||
*
|
||||
*
|
||||
* \sa \ref TutorialSparseDirectSolvers
|
||||
@ -208,6 +209,7 @@ class SparseLU : public internal::SparseLUImpl<typename _MatrixType::Scalar, typ
|
||||
eigen_assert(m_isInitialized && "Decomposition is not initialized.");
|
||||
return m_info;
|
||||
}
|
||||
|
||||
/**
|
||||
* \returns A string describing the type of error
|
||||
*/
|
||||
@ -240,6 +242,7 @@ class SparseLU : public internal::SparseLUImpl<typename _MatrixType::Scalar, typ
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
/**
|
||||
* \returns the absolute value of the determinant of the matrix of which
|
||||
* *this is the QR decomposition.
|
||||
@ -276,8 +279,8 @@ class SparseLU : public internal::SparseLUImpl<typename _MatrixType::Scalar, typ
|
||||
* of which **this is the QR decomposition
|
||||
*
|
||||
* \note This method is useful to work around the risk of overflow/underflow that's
|
||||
* inherent to the determinant computation
|
||||
*a
|
||||
* inherent to the determinant computation.
|
||||
*
|
||||
* \sa absDeterminant(), signDeterminant()
|
||||
*/
|
||||
Scalar logAbsDeterminant() const
|
||||
|
@ -26,7 +26,7 @@
|
||||
// SparseLU solve does not accept column major matrices for the destination.
|
||||
// However, as expected, the generic check_sparse_square_solving routines produces row-major
|
||||
// rhs and destination matrices when compiled with EIGEN_DEFAULT_TO_ROW_MAJOR
|
||||
//
|
||||
|
||||
#ifdef EIGEN_DEFAULT_TO_ROW_MAJOR
|
||||
#undef EIGEN_DEFAULT_TO_ROW_MAJOR
|
||||
#endif
|
||||
@ -37,7 +37,7 @@
|
||||
|
||||
template<typename T> void test_sparselu_T()
|
||||
{
|
||||
SparseLU<SparseMatrix<T, ColMajor>, COLAMDOrdering<int> > sparselu_colamd;
|
||||
SparseLU<SparseMatrix<T, ColMajor> /*, COLAMDOrdering<int>*/ > sparselu_colamd; // COLAMDOrdering is the default
|
||||
SparseLU<SparseMatrix<T, ColMajor>, AMDOrdering<int> > sparselu_amd;
|
||||
SparseLU<SparseMatrix<T, ColMajor, long int>, NaturalOrdering<long int> > sparselu_natural;
|
||||
|
||||
|
Loading…
x
Reference in New Issue
Block a user