blob: c8f2bf23d64fba25409afbebbfc422245e9e5537 [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:
32<table class="tutorial_code">
33<tr>
34 <td>\include TutorialLinAlgExSolveColPivHouseholderQR.cpp </td>
35 <td>output: \verbinclude TutorialLinAlgExSolveColPivHouseholderQR.out </td>
36</tr>
37</table>
38
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040039In this example, the colPivHouseholderQr() method returns an object of class ColPivHouseholderQR. Since here the
40matrix is of type Matrix3f, this line could have been replaced by:
Benoit Jacob76152e92010-06-29 10:02:33 -040041\code
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040042ColPivHouseholderQR<Matrix3f> dec(A);
Benoit Jacob76152e92010-06-29 10:02:33 -040043Vector3f x = dec.solve(b);
44\endcode
45
46Here, ColPivHouseholderQR is a QR decomposition with column pivoting. It's a good compromise for this tutorial, as it
47works for all matrices while being quite fast. Here is a table of some other decompositions that you can choose from,
48depending on your matrix and the trade-off you want to make:
49
50<table border="1">
51
52 <tr>
53 <td>Decomposition</td>
54 <td>Method</td>
55 <td>Requirements on the matrix</td>
56 <td>Speed</td>
57 <td>Accuracy</td>
58 </tr>
59
60 <tr>
61 <td>PartialPivLU</td>
62 <td>partialPivLu()</td>
63 <td>Invertible</td>
64 <td>++</td>
65 <td>+</td>
66 </tr>
67
68 <tr>
69 <td>FullPivLU</td>
70 <td>fullPivLu()</td>
71 <td>None</td>
72 <td>-</td>
73 <td>+++</td>
74 </tr>
75
76 <tr>
77 <td>HouseholderQR</td>
78 <td>householderQr()</td>
79 <td>None</td>
80 <td>++</td>
81 <td>+</td>
82 </tr>
83
84 <tr>
85 <td>ColPivHouseholderQR</td>
86 <td>colPivHouseholderQr()</td>
87 <td>None</td>
88 <td>+</td>
89 <td>++</td>
90 </tr>
91
92 <tr>
93 <td>FullPivHouseholderQR</td>
94 <td>fullPivHouseholderQr()</td>
95 <td>None</td>
96 <td>-</td>
97 <td>+++</td>
98 </tr>
99
100 <tr>
101 <td>LLT</td>
102 <td>llt()</td>
103 <td>Positive definite</td>
104 <td>+++</td>
105 <td>+</td>
106 </tr>
107
108 <tr>
109 <td>LDLT</td>
110 <td>ldlt()</td>
111 <td>Positive or negative semidefinite</td>
112 <td>+++</td>
113 <td>++</td>
114 </tr>
115
116</table>
117
118All of these decompositions offer a solve() method that works as in the above example.
119
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400120For example, if your matrix is positive definite, the above table says that a very good
121choice is then the LDLT decomposition. Here's an example, also demonstrating that using a general
122matrix (not a vector) as right hand side is possible.
123
124<table class="tutorial_code">
125<tr>
126 <td>\include TutorialLinAlgExSolveLDLT.cpp </td>
127 <td>output: \verbinclude TutorialLinAlgExSolveLDLT.out </td>
128</tr>
129</table>
130
131For a \ref TopicLinearAlgebraDecompositions "much more complete table" comparing all decompositions supported by Eigen (notice that Eigen
Benoit Jacob76152e92010-06-29 10:02:33 -0400132supports many other decompositions), see our special page on
133\ref TopicLinearAlgebraDecompositions "this topic".
134
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400135\section TutorialLinAlgSolutionExists Checking if a solution really exists
Benoit Jacob76152e92010-06-29 10:02:33 -0400136
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400137Only you know what error margin you want to allow for a solution to be considered valid.
138So Eigen lets you do this computation for yourself, if you want to, as in this example:
139
140<table class="tutorial_code">
141<tr>
142 <td>\include TutorialLinAlgExComputeSolveError.cpp </td>
143 <td>output: \verbinclude TutorialLinAlgExComputeSolveError.out </td>
144</tr>
145</table>
146
147\section TutorialLinAlgEigensolving Computing eigenvalues and eigenvectors
148
149You need an eigendecomposition here, see available such decompositions on \ref TopicLinearAlgebraDecompositions "this page".
150Make sure to check if your matrix is self-adjoint, as is often the case in these problems. Here's an example using
151SelfAdjointEigenSolver, it could easily be adapted to general matrices using EigenSolver or ComplexEigenSolver.
152
153<table class="tutorial_code">
154<tr>
155 <td>\include TutorialLinAlgSelfAdjointEigenSolver.cpp </td>
156 <td>output: \verbinclude TutorialLinAlgSelfAdjointEigenSolver.out </td>
157</tr>
158</table>
159
Jitse Niesend0f6b1c2010-07-22 21:52:04 +0100160\section TutorialLinAlgInverse Computing inverse and determinant
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400161
162First of all, make sure that you really want this. While inverse and determinant are fundamental mathematical concepts,
163in \em numerical linear algebra they are not as popular as in pure mathematics. Inverse computations are often
164advantageously replaced by solve() operations, and the determinant is often \em not a good way of checking if a matrix
165is invertible.
166
167However, for \em very \em small matrices, the above is not true, and inverse and determinant can be very useful.
168
169While certain decompositions, such as PartialPivLU and FullPivLU, offer inverse() and determinant() methods, you can also
170call inverse() and determinant() directly on a matrix. If your matrix is of a very small fixed size (at most 4x4) this
171allows Eigen to avoid performing a LU decomposition, and instead use formulas that are more efficient on such small matrices.
172
173Here is an example:
174<table class="tutorial_code">
175<tr>
176 <td>\include TutorialLinAlgInverseDeterminant.cpp </td>
177 <td>output: \verbinclude TutorialLinAlgInverseDeterminant.out </td>
178</tr>
179</table>
180
181\section TutorialLinAlgLeastsquares Least squares solving
182
Benoit Jacob26129222010-10-15 09:44:43 -0400183The best way to do least squares solving is with a SVD decomposition. Eigen provides one as the JacobiSVD class, and its solve()
184is doing least-squares solving.
185
186Here is an example:
187<table class="tutorial_code">
188<tr>
189 <td>\include TutorialLinAlgSVDSolve.cpp </td>
190 <td>output: \verbinclude TutorialLinAlgSVDSolve.out </td>
191</tr>
192</table>
193
194Another way, potentially faster but less reliable, is to use a LDLT decomposition
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400195of 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="tutorial_code">
213<tr>
214 <td>\include TutorialLinAlgComputeTwice.cpp </td>
215 <td>output: \verbinclude TutorialLinAlgComputeTwice.out </td>
216</tr>
217</table>
218
219Finally, you can tell the decomposition constructor to preallocate storage for decomposing matrices of a given size,
220so that when you subsequently decompose such matrices, no dynamic memory allocation is performed (of course, if you
221are using fixed-size matrices, no dynamic memory allocation happens at all). This is done by just
222passing the size to the decomposition constructor, as in this example:
223\code
224HouseholderQR<MatrixXf> qr(50,50);
225MatrixXf A = MatrixXf::Random(50,50);
226qr.compute(A); // no dynamic memory allocation
227\endcode
228
229\section TutorialLinAlgRankRevealing Rank-revealing decompositions
230
231Certain decompositions are rank-revealing, i.e. are able to compute the rank of a matrix. These are typically
232also the decompositions that behave best in the face of a non-full-rank matrix (which in the square case means a
233singular matrix). On \ref TopicLinearAlgebraDecompositions "this table" you can see for all our decompositions
234whether they are rank-revealing or not.
235
236Rank-revealing decompositions offer at least a rank() method. They can also offer convenience methods such as isInvertible(),
237and some are also providing methods to compute the kernel (null-space) and image (column-space) of the matrix, as is the
238case with FullPivLU:
239
240<table class="tutorial_code">
241<tr>
242 <td>\include TutorialLinAlgRankRevealing.cpp </td>
243 <td>output: \verbinclude TutorialLinAlgRankRevealing.out </td>
244</tr>
245</table>
246
247Of course, any rank computation depends on the choice of an arbitrary threshold, since practically no
248floating-point matrix is \em exactly rank-deficient. Eigen picks a sensible default threshold, which depends
249on the decomposition but is typically the diagonal size times machine epsilon. While this is the best default we
250could 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 -0400251on your decomposition object before calling rank() or any other method that needs to use such a threshold.
252The decomposition itself, i.e. the compute() method, is independent of the threshold. You don't need to recompute the
253decomposition after you've changed the threshold.
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400254
255<table class="tutorial_code">
256<tr>
257 <td>\include TutorialLinAlgSetThreshold.cpp </td>
258 <td>output: \verbinclude TutorialLinAlgSetThreshold.out </td>
259</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}