cart-elc

Source code for CART-ELC
git clone git://git.laack.co/cart-elc.git
Log | Files | Refs | README | LICENSE

TutorialLinearAlgebra.dox (11860B)


      1 namespace Eigen {
      2 
      3 /** \eigenManualPage TutorialLinearAlgebra Linear algebra and decompositions
      4 
      5 This page explains how to solve linear systems, compute various decompositions such as LU,
      6 QR, %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]
     15 Where \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 the properties of your matrix \a A,
     18 and depending on whether you favor speed or accuracy. However, let's start with an example that works in all cases,
     19 and 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 
     28 In this example, the colPivHouseholderQr() method returns an object of class ColPivHouseholderQR. Since here the
     29 matrix is of type Matrix3f, this line could have been replaced by:
     30 \code
     31 ColPivHouseholderQR<Matrix3f> dec(A);
     32 Vector3f x = dec.solve(b);
     33 \endcode
     34 
     35 Here, ColPivHouseholderQR is a QR decomposition with column pivoting. It's a good compromise for this tutorial, as it
     36 works for all matrices while being quite fast. Here is a table of some other decompositions that you can choose from,
     37 depending on your matrix, the problem you are trying to solve, 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>CompleteOrthogonalDecomposition</td>
     90         <td>completeOrthogonalDecomposition()</td>
     91         <td>None</td>
     92         <td>+</td>
     93         <td>-</td>
     94         <td>+++</td>
     95     </tr>
     96     <tr class="alt">
     97         <td>LLT</td>
     98         <td>llt()</td>
     99         <td>Positive definite</td>
    100         <td>+++</td>
    101         <td>+++</td>
    102         <td>+</td>
    103     </tr>
    104     <tr>
    105         <td>LDLT</td>
    106         <td>ldlt()</td>
    107         <td>Positive or negative<br/> semidefinite</td>
    108         <td>+++</td>
    109         <td>+</td>
    110         <td>++</td>
    111     </tr>
    112     <tr class="alt">
    113         <td>BDCSVD</td>
    114         <td>bdcSvd()</td>
    115         <td>None</td>
    116         <td>-</td>
    117         <td>-</td>
    118         <td>+++</td>
    119     </tr>
    120     <tr class="alt">
    121         <td>JacobiSVD</td>
    122         <td>jacobiSvd()</td>
    123         <td>None</td>
    124         <td>-</td>
    125         <td>- - -</td>
    126         <td>+++</td>
    127     </tr>
    128 </table>
    129 To get an overview of the true relative speed of the different decompositions, check this \link DenseDecompositionBenchmark benchmark \endlink.
    130 
    131 All of these decompositions offer a solve() method that works as in the above example. 
    132 
    133 If you know more about the properties of your matrix, you can use the above table to select the best method.
    134 For example, a good choice for solving linear systems with a non-symmetric matrix of full rank is PartialPivLU.
    135 If you know that your matrix is also symmetric and positive definite, the above table says that
    136 a very good choice is the LLT or LDLT decomposition. Here's an example, also demonstrating that using a general
    137 matrix (not a vector) as right hand side is possible:
    138 
    139 <table class="example">
    140 <tr><th>Example:</th><th>Output:</th></tr>
    141 <tr>
    142   <td>\include TutorialLinAlgExSolveLDLT.cpp </td>
    143   <td>\verbinclude TutorialLinAlgExSolveLDLT.out </td>
    144 </tr>
    145 </table>
    146 
    147 For a \ref TopicLinearAlgebraDecompositions "much more complete table" comparing all decompositions supported by Eigen (notice that Eigen
    148 supports many other decompositions), see our special page on
    149 \ref TopicLinearAlgebraDecompositions "this topic".
    150 
    151 
    152 \section TutorialLinAlgLeastsquares Least squares solving
    153 
    154 The most general and accurate method to solve under- or over-determined linear systems
    155 in the least squares sense, is the SVD decomposition. Eigen provides two implementations.
    156 The recommended one is the BDCSVD class, which scales well for large problems
    157 and automatically falls back to the JacobiSVD class for smaller problems.
    158 For both classes, their solve() method solved the linear system in the least-squares
    159 sense. 
    160 
    161 Here is an example:
    162 <table class="example">
    163 <tr><th>Example:</th><th>Output:</th></tr>
    164 <tr>
    165   <td>\include TutorialLinAlgSVDSolve.cpp </td>
    166   <td>\verbinclude TutorialLinAlgSVDSolve.out </td>
    167 </tr>
    168 </table>
    169 
    170 An alternative to the SVD, which is usually faster and about as accurate, is CompleteOrthogonalDecomposition. 
    171 
    172 Again, if you know more about the problem, the table above contains methods that are potentially faster.
    173 If your matrix is full rank, HouseHolderQR is the method of choice. If your matrix is full rank and well conditioned,
    174 using the Cholesky decomposition (LLT) on the matrix of the normal equations can be faster still.
    175 Our page on \link LeastSquares least squares solving \endlink has more details.
    176 
    177 
    178 \section TutorialLinAlgSolutionExists Checking if a matrix is singular
    179 
    180 Only you know what error margin you want to allow for a solution to be considered valid.
    181 So Eigen lets you do this computation for yourself, if you want to, as in this example:
    182 
    183 <table class="example">
    184 <tr><th>Example:</th><th>Output:</th></tr>
    185 <tr>
    186   <td>\include TutorialLinAlgExComputeSolveError.cpp </td>
    187   <td>\verbinclude TutorialLinAlgExComputeSolveError.out </td>
    188 </tr>
    189 </table>
    190 
    191 \section TutorialLinAlgEigensolving Computing eigenvalues and eigenvectors
    192 
    193 You need an eigendecomposition here, see available such decompositions on \ref TopicLinearAlgebraDecompositions "this page".
    194 Make sure to check if your matrix is self-adjoint, as is often the case in these problems. Here's an example using
    195 SelfAdjointEigenSolver, it could easily be adapted to general matrices using EigenSolver or ComplexEigenSolver.
    196 
    197 The computation of eigenvalues and eigenvectors does not necessarily converge, but such failure to converge is
    198 very rare. The call to info() is to check for this possibility.
    199 
    200 <table class="example">
    201 <tr><th>Example:</th><th>Output:</th></tr>
    202 <tr>
    203   <td>\include TutorialLinAlgSelfAdjointEigenSolver.cpp </td>
    204   <td>\verbinclude TutorialLinAlgSelfAdjointEigenSolver.out </td>
    205 </tr>
    206 </table>
    207 
    208 \section TutorialLinAlgInverse Computing inverse and determinant
    209 
    210 First of all, make sure that you really want this. While inverse and determinant are fundamental mathematical concepts,
    211 in \em numerical linear algebra they are not as useful as in pure mathematics. Inverse computations are often
    212 advantageously replaced by solve() operations, and the determinant is often \em not a good way of checking if a matrix
    213 is invertible.
    214 
    215 However, for \em very \em small matrices, the above may not be true, and inverse and determinant can be very useful.
    216 
    217 While certain decompositions, such as PartialPivLU and FullPivLU, offer inverse() and determinant() methods, you can also
    218 call inverse() and determinant() directly on a matrix. If your matrix is of a very small fixed size (at most 4x4) this
    219 allows Eigen to avoid performing a LU decomposition, and instead use formulas that are more efficient on such small matrices.
    220 
    221 Here is an example:
    222 <table class="example">
    223 <tr><th>Example:</th><th>Output:</th></tr>
    224 <tr>
    225   <td>\include TutorialLinAlgInverseDeterminant.cpp </td>
    226   <td>\verbinclude TutorialLinAlgInverseDeterminant.out </td>
    227 </tr>
    228 </table>
    229 
    230 \section TutorialLinAlgSeparateComputation Separating the computation from the construction
    231 
    232 In the above examples, the decomposition was computed at the same time that the decomposition object was constructed.
    233 There are however situations where you might want to separate these two things, for example if you don't know,
    234 at the time of the construction, the matrix that you will want to decompose; or if you want to reuse an existing
    235 decomposition object.
    236 
    237 What makes this possible is that:
    238 \li all decompositions have a default constructor,
    239 \li all decompositions have a compute(matrix) method that does the computation, and that may be called again
    240     on an already-computed decomposition, reinitializing it.
    241 
    242 For example:
    243 
    244 <table class="example">
    245 <tr><th>Example:</th><th>Output:</th></tr>
    246 <tr>
    247   <td>\include TutorialLinAlgComputeTwice.cpp </td>
    248   <td>\verbinclude TutorialLinAlgComputeTwice.out </td>
    249 </tr>
    250 </table>
    251 
    252 Finally, you can tell the decomposition constructor to preallocate storage for decomposing matrices of a given size,
    253 so that when you subsequently decompose such matrices, no dynamic memory allocation is performed (of course, if you
    254 are using fixed-size matrices, no dynamic memory allocation happens at all). This is done by just
    255 passing the size to the decomposition constructor, as in this example:
    256 \code
    257 HouseholderQR<MatrixXf> qr(50,50);
    258 MatrixXf A = MatrixXf::Random(50,50);
    259 qr.compute(A); // no dynamic memory allocation
    260 \endcode
    261 
    262 \section TutorialLinAlgRankRevealing Rank-revealing decompositions
    263 
    264 Certain decompositions are rank-revealing, i.e. are able to compute the rank of a matrix. These are typically
    265 also the decompositions that behave best in the face of a non-full-rank matrix (which in the square case means a
    266 singular matrix). On \ref TopicLinearAlgebraDecompositions "this table" you can see for all our decompositions
    267 whether they are rank-revealing or not.
    268 
    269 Rank-revealing decompositions offer at least a rank() method. They can also offer convenience methods such as isInvertible(),
    270 and some are also providing methods to compute the kernel (null-space) and image (column-space) of the matrix, as is the
    271 case with FullPivLU:
    272 
    273 <table class="example">
    274 <tr><th>Example:</th><th>Output:</th></tr>
    275 <tr>
    276   <td>\include TutorialLinAlgRankRevealing.cpp </td>
    277   <td>\verbinclude TutorialLinAlgRankRevealing.out </td>
    278 </tr>
    279 </table>
    280 
    281 Of course, any rank computation depends on the choice of an arbitrary threshold, since practically no
    282 floating-point matrix is \em exactly rank-deficient. Eigen picks a sensible default threshold, which depends
    283 on the decomposition but is typically the diagonal size times machine epsilon. While this is the best default we
    284 could pick, only you know what is the right threshold for your application. You can set this by calling setThreshold()
    285 on your decomposition object before calling rank() or any other method that needs to use such a threshold.
    286 The decomposition itself, i.e. the compute() method, is independent of the threshold. You don't need to recompute the
    287 decomposition after you've changed the threshold.
    288 
    289 <table class="example">
    290 <tr><th>Example:</th><th>Output:</th></tr>
    291 <tr>
    292   <td>\include TutorialLinAlgSetThreshold.cpp </td>
    293   <td>\verbinclude TutorialLinAlgSetThreshold.out </td>
    294 </tr>
    295 </table>
    296 
    297 */
    298 
    299 }