blob: f31fcfc0f630fe9e3e2ce1b6fe31c3a53c42f9b0 [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
6\li \b Previous: TODO
7\li \b Next: TODO
8
9This tutorial explains how to solve linear systems, compute various decompositions such as LU,
10QR, SVD, eigendecompositions... for more advanced topics, don't miss our special page on
11\ref TopicLinearAlgebraDecompositions "this topic".
12
13\section TutorialLinAlgBasicSolve How do I solve a system of linear equations?
14
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
29In this example, the colPivHouseholderQr() method returns an object of class ColPivHouseholderQR. This line could
30have been replaced by:
31\code
32ColPivHouseholderQR dec(A);
33Vector3f 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
110For a much more complete table comparing all decompositions supported by Eigen (notice that Eigen
111supports many other decompositions), see our special page on
112\ref TopicLinearAlgebraDecompositions "this topic".
113
114
115
116*/
117
118}