blob: c471a6d04cab2a02d46fa065ea6312e46ce9521d [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
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040013\section TutorialLinAlgBasicSolve Basic linear solving
Benoit Jacob76152e92010-06-29 10:02:33 -040014
15\b The \b problem: You have a system of equations, that you have written as a single matrix equation
16 \f[ Ax \: = \: b \f]
17Where \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.
18
19\b The \b solution: You can choose between various decompositions, depending on what your matrix \a A looks like,
20and depending on whether you favor speed or accuracy. However, let's start with an example that works in all cases,
21and is a good compromise:
22<table class="tutorial_code">
23<tr>
24 <td>\include TutorialLinAlgExSolveColPivHouseholderQR.cpp </td>
25 <td>output: \verbinclude TutorialLinAlgExSolveColPivHouseholderQR.out </td>
26</tr>
27</table>
28
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040029In this example, the colPivHouseholderQr() method returns an object of class ColPivHouseholderQR. Since here the
30matrix is of type Matrix3f, this line could have been replaced by:
Benoit Jacob76152e92010-06-29 10:02:33 -040031\code
Benoit Jacob4d4a23c2010-06-30 10:11:55 -040032ColPivHouseholderQR<Matrix3f> dec(A);
Benoit Jacob76152e92010-06-29 10:02:33 -040033Vector3f x = dec.solve(b);
34\endcode
35
36Here, ColPivHouseholderQR is a QR decomposition with column pivoting. It's a good compromise for this tutorial, as it
37works for all matrices while being quite fast. Here is a table of some other decompositions that you can choose from,
38depending on your matrix and the trade-off you want to make:
39
40<table border="1">
41
42 <tr>
43 <td>Decomposition</td>
44 <td>Method</td>
45 <td>Requirements on the matrix</td>
46 <td>Speed</td>
47 <td>Accuracy</td>
48 </tr>
49
50 <tr>
51 <td>PartialPivLU</td>
52 <td>partialPivLu()</td>
53 <td>Invertible</td>
54 <td>++</td>
55 <td>+</td>
56 </tr>
57
58 <tr>
59 <td>FullPivLU</td>
60 <td>fullPivLu()</td>
61 <td>None</td>
62 <td>-</td>
63 <td>+++</td>
64 </tr>
65
66 <tr>
67 <td>HouseholderQR</td>
68 <td>householderQr()</td>
69 <td>None</td>
70 <td>++</td>
71 <td>+</td>
72 </tr>
73
74 <tr>
75 <td>ColPivHouseholderQR</td>
76 <td>colPivHouseholderQr()</td>
77 <td>None</td>
78 <td>+</td>
79 <td>++</td>
80 </tr>
81
82 <tr>
83 <td>FullPivHouseholderQR</td>
84 <td>fullPivHouseholderQr()</td>
85 <td>None</td>
86 <td>-</td>
87 <td>+++</td>
88 </tr>
89
90 <tr>
91 <td>LLT</td>
92 <td>llt()</td>
93 <td>Positive definite</td>
94 <td>+++</td>
95 <td>+</td>
96 </tr>
97
98 <tr>
99 <td>LDLT</td>
100 <td>ldlt()</td>
101 <td>Positive or negative semidefinite</td>
102 <td>+++</td>
103 <td>++</td>
104 </tr>
105
106</table>
107
108All of these decompositions offer a solve() method that works as in the above example.
109
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400110For example, if your matrix is positive definite, the above table says that a very good
111choice is then the LDLT decomposition. Here's an example, also demonstrating that using a general
112matrix (not a vector) as right hand side is possible.
113
114<table class="tutorial_code">
115<tr>
116 <td>\include TutorialLinAlgExSolveLDLT.cpp </td>
117 <td>output: \verbinclude TutorialLinAlgExSolveLDLT.out </td>
118</tr>
119</table>
120
121For a \ref TopicLinearAlgebraDecompositions "much more complete table" comparing all decompositions supported by Eigen (notice that Eigen
Benoit Jacob76152e92010-06-29 10:02:33 -0400122supports many other decompositions), see our special page on
123\ref TopicLinearAlgebraDecompositions "this topic".
124
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400125\section TutorialLinAlgSolutionExists Checking if a solution really exists
Benoit Jacob76152e92010-06-29 10:02:33 -0400126
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400127Only you know what error margin you want to allow for a solution to be considered valid.
128So Eigen lets you do this computation for yourself, if you want to, as in this example:
129
130<table class="tutorial_code">
131<tr>
132 <td>\include TutorialLinAlgExComputeSolveError.cpp </td>
133 <td>output: \verbinclude TutorialLinAlgExComputeSolveError.out </td>
134</tr>
135</table>
136
137\section TutorialLinAlgEigensolving Computing eigenvalues and eigenvectors
138
139You need an eigendecomposition here, see available such decompositions on \ref TopicLinearAlgebraDecompositions "this page".
140Make sure to check if your matrix is self-adjoint, as is often the case in these problems. Here's an example using
141SelfAdjointEigenSolver, it could easily be adapted to general matrices using EigenSolver or ComplexEigenSolver.
142
143<table class="tutorial_code">
144<tr>
145 <td>\include TutorialLinAlgSelfAdjointEigenSolver.cpp </td>
146 <td>output: \verbinclude TutorialLinAlgSelfAdjointEigenSolver.out </td>
147</tr>
148</table>
149
150\section TutorialLinAlgEigensolving Computing inverse and determinant
151
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:
164<table class="tutorial_code">
165<tr>
166 <td>\include TutorialLinAlgInverseDeterminant.cpp </td>
167 <td>output: \verbinclude TutorialLinAlgInverseDeterminant.out </td>
168</tr>
169</table>
170
171\section TutorialLinAlgLeastsquares Least squares solving
172
173Eigen doesn't currently provide built-in linear least squares solving functions, but you can easily compute that yourself
174from Eigen's decompositions. The most reliable way is to use a SVD (or better yet, JacobiSVD), and in the future
175these classes will offer methods for least squares solving. Another, potentially faster way, is to use a LLT decomposition
176of the normal matrix. In any case, just read any reference text on least squares, and it will be very easy for you
177to implement any linear least squares computation on top of Eigen.
178
179\section TutorialLinAlgSeparateComputation Separating the computation from the construction
180
181In the above examples, the decomposition was computed at the same time that the decomposition object was constructed.
182There are however situations where you might want to separate these two things, for example if you don't know,
183at the time of the construction, the matrix that you will want to decompose; or if you want to reuse an existing
184decomposition object.
185
186What makes this possible is that:
187\li all decompositions have a default constructor,
188\li all decompositions have a compute(matrix) method that does the computation, and that may be called again
189 on an already-computed decomposition, reinitializing it.
190
191For example:
192
193<table class="tutorial_code">
194<tr>
195 <td>\include TutorialLinAlgComputeTwice.cpp </td>
196 <td>output: \verbinclude TutorialLinAlgComputeTwice.out </td>
197</tr>
198</table>
199
200Finally, you can tell the decomposition constructor to preallocate storage for decomposing matrices of a given size,
201so that when you subsequently decompose such matrices, no dynamic memory allocation is performed (of course, if you
202are using fixed-size matrices, no dynamic memory allocation happens at all). This is done by just
203passing the size to the decomposition constructor, as in this example:
204\code
205HouseholderQR<MatrixXf> qr(50,50);
206MatrixXf A = MatrixXf::Random(50,50);
207qr.compute(A); // no dynamic memory allocation
208\endcode
209
210\section TutorialLinAlgRankRevealing Rank-revealing decompositions
211
212Certain decompositions are rank-revealing, i.e. are able to compute the rank of a matrix. These are typically
213also the decompositions that behave best in the face of a non-full-rank matrix (which in the square case means a
214singular matrix). On \ref TopicLinearAlgebraDecompositions "this table" you can see for all our decompositions
215whether they are rank-revealing or not.
216
217Rank-revealing decompositions offer at least a rank() method. They can also offer convenience methods such as isInvertible(),
218and some are also providing methods to compute the kernel (null-space) and image (column-space) of the matrix, as is the
219case with FullPivLU:
220
221<table class="tutorial_code">
222<tr>
223 <td>\include TutorialLinAlgRankRevealing.cpp </td>
224 <td>output: \verbinclude TutorialLinAlgRankRevealing.out </td>
225</tr>
226</table>
227
228Of course, any rank computation depends on the choice of an arbitrary threshold, since practically no
229floating-point matrix is \em exactly rank-deficient. Eigen picks a sensible default threshold, which depends
230on the decomposition but is typically the diagonal size times machine epsilon. While this is the best default we
231could 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 -0400232on your decomposition object before calling rank() or any other method that needs to use such a threshold.
233The decomposition itself, i.e. the compute() method, is independent of the threshold. You don't need to recompute the
234decomposition after you've changed the threshold.
Benoit Jacob4d4a23c2010-06-30 10:11:55 -0400235
236<table class="tutorial_code">
237<tr>
238 <td>\include TutorialLinAlgSetThreshold.cpp </td>
239 <td>output: \verbinclude TutorialLinAlgSetThreshold.out </td>
240</tr>
241</table>
242
Jitse Niesen26cfe5a2010-07-09 11:59:29 +0100243\li \b Next: \ref TutorialReductionsVisitorsBroadcasting
Benoit Jacob76152e92010-06-29 10:02:33 -0400244
245*/
246
247}