blob: 77f13f4a0e534a4b532546fe376b349bbb3f4d70 [file] [log] [blame]
Benoit Jacob76152e92010-06-29 10:02:33 -04001namespace Eigen {
2
3/** \page TutorialLinearAlgebra Tutorial page 6 - Linear algebra and decompositions
4 \ingroup Tutorial
5
Benoit Jacob4d4a23c2010-06-30 10:11:55 -04006\li \b Previous: \ref TutorialAdvancedInitialization
Jitse Niesen26cfe5a2010-07-09 11:59:29 +01007\li \b Next: \ref TutorialReductionsVisitorsBroadcasting
Benoit Jacob76152e92010-06-29 10:02:33 -04008
9This tutorial explains how to solve linear systems, compute various decompositions such as LU,
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040010QR, %SVD, eigendecompositions... for more advanced topics, don't miss our special page on
Benoit Jacob76152e92010-06-29 10:02:33 -040011\ref TopicLinearAlgebraDecompositions "this topic".
12
Jitse Niesend0f6b1c2010-07-22 21:52:04 +010013\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
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040023\section TutorialLinAlgBasicSolve Basic linear solving
Benoit Jacob76152e92010-06-29 10:02:33 -040024
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:
Gael Guennebaudf66fe262010-10-19 11:40:49 +020032<table class="example">
33<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040034<tr>
35 <td>\include TutorialLinAlgExSolveColPivHouseholderQR.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020036 <td>\verbinclude TutorialLinAlgExSolveColPivHouseholderQR.out </td>
Benoit Jacob76152e92010-06-29 10:02:33 -040037</tr>
38</table>
39
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040040In 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:
Benoit Jacob76152e92010-06-29 10:02:33 -040042\code
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040043ColPivHouseholderQR<Matrix3f> dec(A);
Benoit Jacob76152e92010-06-29 10:02:33 -040044Vector3f 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
Gael Guennebaudf66fe262010-10-19 11:40:49 +020051<table class="manual">
Benoit Jacob76152e92010-06-29 10:02:33 -040052 <tr>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020053 <th>Decomposition</th>
54 <th>Method</th>
55 <th>Requirements on the matrix</th>
56 <th>Speed</th>
57 <th>Accuracy</th>
Benoit Jacob76152e92010-06-29 10:02:33 -040058 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040059 <tr>
60 <td>PartialPivLU</td>
61 <td>partialPivLu()</td>
62 <td>Invertible</td>
63 <td>++</td>
64 <td>+</td>
65 </tr>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020066 <tr class="alt">
Benoit Jacob76152e92010-06-29 10:02:33 -040067 <td>FullPivLU</td>
68 <td>fullPivLu()</td>
69 <td>None</td>
70 <td>-</td>
71 <td>+++</td>
72 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040073 <tr>
74 <td>HouseholderQR</td>
75 <td>householderQr()</td>
76 <td>None</td>
77 <td>++</td>
78 <td>+</td>
79 </tr>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020080 <tr class="alt">
Benoit Jacob76152e92010-06-29 10:02:33 -040081 <td>ColPivHouseholderQR</td>
82 <td>colPivHouseholderQr()</td>
83 <td>None</td>
84 <td>+</td>
85 <td>++</td>
86 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040087 <tr>
88 <td>FullPivHouseholderQR</td>
89 <td>fullPivHouseholderQr()</td>
90 <td>None</td>
91 <td>-</td>
92 <td>+++</td>
93 </tr>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020094 <tr class="alt">
Benoit Jacob76152e92010-06-29 10:02:33 -040095 <td>LLT</td>
96 <td>llt()</td>
97 <td>Positive definite</td>
98 <td>+++</td>
99 <td>+</td>
100 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -0400101 <tr>
102 <td>LDLT</td>
103 <td>ldlt()</td>
104 <td>Positive or negative semidefinite</td>
105 <td>+++</td>
106 <td>++</td>
107 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -0400108</table>
109
110All of these decompositions offer a solve() method that works as in the above example.
111
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400112For 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
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200116<table class="example">
117<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400118<tr>
119 <td>\include TutorialLinAlgExSolveLDLT.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200120 <td>\verbinclude TutorialLinAlgExSolveLDLT.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400121</tr>
122</table>
123
124For a \ref TopicLinearAlgebraDecompositions "much more complete table" comparing all decompositions supported by Eigen (notice that Eigen
Benoit Jacob76152e92010-06-29 10:02:33 -0400125supports many other decompositions), see our special page on
126\ref TopicLinearAlgebraDecompositions "this topic".
127
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400128\section TutorialLinAlgSolutionExists Checking if a solution really exists
Benoit Jacob76152e92010-06-29 10:02:33 -0400129
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400130Only 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
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200133<table class="example">
134<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400135<tr>
136 <td>\include TutorialLinAlgExComputeSolveError.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200137 <td>\verbinclude TutorialLinAlgExComputeSolveError.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400138</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
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200147<table class="example">
148<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400149<tr>
150 <td>\include TutorialLinAlgSelfAdjointEigenSolver.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200151 <td>\verbinclude TutorialLinAlgSelfAdjointEigenSolver.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400152</tr>
153</table>
154
Jitse Niesend0f6b1c2010-07-22 21:52:04 +0100155\section TutorialLinAlgInverse Computing inverse and determinant
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400156
157First of all, make sure that you really want this. While inverse and determinant are fundamental mathematical concepts,
158in \em numerical linear algebra they are not as popular as in pure mathematics. Inverse computations are often
159advantageously replaced by solve() operations, and the determinant is often \em not a good way of checking if a matrix
160is invertible.
161
162However, for \em very \em small matrices, the above is not true, and inverse and determinant can be very useful.
163
164While certain decompositions, such as PartialPivLU and FullPivLU, offer inverse() and determinant() methods, you can also
165call inverse() and determinant() directly on a matrix. If your matrix is of a very small fixed size (at most 4x4) this
166allows Eigen to avoid performing a LU decomposition, and instead use formulas that are more efficient on such small matrices.
167
168Here is an example:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200169<table class="example">
170<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400171<tr>
172 <td>\include TutorialLinAlgInverseDeterminant.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200173 <td>\verbinclude TutorialLinAlgInverseDeterminant.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400174</tr>
175</table>
176
177\section TutorialLinAlgLeastsquares Least squares solving
178
Benoit Jacob26129222010-10-15 09:44:43 -0400179The best way to do least squares solving is with a SVD decomposition. Eigen provides one as the JacobiSVD class, and its solve()
180is doing least-squares solving.
181
182Here is an example:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200183<table class="example">
184<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob26129222010-10-15 09:44:43 -0400185<tr>
186 <td>\include TutorialLinAlgSVDSolve.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200187 <td>\verbinclude TutorialLinAlgSVDSolve.out </td>
Benoit Jacob26129222010-10-15 09:44:43 -0400188</tr>
189</table>
190
191Another way, potentially faster but less reliable, is to use a LDLT decomposition
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400192of the normal matrix. In any case, just read any reference text on least squares, and it will be very easy for you
193to implement any linear least squares computation on top of Eigen.
194
195\section TutorialLinAlgSeparateComputation Separating the computation from the construction
196
197In the above examples, the decomposition was computed at the same time that the decomposition object was constructed.
198There are however situations where you might want to separate these two things, for example if you don't know,
199at the time of the construction, the matrix that you will want to decompose; or if you want to reuse an existing
200decomposition object.
201
202What makes this possible is that:
203\li all decompositions have a default constructor,
204\li all decompositions have a compute(matrix) method that does the computation, and that may be called again
205 on an already-computed decomposition, reinitializing it.
206
207For example:
208
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200209<table class="example">
210<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400211<tr>
212 <td>\include TutorialLinAlgComputeTwice.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200213 <td>\verbinclude TutorialLinAlgComputeTwice.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400214</tr>
215</table>
216
217Finally, you can tell the decomposition constructor to preallocate storage for decomposing matrices of a given size,
218so that when you subsequently decompose such matrices, no dynamic memory allocation is performed (of course, if you
219are using fixed-size matrices, no dynamic memory allocation happens at all). This is done by just
220passing the size to the decomposition constructor, as in this example:
221\code
222HouseholderQR<MatrixXf> qr(50,50);
223MatrixXf A = MatrixXf::Random(50,50);
224qr.compute(A); // no dynamic memory allocation
225\endcode
226
227\section TutorialLinAlgRankRevealing Rank-revealing decompositions
228
229Certain decompositions are rank-revealing, i.e. are able to compute the rank of a matrix. These are typically
230also the decompositions that behave best in the face of a non-full-rank matrix (which in the square case means a
231singular matrix). On \ref TopicLinearAlgebraDecompositions "this table" you can see for all our decompositions
232whether they are rank-revealing or not.
233
234Rank-revealing decompositions offer at least a rank() method. They can also offer convenience methods such as isInvertible(),
235and some are also providing methods to compute the kernel (null-space) and image (column-space) of the matrix, as is the
236case with FullPivLU:
237
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200238<table class="example">
239<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400240<tr>
241 <td>\include TutorialLinAlgRankRevealing.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200242 <td>\verbinclude TutorialLinAlgRankRevealing.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400243</tr>
244</table>
245
246Of course, any rank computation depends on the choice of an arbitrary threshold, since practically no
247floating-point matrix is \em exactly rank-deficient. Eigen picks a sensible default threshold, which depends
248on the decomposition but is typically the diagonal size times machine epsilon. While this is the best default we
249could pick, only you know what is the right threshold for your application. You can set this by calling setThreshold()
Benoit Jacob962b30d2010-06-30 19:27:30 -0400250on your decomposition object before calling rank() or any other method that needs to use such a threshold.
251The decomposition itself, i.e. the compute() method, is independent of the threshold. You don't need to recompute the
252decomposition after you've changed the threshold.
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400253
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200254<table class="example">
255<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400256<tr>
257 <td>\include TutorialLinAlgSetThreshold.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200258 <td>\verbinclude TutorialLinAlgSetThreshold.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400259</tr>
260</table>
261
Jitse Niesen26cfe5a2010-07-09 11:59:29 +0100262\li \b Next: \ref TutorialReductionsVisitorsBroadcasting
Benoit Jacob76152e92010-06-29 10:02:33 -0400263
264*/
265
266}