aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2015-06-26 10:49:40 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2015-06-26 10:49:40 +0200
commit555b9c684346129057b14d73db64a75274fdfe0f (patch)
treed35749e34b390790a9cd12c7e9309391baf69309
parent53b930887d118af5204840231f08b3307addce4e (diff)
Doc: explain perf and multithreading issues in sparse iterative solvers
-rw-r--r--Eigen/src/IterativeLinearSolvers/BiCGSTAB.h6
-rw-r--r--Eigen/src/IterativeLinearSolvers/ConjugateGradient.h12
-rw-r--r--doc/SparseLinearSystems.dox2
-rw-r--r--doc/TopicMultithreading.dox8
4 files changed, 21 insertions, 7 deletions
diff --git a/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h b/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h
index f6d54c5e7..a34ee7628 100644
--- a/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h
+++ b/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h
@@ -136,7 +136,11 @@ struct traits<BiCGSTAB<_MatrixType,_Preconditioner> >
* and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
* and NumTraits<Scalar>::epsilon() for the tolerance.
*
- * The tolerance is the relative residual error: |Ax-b|/|b|
+ * The tolerance corresponds to the relative residual error: |Ax-b|/|b|
+ *
+ * \b Performance: when using sparse matrices, best performance is achied for a row-major sparse matrix format.
+ * Moreover, in this case multi-threading can be exploited if the user code is compiled with OpenMP enabled.
+ * See \ref TopicMultiThreading for details.
*
* This class can be used as the direct solver classes. Here is a typical usage example:
* \include BiCGSTAB_simple.cpp
diff --git a/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h b/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
index 5f0159e52..8f33c446d 100644
--- a/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
+++ b/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
@@ -114,14 +114,20 @@ struct traits<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
*
* \tparam _MatrixType the type of the matrix A, can be a dense or a sparse matrix.
* \tparam _UpLo the triangular part that will be used for the computations. It can be Lower,
- * Upper, or Lower|Upper in which the full matrix entries will be considered. Default is Lower.
+ * \c Upper, or \c Lower|Upper in which the full matrix entries will be considered.
+ * Default is \c Lower, best performance is \c Lower|Upper.
* \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
*
* 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.
*
- * The tolerance is the relative residual error: |Ax-b|/|b|
+ * The tolerance corresponds to the relative residual error: |Ax-b|/|b|
+ *
+ * \b Performance: Even though the default value of \c _UpLo is \c Lower, significantly higher performance is
+ * achieved when using a complete matrix and \b Lower|Upper as the \a _UpLo template parameter. Moreover, in this
+ * case multi-threading can be exploited if the user code is compiled with OpenMP enabled.
+ * See \ref TopicMultiThreading for details.
*
* This class can be used as the direct solver classes. Here is a typical usage example:
\code
@@ -129,7 +135,7 @@ struct traits<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
VectorXd x(n), b(n);
SparseMatrix<double> A(n,n);
// fill A and b
- ConjugateGradient<SparseMatrix<double> > cg;
+ ConjugateGradient<SparseMatrix<double>, Lower|Upper> cg;
cg.compute(A);
x = cg.solve(b);
std::cout << "#iterations: " << cg.iterations() << std::endl;
diff --git a/doc/SparseLinearSystems.dox b/doc/SparseLinearSystems.dox
index 13741280a..48c18f46f 100644
--- a/doc/SparseLinearSystems.dox
+++ b/doc/SparseLinearSystems.dox
@@ -21,7 +21,7 @@ They are summarized in the following table:
<tr><td>ConjugateGradient</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>Classic iterative CG</td><td>SPD</td><td>Preconditionning</td>
<td>built-in, MPL2</td>
<td>Recommended for large symmetric problems (e.g., 3D Poisson eq.)</td></tr>
-<tr><td>LSCG</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>CG for rectangular least-square problem</td><td>Rectangular</td><td>Preconditionning</td>
+<tr><td>LeastSquaresConjugateGradient</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>CG for rectangular least-square problem</td><td>Rectangular</td><td>Preconditionning</td>
<td>built-in, MPL2</td>
<td>Solve for min |A'Ax-b|^2 without forming A'A</td></tr>
<tr><td>BiCGSTAB</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>Iterative stabilized bi-conjugate gradient</td><td>Square</td><td>Preconditionning</td>
diff --git a/doc/TopicMultithreading.dox b/doc/TopicMultithreading.dox
index e22e1c613..66028d7a8 100644
--- a/doc/TopicMultithreading.dox
+++ b/doc/TopicMultithreading.dox
@@ -22,8 +22,12 @@ n = Eigen::nbThreads( );
You can disable Eigen's multi threading at compile time by defining the EIGEN_DONT_PARALLELIZE preprocessor token.
Currently, the following algorithms can make use of multi-threading:
- * general matrix - matrix products
- * PartialPivLU
+ - general dense matrix - matrix products
+ - PartialPivLU
+ - row-major-sparse * dense vector/matrix products
+ - ConjugateGradient with \c Lower|Upper as the \c UpLo template parameter.
+ - BiCGSTAB with a row-major sparse matrix format.
+ - LeastSquaresConjugateGradient
\section TopicMultiThreading_UsingEigenWithMT Using Eigen in a multi-threaded application