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.

260 lines
7.9 KiB

  1. namespace Eigen {
  2. /** \page TopicLinearAlgebraDecompositions Linear algebra and decompositions
  3. \section TopicLinAlgBigTable Catalogue of decompositions offered by Eigen
  4. <table class="manual-vl">
  5. <tr>
  6. <th class="meta"></th>
  7. <th class="meta" colspan="5">Generic information, not Eigen-specific</th>
  8. <th class="meta" colspan="3">Eigen-specific</th>
  9. </tr>
  10. <tr>
  11. <th>Decomposition</th>
  12. <th>Requirements on the matrix</th>
  13. <th>Speed</th>
  14. <th>Algorithm reliability and accuracy</th>
  15. <th>Rank-revealing</th>
  16. <th>Allows to compute (besides linear solving)</th>
  17. <th>Linear solver provided by Eigen</th>
  18. <th>Maturity of Eigen's implementation</th>
  19. <th>Optimizations</th>
  20. </tr>
  21. <tr>
  22. <td>PartialPivLU</td>
  23. <td>Invertible</td>
  24. <td>Fast</td>
  25. <td>Depends on condition number</td>
  26. <td>-</td>
  27. <td>-</td>
  28. <td>Yes</td>
  29. <td>Excellent</td>
  30. <td>Blocking, Implicit MT</td>
  31. </tr>
  32. <tr class="alt">
  33. <td>FullPivLU</td>
  34. <td>-</td>
  35. <td>Slow</td>
  36. <td>Proven</td>
  37. <td>Yes</td>
  38. <td>-</td>
  39. <td>Yes</td>
  40. <td>Excellent</td>
  41. <td>-</td>
  42. </tr>
  43. <tr>
  44. <td>HouseholderQR</td>
  45. <td>-</td>
  46. <td>Fast</td>
  47. <td>Depends on condition number</td>
  48. <td>-</td>
  49. <td>Orthogonalization</td>
  50. <td>Yes</td>
  51. <td>Excellent</td>
  52. <td>Blocking</td>
  53. </tr>
  54. <tr class="alt">
  55. <td>ColPivHouseholderQR</td>
  56. <td>-</td>
  57. <td>Fast</td>
  58. <td>Good</td>
  59. <td>Yes</td>
  60. <td>Orthogonalization</td>
  61. <td>Yes</td>
  62. <td>Excellent</td>
  63. <td><em>Soon: blocking</em></td>
  64. </tr>
  65. <tr>
  66. <td>FullPivHouseholderQR</td>
  67. <td>-</td>
  68. <td>Slow</td>
  69. <td>Proven</td>
  70. <td>Yes</td>
  71. <td>Orthogonalization</td>
  72. <td>Yes</td>
  73. <td>Average</td>
  74. <td>-</td>
  75. </tr>
  76. <tr class="alt">
  77. <td>LLT</td>
  78. <td>Positive definite</td>
  79. <td>Very fast</td>
  80. <td>Depends on condition number</td>
  81. <td>-</td>
  82. <td>-</td>
  83. <td>Yes</td>
  84. <td>Excellent</td>
  85. <td>Blocking</td>
  86. </tr>
  87. <tr>
  88. <td>LDLT</td>
  89. <td>Positive or negative semidefinite<sup><a href="#note1">1</a></sup></td>
  90. <td>Very fast</td>
  91. <td>Good</td>
  92. <td>-</td>
  93. <td>-</td>
  94. <td>Yes</td>
  95. <td>Excellent</td>
  96. <td><em>Soon: blocking</em></td>
  97. </tr>
  98. <tr><th class="inter" colspan="9">\n Singular values and eigenvalues decompositions</th></tr>
  99. <tr>
  100. <td>JacobiSVD (two-sided)</td>
  101. <td>-</td>
  102. <td>Slow (but fast for small matrices)</td>
  103. <td>Excellent-Proven<sup><a href="#note3">3</a></sup></td>
  104. <td>Yes</td>
  105. <td>Singular values/vectors, least squares</td>
  106. <td>Yes (and does least squares)</td>
  107. <td>Excellent</td>
  108. <td>R-SVD</td>
  109. </tr>
  110. <tr class="alt">
  111. <td>SelfAdjointEigenSolver</td>
  112. <td>Self-adjoint</td>
  113. <td>Fast-average<sup><a href="#note2">2</a></sup></td>
  114. <td>Good</td>
  115. <td>Yes</td>
  116. <td>Eigenvalues/vectors</td>
  117. <td>-</td>
  118. <td>Good</td>
  119. <td><em>Closed forms for 2x2 and 3x3</em></td>
  120. </tr>
  121. <tr>
  122. <td>ComplexEigenSolver</td>
  123. <td>Square</td>
  124. <td>Slow-very slow<sup><a href="#note2">2</a></sup></td>
  125. <td>Depends on condition number</td>
  126. <td>Yes</td>
  127. <td>Eigenvalues/vectors</td>
  128. <td>-</td>
  129. <td>Average</td>
  130. <td>-</td>
  131. </tr>
  132. <tr class="alt">
  133. <td>EigenSolver</td>
  134. <td>Square and real</td>
  135. <td>Average-slow<sup><a href="#note2">2</a></sup></td>
  136. <td>Depends on condition number</td>
  137. <td>Yes</td>
  138. <td>Eigenvalues/vectors</td>
  139. <td>-</td>
  140. <td>Average</td>
  141. <td>-</td>
  142. </tr>
  143. <tr>
  144. <td>GeneralizedSelfAdjointEigenSolver</td>
  145. <td>Square</td>
  146. <td>Fast-average<sup><a href="#note2">2</a></sup></td>
  147. <td>Depends on condition number</td>
  148. <td>-</td>
  149. <td>Generalized eigenvalues/vectors</td>
  150. <td>-</td>
  151. <td>Good</td>
  152. <td>-</td>
  153. </tr>
  154. <tr><th class="inter" colspan="9">\n Helper decompositions</th></tr>
  155. <tr>
  156. <td>RealSchur</td>
  157. <td>Square and real</td>
  158. <td>Average-slow<sup><a href="#note2">2</a></sup></td>
  159. <td>Depends on condition number</td>
  160. <td>Yes</td>
  161. <td>-</td>
  162. <td>-</td>
  163. <td>Average</td>
  164. <td>-</td>
  165. </tr>
  166. <tr class="alt">
  167. <td>ComplexSchur</td>
  168. <td>Square</td>
  169. <td>Slow-very slow<sup><a href="#note2">2</a></sup></td>
  170. <td>Depends on condition number</td>
  171. <td>Yes</td>
  172. <td>-</td>
  173. <td>-</td>
  174. <td>Average</td>
  175. <td>-</td>
  176. </tr>
  177. <tr class="alt">
  178. <td>Tridiagonalization</td>
  179. <td>Self-adjoint</td>
  180. <td>Fast</td>
  181. <td>Good</td>
  182. <td>-</td>
  183. <td>-</td>
  184. <td>-</td>
  185. <td>Good</td>
  186. <td><em>Soon: blocking</em></td>
  187. </tr>
  188. <tr>
  189. <td>HessenbergDecomposition</td>
  190. <td>Square</td>
  191. <td>Average</td>
  192. <td>Good</td>
  193. <td>-</td>
  194. <td>-</td>
  195. <td>-</td>
  196. <td>Good</td>
  197. <td><em>Soon: blocking</em></td>
  198. </tr>
  199. </table>
  200. \b Notes:
  201. <ul>
  202. <li><a name="note1">\b 1: </a>There exist two variants of the LDLT algorithm. Eigen's one produces a pure diagonal D matrix, and therefore it cannot handle indefinite matrices, unlike Lapack's one which produces a block diagonal D matrix.</li>
  203. <li><a name="note2">\b 2: </a>Eigenvalues, SVD and Schur decompositions rely on iterative algorithms. Their convergence speed depends on how well the eigenvalues are separated.</li>
  204. <li><a name="note3">\b 3: </a>Our JacobiSVD is two-sided, making for proven and optimal precision for square matrices. For non-square matrices, we have to use a QR preconditioner first. The default choice, ColPivHouseholderQR, is already very reliable, but if you want it to be proven, use FullPivHouseholderQR instead.
  205. </ul>
  206. \section TopicLinAlgTerminology Terminology
  207. <dl>
  208. <dt><b>Selfadjoint</b></dt>
  209. <dd>For a real matrix, selfadjoint is a synonym for symmetric. For a complex matrix, selfadjoint is a synonym for \em hermitian.
  210. More generally, a matrix \f$ A \f$ is selfadjoint if and only if it is equal to its adjoint \f$ A^* \f$. The adjoint is also called the \em conjugate \em transpose. </dd>
  211. <dt><b>Positive/negative definite</b></dt>
  212. <dd>A selfadjoint matrix \f$ A \f$ is positive definite if \f$ v^* A v > 0 \f$ for any non zero vector \f$ v \f$.
  213. In the same vein, it is negative definite if \f$ v^* A v < 0 \f$ for any non zero vector \f$ v \f$ </dd>
  214. <dt><b>Positive/negative semidefinite</b></dt>
  215. <dd>A selfadjoint matrix \f$ A \f$ is positive semi-definite if \f$ v^* A v \ge 0 \f$ for any non zero vector \f$ v \f$.
  216. In the same vein, it is negative semi-definite if \f$ v^* A v \le 0 \f$ for any non zero vector \f$ v \f$ </dd>
  217. <dt><b>Blocking</b></dt>
  218. <dd>Means the algorithm can work per block, whence guaranteeing a good scaling of the performance for large matrices.</dd>
  219. <dt><b>Implicit Multi Threading (MT)</b></dt>
  220. <dd>Means the algorithm can take advantage of multicore processors via OpenMP. "Implicit" means the algortihm itself is not parallelized, but that it relies on parallelized matrix-matrix product rountines.</dd>
  221. <dt><b>Explicit Multi Threading (MT)</b></dt>
  222. <dd>Means the algorithm is explicitely parallelized to take advantage of multicore processors via OpenMP.</dd>
  223. <dt><b>Meta-unroller</b></dt>
  224. <dd>Means the algorithm is automatically and explicitly unrolled for very small fixed size matrices.</dd>
  225. <dt><b></b></dt>
  226. <dd></dd>
  227. </dl>
  228. */
  229. }