From 26129229ec12961687b0414c40e10e2880beec79 Mon Sep 17 00:00:00 2001 From: Benoit Jacob Date: Fri, 15 Oct 2010 09:44:43 -0400 Subject: doc updates/improvements --- doc/C00_QuickStartGuide.dox | 5 +++-- doc/C01_TutorialMatrixClass.dox | 2 +- doc/C06_TutorialLinearAlgebra.dox | 15 ++++++++++++--- doc/QuickReference.dox | 4 ++-- doc/TopicLinearAlgebraDecompositions.dox | 21 +++++---------------- doc/examples/TutorialLinAlgSVDSolve.cpp | 15 +++++++++++++++ 6 files changed, 38 insertions(+), 24 deletions(-) create mode 100644 doc/examples/TutorialLinAlgSVDSolve.cpp diff --git a/doc/C00_QuickStartGuide.dox b/doc/C00_QuickStartGuide.dox index 3b7c405ca..e2381f9e3 100644 --- a/doc/C00_QuickStartGuide.dox +++ b/doc/C00_QuickStartGuide.dox @@ -83,8 +83,9 @@ The use of fixed-size matrices and vectors has two advantages. The compiler emit \section GettingStartedConclusion Where to go from here? -You could directly use our \ref QuickRefPage and class documentation, or if you do not yet feel ready for that, you could -read the longer \ref TutorialMatrixClass "Tutorial" in which the Eigen library is explained in more detail. +It's worth taking the time to read the \ref TutorialMatrixClass "long tutorial". + +However if you think you don't need it, you can directly use the classes documentation and our \ref QuickRefPage. \li \b Next: \ref TutorialMatrixClass diff --git a/doc/C01_TutorialMatrixClass.dox b/doc/C01_TutorialMatrixClass.dox index b5d7ec0c7..5e3098361 100644 --- a/doc/C01_TutorialMatrixClass.dox +++ b/doc/C01_TutorialMatrixClass.dox @@ -92,7 +92,7 @@ Matrix \section TutorialMatrixConstructors Constructors -A default constructor is always available, and always has zero runtime cost. You can do: +A default constructor is always available, never performs any dynamic memory allocation, and never initializes the matrix coefficients. You can do: \code Matrix3f a; MatrixXf b; diff --git a/doc/C06_TutorialLinearAlgebra.dox b/doc/C06_TutorialLinearAlgebra.dox index 3e436d393..c8f2bf23d 100644 --- a/doc/C06_TutorialLinearAlgebra.dox +++ b/doc/C06_TutorialLinearAlgebra.dox @@ -180,9 +180,18 @@ Here is an example: \section TutorialLinAlgLeastsquares Least squares solving -Eigen doesn't currently provide built-in linear least squares solving functions, but you can easily compute that yourself -from Eigen's decompositions. The most reliable way is to use a SVD (or better yet, JacobiSVD), and in the future -these classes will offer methods for least squares solving. Another, potentially faster way, is to use a LLT decomposition +The best way to do least squares solving is with a SVD decomposition. Eigen provides one as the JacobiSVD class, and its solve() +is doing least-squares solving. + +Here is an example: + + + + + +
\include TutorialLinAlgSVDSolve.cpp output: \verbinclude TutorialLinAlgSVDSolve.out
+ +Another way, potentially faster but less reliable, is to use a LDLT decomposition of the normal matrix. In any case, just read any reference text on least squares, and it will be very easy for you to implement any linear least squares computation on top of Eigen. diff --git a/doc/QuickReference.dox b/doc/QuickReference.dox index d426b85de..a7d42767c 100644 --- a/doc/QuickReference.dox +++ b/doc/QuickReference.dox @@ -27,11 +27,11 @@ The Eigen library is divided in a Core module and several additional modules. Ea \link LU_Module LU \endlink\code#include \endcodeInverse, determinant, LU decompositions with solver (FullPivLU, PartialPivLU) \link Cholesky_Module Cholesky \endlink\code#include \endcodeLLT and LDLT Cholesky factorization with solver \link Householder_Module Householder \endlink\code#include \endcodeHouseholder transformations; this module is used by several linear algebra modules -\link SVD_Module SVD \endlink\code#include \endcode%SVD decomposition with solver (SVD, JacobiSVD) +\link SVD_Module SVD \endlink\code#include \endcodeSVD decomposition with least-squares solver (JacobiSVD) \link QR_Module QR \endlink\code#include \endcodeQR decomposition with solver (HouseholderQR, ColPivHouseholderQR, FullPivHouseholderQR) \link Eigenvalues_Module Eigenvalues \endlink\code#include \endcodeEigenvalue, eigenvector decompositions (EigenSolver, SelfAdjointEigenSolver, ComplexEigenSolver) \link Sparse_Module Sparse \endlink\code#include \endcode%Sparse matrix storage and related basic linear algebra (SparseMatrix, DynamicSparseMatrix, SparseVector) -\code#include \endcodeIncludes Core, Geometry, LU, Cholesky, %SVD, QR, and Eigenvalues header files +\code#include \endcodeIncludes Core, Geometry, LU, Cholesky, SVD, QR, and Eigenvalues header files \code#include \endcodeIncludes %Dense and %Sparse header files (the whole Eigen library) diff --git a/doc/TopicLinearAlgebraDecompositions.dox b/doc/TopicLinearAlgebraDecompositions.dox index ad8d0abea..203a05dd8 100644 --- a/doc/TopicLinearAlgebraDecompositions.dox +++ b/doc/TopicLinearAlgebraDecompositions.dox @@ -112,27 +112,15 @@ namespace Eigen { \n Singular values and eigenvalues decompositions - SVD - - - Average - Good - Yes - Singular values/vectors, least squares - Yes - Average - - - - - - JacobiSVD + JacobiSVD (two-sided) - Slow (but fast for small matrices) - Proven + Excellent-Proven3 Yes Singular values/vectors, least squares - - + Yes (and does least squares) Excellent - - + R-SVD @@ -251,6 +239,7 @@ namespace Eigen {
  • \b 1: There exist two variants of the LDLT algorithm. Eigen's one produces a pure diagonal D matrix, and therefore it cannot handle indefinite matrices, unlike Lapack's one which produces a block diagonal D matrix.
  • \b 2: Eigenvalues, SVD and Schur decompositions rely on iterative algorithms. Their convergence speed depends on how well the eigenvalues are separated.
  • +
  • \b 3: Our JacobiSVD is two-sided, making for proven and optimal precision for square matrices. For non-square matrices, we have to use a QR preconditioner first. The default choice, ColPivHouseholderQR, is already very reliable, but if you want it to be proven, use FullPivHouseholderQR instead.
\section TopicLinAlgTerminology Terminology diff --git a/doc/examples/TutorialLinAlgSVDSolve.cpp b/doc/examples/TutorialLinAlgSVDSolve.cpp new file mode 100644 index 000000000..c75779d5f --- /dev/null +++ b/doc/examples/TutorialLinAlgSVDSolve.cpp @@ -0,0 +1,15 @@ +#include +#include + +using namespace std; +using namespace Eigen; + +int main() +{ + MatrixXf A = MatrixXf::Random(3, 2); + cout << "Here is the matrix A:\n" << A << endl; + VectorXf b = VectorXf::Random(3); + cout << "Here is the right hand side b:\n" << b << endl; + JacobiSVD svd(A, ComputeThinU | ComputeThinV); + cout << "The least-squares solution is:\n" << svd.solve(b) << endl; +} -- cgit v1.2.3