You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

272 lines
10 KiB

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