aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/src/CholmodSupport/CholmodSupport.h8
-rw-r--r--Eigen/src/IterativeLinearSolvers/BasicPreconditioners.h6
-rw-r--r--Eigen/src/IterativeLinearSolvers/BiCGSTAB.h2
-rw-r--r--Eigen/src/IterativeLinearSolvers/ConjugateGradient.h2
-rw-r--r--Eigen/src/IterativeLinearSolvers/IncompleteLUT.h2
-rw-r--r--Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h2
-rw-r--r--Eigen/src/PaStiXSupport/PaStiXSupport.h12
-rwxr-xr-xEigen/src/PardisoSupport/PardisoSupport.h6
-rw-r--r--Eigen/src/SPQRSupport/SuiteSparseQRSupport.h44
-rw-r--r--Eigen/src/SparseCholesky/SimplicialCholesky.h4
-rw-r--r--Eigen/src/SparseLU/SparseLU.h3
-rw-r--r--Eigen/src/SparseQR/SparseQR.h2
-rw-r--r--Eigen/src/SuperLUSupport/SuperLUSupport.h4
-rw-r--r--doc/Doxyfile.in3
-rw-r--r--doc/SparseLinearSystems.dox4
15 files changed, 77 insertions, 27 deletions
diff --git a/Eigen/src/CholmodSupport/CholmodSupport.h b/Eigen/src/CholmodSupport/CholmodSupport.h
index d2b0fb282..3ff0c6fc9 100644
--- a/Eigen/src/CholmodSupport/CholmodSupport.h
+++ b/Eigen/src/CholmodSupport/CholmodSupport.h
@@ -350,6 +350,8 @@ class CholmodBase : public SparseSolverBase<Derived>
* \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
* or Upper. Default is Lower.
*
+ * \implsparsesolverconcept
+ *
* This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
*
* \sa \ref TutorialSparseDirectSolvers, class CholmodSupernodalLLT, class SimplicialLLT
@@ -397,6 +399,8 @@ class CholmodSimplicialLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimpl
* \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
* or Upper. Default is Lower.
*
+ * \implsparsesolverconcept
+ *
* This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
*
* \sa \ref TutorialSparseDirectSolvers, class CholmodSupernodalLLT, class SimplicialLDLT
@@ -442,6 +446,8 @@ class CholmodSimplicialLDLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimp
* \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
* or Upper. Default is Lower.
*
+ * \implsparsesolverconcept
+ *
* This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
*
* \sa \ref TutorialSparseDirectSolvers
@@ -489,6 +495,8 @@ class CholmodSupernodalLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSuper
* \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
* or Upper. Default is Lower.
*
+ * \implsparsesolverconcept
+ *
* This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
*
* \sa \ref TutorialSparseDirectSolvers
diff --git a/Eigen/src/IterativeLinearSolvers/BasicPreconditioners.h b/Eigen/src/IterativeLinearSolvers/BasicPreconditioners.h
index ff7f08c1c..b850630a3 100644
--- a/Eigen/src/IterativeLinearSolvers/BasicPreconditioners.h
+++ b/Eigen/src/IterativeLinearSolvers/BasicPreconditioners.h
@@ -23,6 +23,8 @@ namespace Eigen {
*
* \tparam _Scalar the type of the scalar.
*
+ * \implsparsesolverconcept
+ *
* This preconditioner is suitable for both selfadjoint and general problems.
* The diagonal entries are pre-inverted and stored into a dense vector.
*
@@ -114,6 +116,8 @@ class DiagonalPreconditioner
*
* \tparam _Scalar the type of the scalar.
*
+ * \implsparsesolverconcept
+ *
* The diagonal entries are pre-inverted and stored into a dense vector.
*
* \sa class LeastSquaresConjugateGradient, class DiagonalPreconditioner
@@ -172,6 +176,8 @@ class LeastSquareDiagonalPreconditioner : public DiagonalPreconditioner<_Scalar>
/** \ingroup IterativeLinearSolvers_Module
* \brief A naive preconditioner which approximates any matrix as the identity matrix
*
+ * \implsparsesolverconcept
+ *
* \sa class DiagonalPreconditioner
*/
class IdentityPreconditioner
diff --git a/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h b/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h
index a34ee7628..76e86a94a 100644
--- a/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h
+++ b/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h
@@ -132,6 +132,8 @@ struct traits<BiCGSTAB<_MatrixType,_Preconditioner> >
* \tparam _MatrixType the type of the sparse matrix A, can be a dense or a sparse matrix.
* \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
*
+ * \implsparsesolverconcept
+ *
* The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
* and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
* and NumTraits<Scalar>::epsilon() for the tolerance.
diff --git a/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h b/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
index 8f33c446d..59092dc18 100644
--- a/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
+++ b/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
@@ -118,6 +118,8 @@ struct traits<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
* Default is \c Lower, best performance is \c Lower|Upper.
* \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
*
+ * \implsparsesolverconcept
+ *
* The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
* and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
* and NumTraits<Scalar>::epsilon() for the tolerance.
diff --git a/Eigen/src/IterativeLinearSolvers/IncompleteLUT.h b/Eigen/src/IterativeLinearSolvers/IncompleteLUT.h
index b644163f1..10b9fcc18 100644
--- a/Eigen/src/IterativeLinearSolvers/IncompleteLUT.h
+++ b/Eigen/src/IterativeLinearSolvers/IncompleteLUT.h
@@ -67,6 +67,8 @@ Index QuickSplit(VectorV &row, VectorI &ind, Index ncut)
* \class IncompleteLUT
* \brief Incomplete LU factorization with dual-threshold strategy
*
+ * \implsparsesolverconcept
+ *
* During the numerical factorization, two dropping rules are used :
* 1) any element whose magnitude is less than some tolerance is dropped.
* This tolerance is obtained by multiplying the input tolerance @p droptol
diff --git a/Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h b/Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h
index 1d819927e..b578b2a7f 100644
--- a/Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h
+++ b/Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h
@@ -119,6 +119,8 @@ struct traits<LeastSquaresConjugateGradient<_MatrixType,_Preconditioner> >
* \tparam _MatrixType the type of the matrix A, can be a dense or a sparse matrix.
* \tparam _Preconditioner the type of the preconditioner. Default is LeastSquareDiagonalPreconditioner
*
+ * \implsparsesolverconcept
+ *
* The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
* and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
* and NumTraits<Scalar>::epsilon() for the tolerance.
diff --git a/Eigen/src/PaStiXSupport/PaStiXSupport.h b/Eigen/src/PaStiXSupport/PaStiXSupport.h
index 4e73edf5b..cec4149e7 100644
--- a/Eigen/src/PaStiXSupport/PaStiXSupport.h
+++ b/Eigen/src/PaStiXSupport/PaStiXSupport.h
@@ -398,7 +398,9 @@ bool PastixBase<Base>::_solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &x
* NOTE : Note that if the analysis and factorization phase are called separately,
* the input matrix will be symmetrized at each call, hence it is advised to
* symmetrize the matrix in a end-user program and set \p IsStrSym to true
- *
+ *
+ * \implsparsesolverconcept
+ *
* \sa \ref TutorialSparseDirectSolvers
*
*/
@@ -509,7 +511,9 @@ class PastixLU : public PastixBase< PastixLU<_MatrixType> >
*
* \tparam MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
* \tparam UpLo The part of the matrix to use : Lower or Upper. The default is Lower as required by PaStiX
- *
+ *
+ * \implsparsesolverconcept
+ *
* \sa \ref TutorialSparseDirectSolvers
*/
template<typename _MatrixType, int _UpLo>
@@ -590,7 +594,9 @@ class PastixLLT : public PastixBase< PastixLLT<_MatrixType, _UpLo> >
*
* \tparam MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
* \tparam UpLo The part of the matrix to use : Lower or Upper. The default is Lower as required by PaStiX
- *
+ *
+ * \implsparsesolverconcept
+ *
* \sa \ref TutorialSparseDirectSolvers
*/
template<typename _MatrixType, int _UpLo>
diff --git a/Eigen/src/PardisoSupport/PardisoSupport.h b/Eigen/src/PardisoSupport/PardisoSupport.h
index 234e3213b..9c18eb9b9 100755
--- a/Eigen/src/PardisoSupport/PardisoSupport.h
+++ b/Eigen/src/PardisoSupport/PardisoSupport.h
@@ -371,6 +371,8 @@ void PardisoImpl<Derived>::_solve_impl(const MatrixBase<BDerived> &b, MatrixBase
*
* \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
*
+ * \implsparsesolverconcept
+ *
* \sa \ref TutorialSparseDirectSolvers
*/
template<typename MatrixType>
@@ -421,6 +423,8 @@ class PardisoLU : public PardisoImpl< PardisoLU<MatrixType> >
* \tparam UpLo can be any bitwise combination of Upper, Lower. The default is Upper, meaning only the upper triangular part has to be used.
* Upper|Lower can be used to tell both triangular parts can be used as input.
*
+ * \implsparsesolverconcept
+ *
* \sa \ref TutorialSparseDirectSolvers
*/
template<typename MatrixType, int _UpLo>
@@ -479,6 +483,8 @@ class PardisoLLT : public PardisoImpl< PardisoLLT<MatrixType,_UpLo> >
* Symmetric can be used for symmetric, non-selfadjoint complex matrices, the default being to assume a selfadjoint matrix.
* Upper|Lower can be used to tell both triangular parts can be used as input.
*
+ * \implsparsesolverconcept
+ *
* \sa \ref TutorialSparseDirectSolvers
*/
template<typename MatrixType, int Options>
diff --git a/Eigen/src/SPQRSupport/SuiteSparseQRSupport.h b/Eigen/src/SPQRSupport/SuiteSparseQRSupport.h
index 4ad22f8b4..ac2de9b04 100644
--- a/Eigen/src/SPQRSupport/SuiteSparseQRSupport.h
+++ b/Eigen/src/SPQRSupport/SuiteSparseQRSupport.h
@@ -33,27 +33,29 @@ namespace Eigen {
} // End namespace internal
/**
- * \ingroup SPQRSupport_Module
- * \class SPQR
- * \brief Sparse QR factorization based on SuiteSparseQR library
- *
- * This class is used to perform a multithreaded and multifrontal rank-revealing QR decomposition
- * of sparse matrices. The result is then used to solve linear leasts_square systems.
- * Clearly, a QR factorization is returned such that A*P = Q*R where :
- *
- * P is the column permutation. Use colsPermutation() to get it.
- *
- * Q is the orthogonal matrix represented as Householder reflectors.
- * Use matrixQ() to get an expression and matrixQ().transpose() to get the transpose.
- * You can then apply it to a vector.
- *
- * R is the sparse triangular factor. Use matrixQR() to get it as SparseMatrix.
- * NOTE : The Index type of R is always UF_long. You can get it with SPQR::Index
- *
- * \tparam _MatrixType The type of the sparse matrix A, must be a column-major SparseMatrix<>
- * NOTE
- *
- */
+ * \ingroup SPQRSupport_Module
+ * \class SPQR
+ * \brief Sparse QR factorization based on SuiteSparseQR library
+ *
+ * This class is used to perform a multithreaded and multifrontal rank-revealing QR decomposition
+ * of sparse matrices. The result is then used to solve linear leasts_square systems.
+ * Clearly, a QR factorization is returned such that A*P = Q*R where :
+ *
+ * P is the column permutation. Use colsPermutation() to get it.
+ *
+ * Q is the orthogonal matrix represented as Householder reflectors.
+ * Use matrixQ() to get an expression and matrixQ().transpose() to get the transpose.
+ * You can then apply it to a vector.
+ *
+ * R is the sparse triangular factor. Use matrixQR() to get it as SparseMatrix.
+ * NOTE : The Index type of R is always UF_long. You can get it with SPQR::Index
+ *
+ * \tparam _MatrixType The type of the sparse matrix A, must be a column-major SparseMatrix<>
+ *
+ * \implsparsesolverconcept
+ *
+ *
+ */
template<typename _MatrixType>
class SPQR : public SparseSolverBase<SPQR<_MatrixType> >
{
diff --git a/Eigen/src/SparseCholesky/SimplicialCholesky.h b/Eigen/src/SparseCholesky/SimplicialCholesky.h
index f56298e8c..ef612cf45 100644
--- a/Eigen/src/SparseCholesky/SimplicialCholesky.h
+++ b/Eigen/src/SparseCholesky/SimplicialCholesky.h
@@ -319,6 +319,8 @@ template<typename _MatrixType, int _UpLo, typename _Ordering> struct traits<Simp
* or Upper. Default is Lower.
* \tparam _Ordering The ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<>
*
+ * \implsparsesolverconcept
+ *
* \sa class SimplicialLDLT, class AMDOrdering, class NaturalOrdering
*/
template<typename _MatrixType, int _UpLo, typename _Ordering>
@@ -408,6 +410,8 @@ public:
* or Upper. Default is Lower.
* \tparam _Ordering The ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<>
*
+ * \implsparsesolverconcept
+ *
* \sa class SimplicialLLT, class AMDOrdering, class NaturalOrdering
*/
template<typename _MatrixType, int _UpLo, typename _Ordering>
diff --git a/Eigen/src/SparseLU/SparseLU.h b/Eigen/src/SparseLU/SparseLU.h
index 8cdd29c7b..73368cba4 100644
--- a/Eigen/src/SparseLU/SparseLU.h
+++ b/Eigen/src/SparseLU/SparseLU.h
@@ -64,7 +64,8 @@ template <typename MatrixLType, typename MatrixUType> struct SparseLUMatrixURetu
*
* \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. Default is COLMAD
- *
+ *
+ * \implsparsesolverconcept
*
* \sa \ref TutorialSparseDirectSolvers
* \sa \ref OrderingMethods_Module
diff --git a/Eigen/src/SparseQR/SparseQR.h b/Eigen/src/SparseQR/SparseQR.h
index 548b3f9b0..bbd337c40 100644
--- a/Eigen/src/SparseQR/SparseQR.h
+++ b/Eigen/src/SparseQR/SparseQR.h
@@ -62,6 +62,8 @@ namespace internal {
* \tparam _OrderingType The fill-reducing ordering method. See the \link OrderingMethods_Module
* OrderingMethods \endlink module for the list of built-in and external ordering methods.
*
+ * \implsparsesolverconcept
+ *
* \warning The input sparse matrix A must be in compressed mode (see SparseMatrix::makeCompressed()).
*
*/
diff --git a/Eigen/src/SuperLUSupport/SuperLUSupport.h b/Eigen/src/SuperLUSupport/SuperLUSupport.h
index d067d8fdf..7c644eef6 100644
--- a/Eigen/src/SuperLUSupport/SuperLUSupport.h
+++ b/Eigen/src/SuperLUSupport/SuperLUSupport.h
@@ -449,6 +449,8 @@ class SuperLUBase : public SparseSolverBase<Derived>
*
* \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
*
+ * \implsparsesolverconcept
+ *
* \sa \ref TutorialSparseDirectSolvers
*/
template<typename _MatrixType>
@@ -800,6 +802,8 @@ typename SuperLU<MatrixType>::Scalar SuperLU<MatrixType>::determinant() const
*
* \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
*
+ * \implsparsesolverconcept
+ *
* \sa \ref TutorialSparseDirectSolvers, class ConjugateGradient, class BiCGSTAB
*/
diff --git a/doc/Doxyfile.in b/doc/Doxyfile.in
index 800bb30ee..e15ba84bd 100644
--- a/doc/Doxyfile.in
+++ b/doc/Doxyfile.in
@@ -223,7 +223,8 @@ ALIASES = "only_for_vectors=This is only for vectors (either row-
"note_about_using_kernel_to_study_multiple_solutions=If you need a complete analysis of the space of solutions, take the one solution obtained by this method and add to it elements of the kernel, as determined by kernel()." \
"note_about_checking_solutions=This method just tries to find as good a solution as possible. If you want to check whether a solution exists or if it is accurate, just call this function to get a result and then compute the error of this result, or use MatrixBase::isApprox() directly, for instance like this: \code bool a_solution_exists = (A*result).isApprox(b, precision); \endcode This method avoids dividing by zero, so that the non-existence of a solution doesn't by itself mean that you'll get \c inf or \c nan values." \
"note_try_to_help_rvo=This function returns the result by value. In order to make that efficient, it is implemented as just a return statement using a special constructor, hopefully allowing the compiler to perform a RVO (return value optimization)." \
- "nonstableyet=\warning This is not considered to be part of the stable public API yet. Changes may happen in future releases. See \ref Experimental \"Experimental parts of Eigen\""
+ "nonstableyet=\warning This is not considered to be part of the stable public API yet. Changes may happen in future releases. See \ref Experimental \"Experimental parts of Eigen\"" \
+ "implsparsesolverconcept=This class follows the \link TutorialSparseSolverConcept sparse solver concept \endlink."
ALIASES += "eigenAutoToc= "
diff --git a/doc/SparseLinearSystems.dox b/doc/SparseLinearSystems.dox
index 48c18f46f..b7f5c600b 100644
--- a/doc/SparseLinearSystems.dox
+++ b/doc/SparseLinearSystems.dox
@@ -4,7 +4,7 @@ In Eigen, there are several methods available to solve linear systems when the c
\eigenAutoToc
-\section TutorialSparseDirectSolvers Sparse solvers
+\section TutorialSparseSolverList List of sparse solvers
%Eigen currently provides a limited set of built-in solvers, as well as wrappers to external solver libraries.
They are summarized in the following table:
@@ -53,6 +53,8 @@ They are summarized in the following table:
Here \c SPD means symmetric positive definite.
+\section TutorialSparseSolverConcept Sparse solver concept
+
All these solvers follow the same general concept.
Here is a typical and general example:
\code