1namespace Eigen {
2
3/** \eigenManualPage TutorialLinearAlgebra Linear algebra and decompositions
4
5This page explains how to solve linear systems, compute various decompositions such as LU,
6QR, %SVD, eigendecompositions... After reading this page, don't miss our
7\link TopicLinearAlgebraDecompositions catalogue \endlink of dense matrix decompositions.
8
9\eigenAutoToc
10
11\section TutorialLinAlgBasicSolve Basic linear solving
12
13\b The \b problem: You have a system of equations, that you have written as a single matrix equation
14    \f[ Ax \: = \: b \f]
15Where \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.
16
17\b The \b solution: You can choose between various decompositions, depending on what your matrix \a A looks like,
18and depending on whether you favor speed or accuracy. However, let's start with an example that works in all cases,
19and is a good compromise:
20<table class="example">
21<tr><th>Example:</th><th>Output:</th></tr>
22<tr>
23  <td>\include TutorialLinAlgExSolveColPivHouseholderQR.cpp </td>
24  <td>\verbinclude TutorialLinAlgExSolveColPivHouseholderQR.out </td>
25</tr>
26</table>
27
28In this example, the colPivHouseholderQr() method returns an object of class ColPivHouseholderQR. Since here the
29matrix is of type Matrix3f, this line could have been replaced by:
30\code
31ColPivHouseholderQR<Matrix3f> dec(A);
32Vector3f x = dec.solve(b);
33\endcode
34
35Here, ColPivHouseholderQR is a QR decomposition with column pivoting. It's a good compromise for this tutorial, as it
36works for all matrices while being quite fast. Here is a table of some other decompositions that you can choose from,
37depending on your matrix and the trade-off you want to make:
38
39<table class="manual">
40    <tr>
41        <th>Decomposition</th>
42        <th>Method</th>
43        <th>Requirements<br/>on the matrix</th>
44        <th>Speed<br/> (small-to-medium)</th>
45        <th>Speed<br/> (large)</th>
46        <th>Accuracy</th>
47    </tr>
48    <tr>
49        <td>PartialPivLU</td>
50        <td>partialPivLu()</td>
51        <td>Invertible</td>
52        <td>++</td>
53        <td>++</td>
54        <td>+</td>
55    </tr>
56    <tr class="alt">
57        <td>FullPivLU</td>
58        <td>fullPivLu()</td>
59        <td>None</td>
60        <td>-</td>
61        <td>- -</td>
62        <td>+++</td>
63    </tr>
64    <tr>
65        <td>HouseholderQR</td>
66        <td>householderQr()</td>
67        <td>None</td>
68        <td>++</td>
69        <td>++</td>
70        <td>+</td>
71    </tr>
72    <tr class="alt">
73        <td>ColPivHouseholderQR</td>
74        <td>colPivHouseholderQr()</td>
75        <td>None</td>
76        <td>++</td>
77        <td>-</td>
78        <td>+++</td>
79    </tr>
80    <tr>
81        <td>FullPivHouseholderQR</td>
82        <td>fullPivHouseholderQr()</td>
83        <td>None</td>
84        <td>-</td>
85        <td>- -</td>
86        <td>+++</td>
87    </tr>
88    <tr class="alt">
89        <td>LLT</td>
90        <td>llt()</td>
91        <td>Positive definite</td>
92        <td>+++</td>
93        <td>+++</td>
94        <td>+</td>
95    </tr>
96    <tr>
97        <td>LDLT</td>
98        <td>ldlt()</td>
99        <td>Positive or negative<br/> semidefinite</td>
100        <td>+++</td>
101        <td>+</td>
102        <td>++</td>
103    </tr>
104    <tr class="alt">
105        <td>JacobiSVD</td>
106        <td>jacobiSvd()</td>
107        <td>None</td>
108        <td>- -</td>
109        <td>- - -</td>
110        <td>+++</td>
111    </tr>
112</table>
113
114All of these decompositions offer a solve() method that works as in the above example.
115
116For example, if your matrix is positive definite, the above table says that a very good
117choice is then the LLT or LDLT decomposition. Here's an example, also demonstrating that using a general
118matrix (not a vector) as right hand side is possible.
119
120<table class="example">
121<tr><th>Example:</th><th>Output:</th></tr>
122<tr>
123  <td>\include TutorialLinAlgExSolveLDLT.cpp </td>
124  <td>\verbinclude TutorialLinAlgExSolveLDLT.out </td>
125</tr>
126</table>
127
128For a \ref TopicLinearAlgebraDecompositions "much more complete table" comparing all decompositions supported by Eigen (notice that Eigen
129supports many other decompositions), see our special page on
130\ref TopicLinearAlgebraDecompositions "this topic".
131
132\section TutorialLinAlgSolutionExists Checking if a solution really exists
133
134Only you know what error margin you want to allow for a solution to be considered valid.
135So Eigen lets you do this computation for yourself, if you want to, as in this example:
136
137<table class="example">
138<tr><th>Example:</th><th>Output:</th></tr>
139<tr>
140  <td>\include TutorialLinAlgExComputeSolveError.cpp </td>
141  <td>\verbinclude TutorialLinAlgExComputeSolveError.out </td>
142</tr>
143</table>
144
145\section TutorialLinAlgEigensolving Computing eigenvalues and eigenvectors
146
147You need an eigendecomposition here, see available such decompositions on \ref TopicLinearAlgebraDecompositions "this page".
148Make sure to check if your matrix is self-adjoint, as is often the case in these problems. Here's an example using
149SelfAdjointEigenSolver, it could easily be adapted to general matrices using EigenSolver or ComplexEigenSolver.
150
151The computation of eigenvalues and eigenvectors does not necessarily converge, but such failure to converge is
152very rare. The call to info() is to check for this possibility.
153
154<table class="example">
155<tr><th>Example:</th><th>Output:</th></tr>
156<tr>
157  <td>\include TutorialLinAlgSelfAdjointEigenSolver.cpp </td>
158  <td>\verbinclude TutorialLinAlgSelfAdjointEigenSolver.out </td>
159</tr>
160</table>
161
162\section TutorialLinAlgInverse Computing inverse and determinant
163
164First of all, make sure that you really want this. While inverse and determinant are fundamental mathematical concepts,
165in \em numerical linear algebra they are not as popular as in pure mathematics. Inverse computations are often
166advantageously replaced by solve() operations, and the determinant is often \em not a good way of checking if a matrix
167is invertible.
168
169However, for \em very \em small matrices, the above is not true, and inverse and determinant can be very useful.
170
171While certain decompositions, such as PartialPivLU and FullPivLU, offer inverse() and determinant() methods, you can also
172call inverse() and determinant() directly on a matrix. If your matrix is of a very small fixed size (at most 4x4) this
173allows Eigen to avoid performing a LU decomposition, and instead use formulas that are more efficient on such small matrices.
174
175Here is an example:
176<table class="example">
177<tr><th>Example:</th><th>Output:</th></tr>
178<tr>
179  <td>\include TutorialLinAlgInverseDeterminant.cpp </td>
180  <td>\verbinclude TutorialLinAlgInverseDeterminant.out </td>
181</tr>
182</table>
183
184\section TutorialLinAlgLeastsquares Least squares solving
185
186The most accurate method to do least squares solving is with a SVD decomposition. Eigen provides one
187as the JacobiSVD class, and its solve() is doing least-squares solving.
188
189Here is an example:
190<table class="example">
191<tr><th>Example:</th><th>Output:</th></tr>
192<tr>
193  <td>\include TutorialLinAlgSVDSolve.cpp </td>
194  <td>\verbinclude TutorialLinAlgSVDSolve.out </td>
195</tr>
196</table>
197
198Another methods, potentially faster but less reliable, are to use a Cholesky decomposition of the
199normal matrix or a QR decomposition. Our page on \link LeastSquares least squares solving \endlink
200has more details.
201
202
203\section TutorialLinAlgSeparateComputation Separating the computation from the construction
204
205In the above examples, the decomposition was computed at the same time that the decomposition object was constructed.
206There are however situations where you might want to separate these two things, for example if you don't know,
207at the time of the construction, the matrix that you will want to decompose; or if you want to reuse an existing
208decomposition object.
209
210What makes this possible is that:
211\li all decompositions have a default constructor,
212\li all decompositions have a compute(matrix) method that does the computation, and that may be called again
213    on an already-computed decomposition, reinitializing it.
214
215For example:
216
217<table class="example">
218<tr><th>Example:</th><th>Output:</th></tr>
219<tr>
220  <td>\include TutorialLinAlgComputeTwice.cpp </td>
221  <td>\verbinclude TutorialLinAlgComputeTwice.out </td>
222</tr>
223</table>
224
225Finally, you can tell the decomposition constructor to preallocate storage for decomposing matrices of a given size,
226so that when you subsequently decompose such matrices, no dynamic memory allocation is performed (of course, if you
227are using fixed-size matrices, no dynamic memory allocation happens at all). This is done by just
228passing the size to the decomposition constructor, as in this example:
229\code
230HouseholderQR<MatrixXf> qr(50,50);
231MatrixXf A = MatrixXf::Random(50,50);
232qr.compute(A); // no dynamic memory allocation
233\endcode
234
235\section TutorialLinAlgRankRevealing Rank-revealing decompositions
236
237Certain decompositions are rank-revealing, i.e. are able to compute the rank of a matrix. These are typically
238also the decompositions that behave best in the face of a non-full-rank matrix (which in the square case means a
239singular matrix). On \ref TopicLinearAlgebraDecompositions "this table" you can see for all our decompositions
240whether they are rank-revealing or not.
241
242Rank-revealing decompositions offer at least a rank() method. They can also offer convenience methods such as isInvertible(),
243and some are also providing methods to compute the kernel (null-space) and image (column-space) of the matrix, as is the
244case with FullPivLU:
245
246<table class="example">
247<tr><th>Example:</th><th>Output:</th></tr>
248<tr>
249  <td>\include TutorialLinAlgRankRevealing.cpp </td>
250  <td>\verbinclude TutorialLinAlgRankRevealing.out </td>
251</tr>
252</table>
253
254Of course, any rank computation depends on the choice of an arbitrary threshold, since practically no
255floating-point matrix is \em exactly rank-deficient. Eigen picks a sensible default threshold, which depends
256on the decomposition but is typically the diagonal size times machine epsilon. While this is the best default we
257could pick, only you know what is the right threshold for your application. You can set this by calling setThreshold()
258on your decomposition object before calling rank() or any other method that needs to use such a threshold.
259The decomposition itself, i.e. the compute() method, is independent of the threshold. You don't need to recompute the
260decomposition after you've changed the threshold.
261
262<table class="example">
263<tr><th>Example:</th><th>Output:</th></tr>
264<tr>
265  <td>\include TutorialLinAlgSetThreshold.cpp </td>
266  <td>\verbinclude TutorialLinAlgSetThreshold.out </td>
267</tr>
268</table>
269
270*/
271
272}
273