blob: c6235f2687f0499c27738fd7d236d07672fb04a6 [file] [log] [blame]
Benoit Jacob76152e92010-06-29 10:02:33 -04001namespace Eigen {
2
Gael Guennebaud93ee82b2013-01-05 16:37:11 +01003/** \eigenManualPage TutorialLinearAlgebra Linear algebra and decompositions
Benoit Jacob76152e92010-06-29 10:02:33 -04004
Benoit Jacob4d4a23c2010-06-30 10:11:55 -04005\li \b Previous: \ref TutorialAdvancedInitialization
Jitse Niesen26cfe5a2010-07-09 11:59:29 +01006\li \b Next: \ref TutorialReductionsVisitorsBroadcasting
Benoit Jacob76152e92010-06-29 10:02:33 -04007
8This tutorial explains how to solve linear systems, compute various decompositions such as LU,
Benoit Jacob4d4a23c2010-06-30 10:11:55 -04009QR, %SVD, eigendecompositions... for more advanced topics, don't miss our special page on
Benoit Jacob76152e92010-06-29 10:02:33 -040010\ref TopicLinearAlgebraDecompositions "this topic".
11
Gael Guennebaud93ee82b2013-01-05 16:37:11 +010012\eigenAutoToc
Jitse Niesend0f6b1c2010-07-22 21:52:04 +010013
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040014\section TutorialLinAlgBasicSolve Basic linear solving
Benoit Jacob76152e92010-06-29 10:02:33 -040015
16\b The \b problem: You have a system of equations, that you have written as a single matrix equation
17 \f[ Ax \: = \: b \f]
18Where \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.
19
20\b The \b solution: You can choose between various decompositions, depending on what your matrix \a A looks like,
21and depending on whether you favor speed or accuracy. However, let's start with an example that works in all cases,
22and is a good compromise:
Gael Guennebaudf66fe262010-10-19 11:40:49 +020023<table class="example">
24<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040025<tr>
26 <td>\include TutorialLinAlgExSolveColPivHouseholderQR.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020027 <td>\verbinclude TutorialLinAlgExSolveColPivHouseholderQR.out </td>
Benoit Jacob76152e92010-06-29 10:02:33 -040028</tr>
29</table>
30
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040031In this example, the colPivHouseholderQr() method returns an object of class ColPivHouseholderQR. Since here the
32matrix is of type Matrix3f, this line could have been replaced by:
Benoit Jacob76152e92010-06-29 10:02:33 -040033\code
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040034ColPivHouseholderQR<Matrix3f> dec(A);
Benoit Jacob76152e92010-06-29 10:02:33 -040035Vector3f x = dec.solve(b);
36\endcode
37
38Here, ColPivHouseholderQR is a QR decomposition with column pivoting. It's a good compromise for this tutorial, as it
39works for all matrices while being quite fast. Here is a table of some other decompositions that you can choose from,
40depending on your matrix and the trade-off you want to make:
41
Gael Guennebaudf66fe262010-10-19 11:40:49 +020042<table class="manual">
Benoit Jacob76152e92010-06-29 10:02:33 -040043 <tr>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020044 <th>Decomposition</th>
45 <th>Method</th>
46 <th>Requirements on the matrix</th>
47 <th>Speed</th>
48 <th>Accuracy</th>
Benoit Jacob76152e92010-06-29 10:02:33 -040049 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040050 <tr>
51 <td>PartialPivLU</td>
52 <td>partialPivLu()</td>
53 <td>Invertible</td>
54 <td>++</td>
55 <td>+</td>
56 </tr>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020057 <tr class="alt">
Benoit Jacob76152e92010-06-29 10:02:33 -040058 <td>FullPivLU</td>
59 <td>fullPivLu()</td>
60 <td>None</td>
61 <td>-</td>
62 <td>+++</td>
63 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040064 <tr>
65 <td>HouseholderQR</td>
66 <td>householderQr()</td>
67 <td>None</td>
68 <td>++</td>
69 <td>+</td>
70 </tr>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020071 <tr class="alt">
Benoit Jacob76152e92010-06-29 10:02:33 -040072 <td>ColPivHouseholderQR</td>
73 <td>colPivHouseholderQr()</td>
74 <td>None</td>
75 <td>+</td>
76 <td>++</td>
77 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040078 <tr>
79 <td>FullPivHouseholderQR</td>
80 <td>fullPivHouseholderQr()</td>
81 <td>None</td>
82 <td>-</td>
83 <td>+++</td>
84 </tr>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020085 <tr class="alt">
Benoit Jacob76152e92010-06-29 10:02:33 -040086 <td>LLT</td>
87 <td>llt()</td>
88 <td>Positive definite</td>
89 <td>+++</td>
90 <td>+</td>
91 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040092 <tr>
93 <td>LDLT</td>
94 <td>ldlt()</td>
95 <td>Positive or negative semidefinite</td>
96 <td>+++</td>
97 <td>++</td>
98 </tr>
Benoit Jacob76152e92010-06-29 10:02:33 -040099</table>
100
101All of these decompositions offer a solve() method that works as in the above example.
102
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400103For example, if your matrix is positive definite, the above table says that a very good
104choice is then the LDLT decomposition. Here's an example, also demonstrating that using a general
105matrix (not a vector) as right hand side is possible.
106
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200107<table class="example">
108<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400109<tr>
110 <td>\include TutorialLinAlgExSolveLDLT.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200111 <td>\verbinclude TutorialLinAlgExSolveLDLT.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400112</tr>
113</table>
114
115For a \ref TopicLinearAlgebraDecompositions "much more complete table" comparing all decompositions supported by Eigen (notice that Eigen
Benoit Jacob76152e92010-06-29 10:02:33 -0400116supports many other decompositions), see our special page on
117\ref TopicLinearAlgebraDecompositions "this topic".
118
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400119\section TutorialLinAlgSolutionExists Checking if a solution really exists
Benoit Jacob76152e92010-06-29 10:02:33 -0400120
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400121Only you know what error margin you want to allow for a solution to be considered valid.
122So Eigen lets you do this computation for yourself, if you want to, as in this example:
123
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200124<table class="example">
125<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400126<tr>
127 <td>\include TutorialLinAlgExComputeSolveError.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200128 <td>\verbinclude TutorialLinAlgExComputeSolveError.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400129</tr>
130</table>
131
132\section TutorialLinAlgEigensolving Computing eigenvalues and eigenvectors
133
134You need an eigendecomposition here, see available such decompositions on \ref TopicLinearAlgebraDecompositions "this page".
135Make sure to check if your matrix is self-adjoint, as is often the case in these problems. Here's an example using
136SelfAdjointEigenSolver, it could easily be adapted to general matrices using EigenSolver or ComplexEigenSolver.
137
Jitse Niesen45a6bb32011-11-07 17:14:06 +0000138The computation of eigenvalues and eigenvectors does not necessarily converge, but such failure to converge is
139very rare. The call to info() is to check for this possibility.
140
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200141<table class="example">
142<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400143<tr>
144 <td>\include TutorialLinAlgSelfAdjointEigenSolver.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200145 <td>\verbinclude TutorialLinAlgSelfAdjointEigenSolver.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400146</tr>
147</table>
148
Jitse Niesend0f6b1c2010-07-22 21:52:04 +0100149\section TutorialLinAlgInverse Computing inverse and determinant
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400150
151First of all, make sure that you really want this. While inverse and determinant are fundamental mathematical concepts,
152in \em numerical linear algebra they are not as popular as in pure mathematics. Inverse computations are often
153advantageously replaced by solve() operations, and the determinant is often \em not a good way of checking if a matrix
154is invertible.
155
156However, for \em very \em small matrices, the above is not true, and inverse and determinant can be very useful.
157
158While certain decompositions, such as PartialPivLU and FullPivLU, offer inverse() and determinant() methods, you can also
159call inverse() and determinant() directly on a matrix. If your matrix is of a very small fixed size (at most 4x4) this
160allows Eigen to avoid performing a LU decomposition, and instead use formulas that are more efficient on such small matrices.
161
162Here is an example:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200163<table class="example">
164<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400165<tr>
166 <td>\include TutorialLinAlgInverseDeterminant.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200167 <td>\verbinclude TutorialLinAlgInverseDeterminant.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400168</tr>
169</table>
170
171\section TutorialLinAlgLeastsquares Least squares solving
172
Benoit Jacob26129222010-10-15 09:44:43 -0400173The best way to do least squares solving is with a SVD decomposition. Eigen provides one as the JacobiSVD class, and its solve()
174is doing least-squares solving.
175
176Here is an example:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200177<table class="example">
178<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob26129222010-10-15 09:44:43 -0400179<tr>
180 <td>\include TutorialLinAlgSVDSolve.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200181 <td>\verbinclude TutorialLinAlgSVDSolve.out </td>
Benoit Jacob26129222010-10-15 09:44:43 -0400182</tr>
183</table>
184
185Another way, potentially faster but less reliable, is to use a LDLT decomposition
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400186of the normal matrix. In any case, just read any reference text on least squares, and it will be very easy for you
187to implement any linear least squares computation on top of Eigen.
188
189\section TutorialLinAlgSeparateComputation Separating the computation from the construction
190
191In the above examples, the decomposition was computed at the same time that the decomposition object was constructed.
192There are however situations where you might want to separate these two things, for example if you don't know,
193at the time of the construction, the matrix that you will want to decompose; or if you want to reuse an existing
194decomposition object.
195
196What makes this possible is that:
197\li all decompositions have a default constructor,
198\li all decompositions have a compute(matrix) method that does the computation, and that may be called again
199 on an already-computed decomposition, reinitializing it.
200
201For example:
202
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200203<table class="example">
204<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400205<tr>
206 <td>\include TutorialLinAlgComputeTwice.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200207 <td>\verbinclude TutorialLinAlgComputeTwice.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400208</tr>
209</table>
210
211Finally, you can tell the decomposition constructor to preallocate storage for decomposing matrices of a given size,
212so that when you subsequently decompose such matrices, no dynamic memory allocation is performed (of course, if you
213are using fixed-size matrices, no dynamic memory allocation happens at all). This is done by just
214passing the size to the decomposition constructor, as in this example:
215\code
216HouseholderQR<MatrixXf> qr(50,50);
217MatrixXf A = MatrixXf::Random(50,50);
218qr.compute(A); // no dynamic memory allocation
219\endcode
220
221\section TutorialLinAlgRankRevealing Rank-revealing decompositions
222
223Certain decompositions are rank-revealing, i.e. are able to compute the rank of a matrix. These are typically
224also the decompositions that behave best in the face of a non-full-rank matrix (which in the square case means a
225singular matrix). On \ref TopicLinearAlgebraDecompositions "this table" you can see for all our decompositions
226whether they are rank-revealing or not.
227
228Rank-revealing decompositions offer at least a rank() method. They can also offer convenience methods such as isInvertible(),
229and some are also providing methods to compute the kernel (null-space) and image (column-space) of the matrix, as is the
230case with FullPivLU:
231
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200232<table class="example">
233<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400234<tr>
235 <td>\include TutorialLinAlgRankRevealing.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200236 <td>\verbinclude TutorialLinAlgRankRevealing.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400237</tr>
238</table>
239
240Of course, any rank computation depends on the choice of an arbitrary threshold, since practically no
241floating-point matrix is \em exactly rank-deficient. Eigen picks a sensible default threshold, which depends
242on the decomposition but is typically the diagonal size times machine epsilon. While this is the best default we
243could 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 -0400244on your decomposition object before calling rank() or any other method that needs to use such a threshold.
245The decomposition itself, i.e. the compute() method, is independent of the threshold. You don't need to recompute the
246decomposition after you've changed the threshold.
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400247
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200248<table class="example">
249<tr><th>Example:</th><th>Output:</th></tr>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400250<tr>
251 <td>\include TutorialLinAlgSetThreshold.cpp </td>
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200252 <td>\verbinclude TutorialLinAlgSetThreshold.out </td>
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400253</tr>
254</table>
255
Jitse Niesen26cfe5a2010-07-09 11:59:29 +0100256\li \b Next: \ref TutorialReductionsVisitorsBroadcasting
Benoit Jacob76152e92010-06-29 10:02:33 -0400257
258*/
259
260}