blob: c1218a562bddf8b0dbf3141f7f330bf44d9ce989 [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
Gael Guennebaud2ea1e492012-12-28 18:58:07 +010013\tableofcontents
Jitse Niesend0f6b1c2010-07-22 21:52:04 +010014
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040015\section TutorialLinAlgBasicSolve Basic linear solving
Benoit Jacob76152e92010-06-29 10:02:33 -040016
17\b The \b problem: You have a system of equations, that you have written as a single matrix equation
18 \f[ Ax \: = \: b \f]
19Where \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.
20
21\b The \b solution: You can choose between various decompositions, depending on what your matrix \a A looks like,
22and depending on whether you favor speed or accuracy. However, let's start with an example that works in all cases,
23and is a good compromise:
Gael Guennebaudf66fe262010-10-19 11:40:49 +020024<table class="example">
25<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040026<tr>
27 <td>\include TutorialLinAlgExSolveColPivHouseholderQR.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020028 <td>\verbinclude TutorialLinAlgExSolveColPivHouseholderQR.out </td>
Benoit Jacob76152e92010-06-29 10:02:33 -040029</tr>
30</table>
31
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040032In this example, the colPivHouseholderQr() method returns an object of class ColPivHouseholderQR. Since here the
33matrix is of type Matrix3f, this line could have been replaced by:
Benoit Jacob76152e92010-06-29 10:02:33 -040034\code
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040035ColPivHouseholderQR<Matrix3f> dec(A);
Benoit Jacob76152e92010-06-29 10:02:33 -040036Vector3f x = dec.solve(b);
37\endcode
38
39Here, ColPivHouseholderQR is a QR decomposition with column pivoting. It's a good compromise for this tutorial, as it
40works for all matrices while being quite fast. Here is a table of some other decompositions that you can choose from,
41depending on your matrix and the trade-off you want to make:
42
Gael Guennebaudf66fe262010-10-19 11:40:49 +020043<table class="manual">
Benoit Jacob76152e92010-06-29 10:02:33 -040044 <tr>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020045 <th>Decomposition</th>
46 <th>Method</th>
47 <th>Requirements on the matrix</th>
48 <th>Speed</th>
49 <th>Accuracy</th>
Benoit Jacob76152e92010-06-29 10:02:33 -040050 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040051 <tr>
52 <td>PartialPivLU</td>
53 <td>partialPivLu()</td>
54 <td>Invertible</td>
55 <td>++</td>
56 <td>+</td>
57 </tr>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020058 <tr class="alt">
Benoit Jacob76152e92010-06-29 10:02:33 -040059 <td>FullPivLU</td>
60 <td>fullPivLu()</td>
61 <td>None</td>
62 <td>-</td>
63 <td>+++</td>
64 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040065 <tr>
66 <td>HouseholderQR</td>
67 <td>householderQr()</td>
68 <td>None</td>
69 <td>++</td>
70 <td>+</td>
71 </tr>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020072 <tr class="alt">
Benoit Jacob76152e92010-06-29 10:02:33 -040073 <td>ColPivHouseholderQR</td>
74 <td>colPivHouseholderQr()</td>
75 <td>None</td>
76 <td>+</td>
77 <td>++</td>
78 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040079 <tr>
80 <td>FullPivHouseholderQR</td>
81 <td>fullPivHouseholderQr()</td>
82 <td>None</td>
83 <td>-</td>
84 <td>+++</td>
85 </tr>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020086 <tr class="alt">
Benoit Jacob76152e92010-06-29 10:02:33 -040087 <td>LLT</td>
88 <td>llt()</td>
89 <td>Positive definite</td>
90 <td>+++</td>
91 <td>+</td>
92 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040093 <tr>
94 <td>LDLT</td>
95 <td>ldlt()</td>
96 <td>Positive or negative semidefinite</td>
97 <td>+++</td>
98 <td>++</td>
99 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -0400100</table>
101
102All of these decompositions offer a solve() method that works as in the above example.
103
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400104For example, if your matrix is positive definite, the above table says that a very good
105choice is then the LDLT decomposition. Here's an example, also demonstrating that using a general
106matrix (not a vector) as right hand side is possible.
107
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200108<table class="example">
109<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400110<tr>
111 <td>\include TutorialLinAlgExSolveLDLT.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200112 <td>\verbinclude TutorialLinAlgExSolveLDLT.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400113</tr>
114</table>
115
116For a \ref TopicLinearAlgebraDecompositions "much more complete table" comparing all decompositions supported by Eigen (notice that Eigen
Benoit Jacob76152e92010-06-29 10:02:33 -0400117supports many other decompositions), see our special page on
118\ref TopicLinearAlgebraDecompositions "this topic".
119
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400120\section TutorialLinAlgSolutionExists Checking if a solution really exists
Benoit Jacob76152e92010-06-29 10:02:33 -0400121
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400122Only you know what error margin you want to allow for a solution to be considered valid.
123So Eigen lets you do this computation for yourself, if you want to, as in this example:
124
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200125<table class="example">
126<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400127<tr>
128 <td>\include TutorialLinAlgExComputeSolveError.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200129 <td>\verbinclude TutorialLinAlgExComputeSolveError.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400130</tr>
131</table>
132
133\section TutorialLinAlgEigensolving Computing eigenvalues and eigenvectors
134
135You need an eigendecomposition here, see available such decompositions on \ref TopicLinearAlgebraDecompositions "this page".
136Make sure to check if your matrix is self-adjoint, as is often the case in these problems. Here's an example using
137SelfAdjointEigenSolver, it could easily be adapted to general matrices using EigenSolver or ComplexEigenSolver.
138
Jitse Niesen45a6bb32011-11-07 17:14:06 +0000139The computation of eigenvalues and eigenvectors does not necessarily converge, but such failure to converge is
140very rare. The call to info() is to check for this possibility.
141
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200142<table class="example">
143<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400144<tr>
145 <td>\include TutorialLinAlgSelfAdjointEigenSolver.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200146 <td>\verbinclude TutorialLinAlgSelfAdjointEigenSolver.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400147</tr>
148</table>
149
Jitse Niesend0f6b1c2010-07-22 21:52:04 +0100150\section TutorialLinAlgInverse Computing inverse and determinant
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400151
152First of all, make sure that you really want this. While inverse and determinant are fundamental mathematical concepts,
153in \em numerical linear algebra they are not as popular as in pure mathematics. Inverse computations are often
154advantageously replaced by solve() operations, and the determinant is often \em not a good way of checking if a matrix
155is invertible.
156
157However, for \em very \em small matrices, the above is not true, and inverse and determinant can be very useful.
158
159While certain decompositions, such as PartialPivLU and FullPivLU, offer inverse() and determinant() methods, you can also
160call inverse() and determinant() directly on a matrix. If your matrix is of a very small fixed size (at most 4x4) this
161allows Eigen to avoid performing a LU decomposition, and instead use formulas that are more efficient on such small matrices.
162
163Here is an example:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200164<table class="example">
165<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400166<tr>
167 <td>\include TutorialLinAlgInverseDeterminant.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200168 <td>\verbinclude TutorialLinAlgInverseDeterminant.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400169</tr>
170</table>
171
172\section TutorialLinAlgLeastsquares Least squares solving
173
Benoit Jacob26129222010-10-15 09:44:43 -0400174The best way to do least squares solving is with a SVD decomposition. Eigen provides one as the JacobiSVD class, and its solve()
175is doing least-squares solving.
176
177Here is an example:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200178<table class="example">
179<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob26129222010-10-15 09:44:43 -0400180<tr>
181 <td>\include TutorialLinAlgSVDSolve.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200182 <td>\verbinclude TutorialLinAlgSVDSolve.out </td>
Benoit Jacob26129222010-10-15 09:44:43 -0400183</tr>
184</table>
185
186Another way, potentially faster but less reliable, is to use a LDLT decomposition
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400187of the normal matrix. In any case, just read any reference text on least squares, and it will be very easy for you
188to implement any linear least squares computation on top of Eigen.
189
190\section TutorialLinAlgSeparateComputation Separating the computation from the construction
191
192In the above examples, the decomposition was computed at the same time that the decomposition object was constructed.
193There are however situations where you might want to separate these two things, for example if you don't know,
194at the time of the construction, the matrix that you will want to decompose; or if you want to reuse an existing
195decomposition object.
196
197What makes this possible is that:
198\li all decompositions have a default constructor,
199\li all decompositions have a compute(matrix) method that does the computation, and that may be called again
200 on an already-computed decomposition, reinitializing it.
201
202For example:
203
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200204<table class="example">
205<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400206<tr>
207 <td>\include TutorialLinAlgComputeTwice.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200208 <td>\verbinclude TutorialLinAlgComputeTwice.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400209</tr>
210</table>
211
212Finally, you can tell the decomposition constructor to preallocate storage for decomposing matrices of a given size,
213so that when you subsequently decompose such matrices, no dynamic memory allocation is performed (of course, if you
214are using fixed-size matrices, no dynamic memory allocation happens at all). This is done by just
215passing the size to the decomposition constructor, as in this example:
216\code
217HouseholderQR<MatrixXf> qr(50,50);
218MatrixXf A = MatrixXf::Random(50,50);
219qr.compute(A); // no dynamic memory allocation
220\endcode
221
222\section TutorialLinAlgRankRevealing Rank-revealing decompositions
223
224Certain decompositions are rank-revealing, i.e. are able to compute the rank of a matrix. These are typically
225also the decompositions that behave best in the face of a non-full-rank matrix (which in the square case means a
226singular matrix). On \ref TopicLinearAlgebraDecompositions "this table" you can see for all our decompositions
227whether they are rank-revealing or not.
228
229Rank-revealing decompositions offer at least a rank() method. They can also offer convenience methods such as isInvertible(),
230and some are also providing methods to compute the kernel (null-space) and image (column-space) of the matrix, as is the
231case with FullPivLU:
232
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200233<table class="example">
234<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400235<tr>
236 <td>\include TutorialLinAlgRankRevealing.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200237 <td>\verbinclude TutorialLinAlgRankRevealing.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400238</tr>
239</table>
240
241Of course, any rank computation depends on the choice of an arbitrary threshold, since practically no
242floating-point matrix is \em exactly rank-deficient. Eigen picks a sensible default threshold, which depends
243on the decomposition but is typically the diagonal size times machine epsilon. While this is the best default we
244could 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 -0400245on your decomposition object before calling rank() or any other method that needs to use such a threshold.
246The decomposition itself, i.e. the compute() method, is independent of the threshold. You don't need to recompute the
247decomposition after you've changed the threshold.
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400248
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200249<table class="example">
250<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400251<tr>
252 <td>\include TutorialLinAlgSetThreshold.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200253 <td>\verbinclude TutorialLinAlgSetThreshold.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400254</tr>
255</table>
256
Jitse Niesen26cfe5a2010-07-09 11:59:29 +0100257\li \b Next: \ref TutorialReductionsVisitorsBroadcasting
Benoit Jacob76152e92010-06-29 10:02:33 -0400258
259*/
260
261}