3 /** \eigenManualPage TutorialLinearAlgebra Linear algebra and decompositions
5 This page explains how to solve linear systems, compute various decompositions such as LU,
6 QR, %SVD, eigendecompositions... After reading this page, don't miss our
7 \link TopicLinearAlgebraDecompositions catalogue \endlink of dense matrix decompositions.
11 \section TutorialLinAlgBasicSolve Basic linear solving
13 \b The \b problem: You have a system of equations, that you have written as a single matrix equation
15 Where \a A and \a b are matrices (\a b could be a vector, as a special case). You want to find a solution \a x.
17 \b The \b solution: You can choose between various decompositions, depending on what your matrix \a A looks like,
18 and depending on whether you favor speed or accuracy. However, let's start with an example that works in all cases,
19 and is a good compromise:
20 <table class="example">
21 <tr><th>Example:</th><th>Output:</th></tr>
23 <td>\include TutorialLinAlgExSolveColPivHouseholderQR.cpp </td>
24 <td>\verbinclude TutorialLinAlgExSolveColPivHouseholderQR.out </td>
28 In this example, the colPivHouseholderQr() method returns an object of class ColPivHouseholderQR. Since here the
29 matrix is of type Matrix3f, this line could have been replaced by:
31 ColPivHouseholderQR<Matrix3f> dec(A);
32 Vector3f x = dec.solve(b);
35 Here, ColPivHouseholderQR is a QR decomposition with column pivoting. It's a good compromise for this tutorial, as it
36 works for all matrices while being quite fast. Here is a table of some other decompositions that you can choose from,
37 depending on your matrix and the trade-off you want to make:
39 <table class="manual">
41 <th>Decomposition</th>
43 <th>Requirements<br/>on the matrix</th>
44 <th>Speed<br/> (small-to-medium)</th>
45 <th>Speed<br/> (large)</th>
50 <td>partialPivLu()</td>
65 <td>HouseholderQR</td>
66 <td>householderQr()</td>
73 <td>ColPivHouseholderQR</td>
74 <td>colPivHouseholderQr()</td>
81 <td>FullPivHouseholderQR</td>
82 <td>fullPivHouseholderQr()</td>
91 <td>Positive definite</td>
99 <td>Positive or negative<br/> semidefinite</td>
114 All of these decompositions offer a solve() method that works as in the above example.
116 For example, if your matrix is positive definite, the above table says that a very good
117 choice is then the LLT or LDLT decomposition. Here's an example, also demonstrating that using a general
118 matrix (not a vector) as right hand side is possible.
120 <table class="example">
121 <tr><th>Example:</th><th>Output:</th></tr>
123 <td>\include TutorialLinAlgExSolveLDLT.cpp </td>
124 <td>\verbinclude TutorialLinAlgExSolveLDLT.out </td>
128 For a \ref TopicLinearAlgebraDecompositions "much more complete table" comparing all decompositions supported by Eigen (notice that Eigen
129 supports many other decompositions), see our special page on
130 \ref TopicLinearAlgebraDecompositions "this topic".
132 \section TutorialLinAlgSolutionExists Checking if a solution really exists
134 Only you know what error margin you want to allow for a solution to be considered valid.
135 So Eigen lets you do this computation for yourself, if you want to, as in this example:
137 <table class="example">
138 <tr><th>Example:</th><th>Output:</th></tr>
140 <td>\include TutorialLinAlgExComputeSolveError.cpp </td>
141 <td>\verbinclude TutorialLinAlgExComputeSolveError.out </td>
145 \section TutorialLinAlgEigensolving Computing eigenvalues and eigenvectors
147 You need an eigendecomposition here, see available such decompositions on \ref TopicLinearAlgebraDecompositions "this page".
148 Make sure to check if your matrix is self-adjoint, as is often the case in these problems. Here's an example using
149 SelfAdjointEigenSolver, it could easily be adapted to general matrices using EigenSolver or ComplexEigenSolver.
151 The computation of eigenvalues and eigenvectors does not necessarily converge, but such failure to converge is
152 very rare. The call to info() is to check for this possibility.
154 <table class="example">
155 <tr><th>Example:</th><th>Output:</th></tr>
157 <td>\include TutorialLinAlgSelfAdjointEigenSolver.cpp </td>
158 <td>\verbinclude TutorialLinAlgSelfAdjointEigenSolver.out </td>
162 \section TutorialLinAlgInverse Computing inverse and determinant
164 First of all, make sure that you really want this. While inverse and determinant are fundamental mathematical concepts,
165 in \em numerical linear algebra they are not as popular as in pure mathematics. Inverse computations are often
166 advantageously replaced by solve() operations, and the determinant is often \em not a good way of checking if a matrix
169 However, for \em very \em small matrices, the above is not true, and inverse and determinant can be very useful.
171 While certain decompositions, such as PartialPivLU and FullPivLU, offer inverse() and determinant() methods, you can also
172 call inverse() and determinant() directly on a matrix. If your matrix is of a very small fixed size (at most 4x4) this
173 allows Eigen to avoid performing a LU decomposition, and instead use formulas that are more efficient on such small matrices.
176 <table class="example">
177 <tr><th>Example:</th><th>Output:</th></tr>
179 <td>\include TutorialLinAlgInverseDeterminant.cpp </td>
180 <td>\verbinclude TutorialLinAlgInverseDeterminant.out </td>
184 \section TutorialLinAlgLeastsquares Least squares solving
186 The most accurate method to do least squares solving is with a SVD decomposition. Eigen provides one
187 as the JacobiSVD class, and its solve() is doing least-squares solving.
190 <table class="example">
191 <tr><th>Example:</th><th>Output:</th></tr>
193 <td>\include TutorialLinAlgSVDSolve.cpp </td>
194 <td>\verbinclude TutorialLinAlgSVDSolve.out </td>
198 Another methods, potentially faster but less reliable, are to use a Cholesky decomposition of the
199 normal matrix or a QR decomposition. Our page on \link LeastSquares least squares solving \endlink
203 \section TutorialLinAlgSeparateComputation Separating the computation from the construction
205 In the above examples, the decomposition was computed at the same time that the decomposition object was constructed.
206 There are however situations where you might want to separate these two things, for example if you don't know,
207 at the time of the construction, the matrix that you will want to decompose; or if you want to reuse an existing
208 decomposition object.
210 What makes this possible is that:
211 \li all decompositions have a default constructor,
212 \li all decompositions have a compute(matrix) method that does the computation, and that may be called again
213 on an already-computed decomposition, reinitializing it.
217 <table class="example">
218 <tr><th>Example:</th><th>Output:</th></tr>
220 <td>\include TutorialLinAlgComputeTwice.cpp </td>
221 <td>\verbinclude TutorialLinAlgComputeTwice.out </td>
225 Finally, you can tell the decomposition constructor to preallocate storage for decomposing matrices of a given size,
226 so that when you subsequently decompose such matrices, no dynamic memory allocation is performed (of course, if you
227 are using fixed-size matrices, no dynamic memory allocation happens at all). This is done by just
228 passing the size to the decomposition constructor, as in this example:
230 HouseholderQR<MatrixXf> qr(50,50);
231 MatrixXf A = MatrixXf::Random(50,50);
232 qr.compute(A); // no dynamic memory allocation
235 \section TutorialLinAlgRankRevealing Rank-revealing decompositions
237 Certain decompositions are rank-revealing, i.e. are able to compute the rank of a matrix. These are typically
238 also the decompositions that behave best in the face of a non-full-rank matrix (which in the square case means a
239 singular matrix). On \ref TopicLinearAlgebraDecompositions "this table" you can see for all our decompositions
240 whether they are rank-revealing or not.
242 Rank-revealing decompositions offer at least a rank() method. They can also offer convenience methods such as isInvertible(),
243 and some are also providing methods to compute the kernel (null-space) and image (column-space) of the matrix, as is the
246 <table class="example">
247 <tr><th>Example:</th><th>Output:</th></tr>
249 <td>\include TutorialLinAlgRankRevealing.cpp </td>
250 <td>\verbinclude TutorialLinAlgRankRevealing.out </td>
254 Of course, any rank computation depends on the choice of an arbitrary threshold, since practically no
255 floating-point matrix is \em exactly rank-deficient. Eigen picks a sensible default threshold, which depends
256 on the decomposition but is typically the diagonal size times machine epsilon. While this is the best default we
257 could pick, only you know what is the right threshold for your application. You can set this by calling setThreshold()
258 on your decomposition object before calling rank() or any other method that needs to use such a threshold.
259 The decomposition itself, i.e. the compute() method, is independent of the threshold. You don't need to recompute the
260 decomposition after you've changed the threshold.
262 <table class="example">
263 <tr><th>Example:</th><th>Output:</th></tr>
265 <td>\include TutorialLinAlgSetThreshold.cpp </td>
266 <td>\verbinclude TutorialLinAlgSetThreshold.out </td>