• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1namespace Eigen {
2
3/** \page TutorialLinearAlgebra Tutorial page 6 - Linear algebra and decompositions
4    \ingroup Tutorial
5
6\li \b Previous: \ref TutorialAdvancedInitialization
7\li \b Next: \ref TutorialReductionsVisitorsBroadcasting
8
9This tutorial explains how to solve linear systems, compute various decompositions such as LU,
10QR, %SVD, eigendecompositions... for more advanced topics, don't miss our special page on
11\ref TopicLinearAlgebraDecompositions "this topic".
12
13\b Table \b of \b contents
14  - \ref TutorialLinAlgBasicSolve
15  - \ref TutorialLinAlgSolutionExists
16  - \ref TutorialLinAlgEigensolving
17  - \ref TutorialLinAlgInverse
18  - \ref TutorialLinAlgLeastsquares
19  - \ref TutorialLinAlgSeparateComputation
20  - \ref TutorialLinAlgRankRevealing
21
22
23\section TutorialLinAlgBasicSolve Basic linear solving
24
25\b The \b problem: You have a system of equations, that you have written as a single matrix equation
26    \f[ Ax \: = \: b \f]
27Where \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.
28
29\b The \b solution: You can choose between various decompositions, depending on what your matrix \a A looks like,
30and depending on whether you favor speed or accuracy. However, let's start with an example that works in all cases,
31and is a good compromise:
32<table class="example">
33<tr><th>Example:</th><th>Output:</th></tr>
34<tr>
35  <td>\include TutorialLinAlgExSolveColPivHouseholderQR.cpp </td>
36  <td>\verbinclude TutorialLinAlgExSolveColPivHouseholderQR.out </td>
37</tr>
38</table>
39
40In this example, the colPivHouseholderQr() method returns an object of class ColPivHouseholderQR. Since here the
41matrix is of type Matrix3f, this line could have been replaced by:
42\code
43ColPivHouseholderQR<Matrix3f> dec(A);
44Vector3f x = dec.solve(b);
45\endcode
46
47Here, ColPivHouseholderQR is a QR decomposition with column pivoting. It's a good compromise for this tutorial, as it
48works for all matrices while being quite fast. Here is a table of some other decompositions that you can choose from,
49depending on your matrix and the trade-off you want to make:
50
51<table class="manual">
52    <tr>
53        <th>Decomposition</th>
54        <th>Method</th>
55        <th>Requirements on the matrix</th>
56        <th>Speed</th>
57        <th>Accuracy</th>
58    </tr>
59    <tr>
60        <td>PartialPivLU</td>
61        <td>partialPivLu()</td>
62        <td>Invertible</td>
63        <td>++</td>
64        <td>+</td>
65    </tr>
66    <tr class="alt">
67        <td>FullPivLU</td>
68        <td>fullPivLu()</td>
69        <td>None</td>
70        <td>-</td>
71        <td>+++</td>
72    </tr>
73    <tr>
74        <td>HouseholderQR</td>
75        <td>householderQr()</td>
76        <td>None</td>
77        <td>++</td>
78        <td>+</td>
79    </tr>
80    <tr class="alt">
81        <td>ColPivHouseholderQR</td>
82        <td>colPivHouseholderQr()</td>
83        <td>None</td>
84        <td>+</td>
85        <td>++</td>
86    </tr>
87    <tr>
88        <td>FullPivHouseholderQR</td>
89        <td>fullPivHouseholderQr()</td>
90        <td>None</td>
91        <td>-</td>
92        <td>+++</td>
93    </tr>
94    <tr class="alt">
95        <td>LLT</td>
96        <td>llt()</td>
97        <td>Positive definite</td>
98        <td>+++</td>
99        <td>+</td>
100    </tr>
101    <tr>
102        <td>LDLT</td>
103        <td>ldlt()</td>
104        <td>Positive or negative semidefinite</td>
105        <td>+++</td>
106        <td>++</td>
107    </tr>
108</table>
109
110All of these decompositions offer a solve() method that works as in the above example.
111
112For example, if your matrix is positive definite, the above table says that a very good
113choice is then the LDLT decomposition. Here's an example, also demonstrating that using a general
114matrix (not a vector) as right hand side is possible.
115
116<table class="example">
117<tr><th>Example:</th><th>Output:</th></tr>
118<tr>
119  <td>\include TutorialLinAlgExSolveLDLT.cpp </td>
120  <td>\verbinclude TutorialLinAlgExSolveLDLT.out </td>
121</tr>
122</table>
123
124For a \ref TopicLinearAlgebraDecompositions "much more complete table" comparing all decompositions supported by Eigen (notice that Eigen
125supports many other decompositions), see our special page on
126\ref TopicLinearAlgebraDecompositions "this topic".
127
128\section TutorialLinAlgSolutionExists Checking if a solution really exists
129
130Only you know what error margin you want to allow for a solution to be considered valid.
131So Eigen lets you do this computation for yourself, if you want to, as in this example:
132
133<table class="example">
134<tr><th>Example:</th><th>Output:</th></tr>
135<tr>
136  <td>\include TutorialLinAlgExComputeSolveError.cpp </td>
137  <td>\verbinclude TutorialLinAlgExComputeSolveError.out </td>
138</tr>
139</table>
140
141\section TutorialLinAlgEigensolving Computing eigenvalues and eigenvectors
142
143You need an eigendecomposition here, see available such decompositions on \ref TopicLinearAlgebraDecompositions "this page".
144Make sure to check if your matrix is self-adjoint, as is often the case in these problems. Here's an example using
145SelfAdjointEigenSolver, it could easily be adapted to general matrices using EigenSolver or ComplexEigenSolver.
146
147The computation of eigenvalues and eigenvectors does not necessarily converge, but such failure to converge is
148very rare. The call to info() is to check for this possibility.
149
150<table class="example">
151<tr><th>Example:</th><th>Output:</th></tr>
152<tr>
153  <td>\include TutorialLinAlgSelfAdjointEigenSolver.cpp </td>
154  <td>\verbinclude TutorialLinAlgSelfAdjointEigenSolver.out </td>
155</tr>
156</table>
157
158\section TutorialLinAlgInverse Computing inverse and determinant
159
160First of all, make sure that you really want this. While inverse and determinant are fundamental mathematical concepts,
161in \em numerical linear algebra they are not as popular as in pure mathematics. Inverse computations are often
162advantageously replaced by solve() operations, and the determinant is often \em not a good way of checking if a matrix
163is invertible.
164
165However, for \em very \em small matrices, the above is not true, and inverse and determinant can be very useful.
166
167While certain decompositions, such as PartialPivLU and FullPivLU, offer inverse() and determinant() methods, you can also
168call inverse() and determinant() directly on a matrix. If your matrix is of a very small fixed size (at most 4x4) this
169allows Eigen to avoid performing a LU decomposition, and instead use formulas that are more efficient on such small matrices.
170
171Here is an example:
172<table class="example">
173<tr><th>Example:</th><th>Output:</th></tr>
174<tr>
175  <td>\include TutorialLinAlgInverseDeterminant.cpp </td>
176  <td>\verbinclude TutorialLinAlgInverseDeterminant.out </td>
177</tr>
178</table>
179
180\section TutorialLinAlgLeastsquares Least squares solving
181
182The best way to do least squares solving is with a SVD decomposition. Eigen provides one as the JacobiSVD class, and its solve()
183is doing least-squares solving.
184
185Here is an example:
186<table class="example">
187<tr><th>Example:</th><th>Output:</th></tr>
188<tr>
189  <td>\include TutorialLinAlgSVDSolve.cpp </td>
190  <td>\verbinclude TutorialLinAlgSVDSolve.out </td>
191</tr>
192</table>
193
194Another way, potentially faster but less reliable, is to use a LDLT decomposition
195of the normal matrix. In any case, just read any reference text on least squares, and it will be very easy for you
196to implement any linear least squares computation on top of Eigen.
197
198\section TutorialLinAlgSeparateComputation Separating the computation from the construction
199
200In the above examples, the decomposition was computed at the same time that the decomposition object was constructed.
201There are however situations where you might want to separate these two things, for example if you don't know,
202at the time of the construction, the matrix that you will want to decompose; or if you want to reuse an existing
203decomposition object.
204
205What makes this possible is that:
206\li all decompositions have a default constructor,
207\li all decompositions have a compute(matrix) method that does the computation, and that may be called again
208    on an already-computed decomposition, reinitializing it.
209
210For example:
211
212<table class="example">
213<tr><th>Example:</th><th>Output:</th></tr>
214<tr>
215  <td>\include TutorialLinAlgComputeTwice.cpp </td>
216  <td>\verbinclude TutorialLinAlgComputeTwice.out </td>
217</tr>
218</table>
219
220Finally, you can tell the decomposition constructor to preallocate storage for decomposing matrices of a given size,
221so that when you subsequently decompose such matrices, no dynamic memory allocation is performed (of course, if you
222are using fixed-size matrices, no dynamic memory allocation happens at all). This is done by just
223passing the size to the decomposition constructor, as in this example:
224\code
225HouseholderQR<MatrixXf> qr(50,50);
226MatrixXf A = MatrixXf::Random(50,50);
227qr.compute(A); // no dynamic memory allocation
228\endcode
229
230\section TutorialLinAlgRankRevealing Rank-revealing decompositions
231
232Certain decompositions are rank-revealing, i.e. are able to compute the rank of a matrix. These are typically
233also the decompositions that behave best in the face of a non-full-rank matrix (which in the square case means a
234singular matrix). On \ref TopicLinearAlgebraDecompositions "this table" you can see for all our decompositions
235whether they are rank-revealing or not.
236
237Rank-revealing decompositions offer at least a rank() method. They can also offer convenience methods such as isInvertible(),
238and some are also providing methods to compute the kernel (null-space) and image (column-space) of the matrix, as is the
239case with FullPivLU:
240
241<table class="example">
242<tr><th>Example:</th><th>Output:</th></tr>
243<tr>
244  <td>\include TutorialLinAlgRankRevealing.cpp </td>
245  <td>\verbinclude TutorialLinAlgRankRevealing.out </td>
246</tr>
247</table>
248
249Of course, any rank computation depends on the choice of an arbitrary threshold, since practically no
250floating-point matrix is \em exactly rank-deficient. Eigen picks a sensible default threshold, which depends
251on the decomposition but is typically the diagonal size times machine epsilon. While this is the best default we
252could pick, only you know what is the right threshold for your application. You can set this by calling setThreshold()
253on your decomposition object before calling rank() or any other method that needs to use such a threshold.
254The decomposition itself, i.e. the compute() method, is independent of the threshold. You don't need to recompute the
255decomposition after you've changed the threshold.
256
257<table class="example">
258<tr><th>Example:</th><th>Output:</th></tr>
259<tr>
260  <td>\include TutorialLinAlgSetThreshold.cpp </td>
261  <td>\verbinclude TutorialLinAlgSetThreshold.out </td>
262</tr>
263</table>
264
265\li \b Next: \ref TutorialReductionsVisitorsBroadcasting
266
267*/
268
269}
270