• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1namespace Eigen {
2
3/** \eigenManualPage TopicLinearAlgebraDecompositions Catalogue of dense decompositions
4
5This page presents a catalogue of the dense matrix decompositions offered by Eigen.
6For an introduction on linear solvers and decompositions, check this \link TutorialLinearAlgebra page \endlink.
7To get an overview of the true relative speed of the different decomposition, check this \link DenseDecompositionBenchmark benchmark \endlink.
8
9\section TopicLinAlgBigTable Catalogue of decompositions offered by Eigen
10
11<table class="manual-vl">
12    <tr>
13        <th class="meta"></th>
14        <th class="meta" colspan="5">Generic information, not Eigen-specific</th>
15        <th class="meta" colspan="3">Eigen-specific</th>
16    </tr>
17
18    <tr>
19        <th>Decomposition</th>
20        <th>Requirements on the matrix</th>
21        <th>Speed</th>
22        <th>Algorithm reliability and accuracy</th>
23        <th>Rank-revealing</th>
24        <th>Allows to compute (besides linear solving)</th>
25        <th>Linear solver provided by Eigen</th>
26        <th>Maturity of Eigen's implementation</th>
27        <th>Optimizations</th>
28    </tr>
29
30    <tr>
31        <td>PartialPivLU</td>
32        <td>Invertible</td>
33        <td>Fast</td>
34        <td>Depends on condition number</td>
35        <td>-</td>
36        <td>-</td>
37        <td>Yes</td>
38        <td>Excellent</td>
39        <td>Blocking, Implicit MT</td>
40    </tr>
41
42    <tr class="alt">
43        <td>FullPivLU</td>
44        <td>-</td>
45        <td>Slow</td>
46        <td>Proven</td>
47        <td>Yes</td>
48        <td>-</td>
49        <td>Yes</td>
50        <td>Excellent</td>
51        <td>-</td>
52    </tr>
53
54    <tr>
55        <td>HouseholderQR</td>
56        <td>-</td>
57        <td>Fast</td>
58        <td>Depends on condition number</td>
59        <td>-</td>
60        <td>Orthogonalization</td>
61        <td>Yes</td>
62        <td>Excellent</td>
63        <td>Blocking</td>
64    </tr>
65
66    <tr class="alt">
67        <td>ColPivHouseholderQR</td>
68        <td>-</td>
69        <td>Fast</td>
70        <td>Good</td>
71        <td>Yes</td>
72        <td>Orthogonalization</td>
73        <td>Yes</td>
74        <td>Excellent</td>
75        <td><em>Soon: blocking</em></td>
76    </tr>
77
78    <tr>
79        <td>FullPivHouseholderQR</td>
80        <td>-</td>
81        <td>Slow</td>
82        <td>Proven</td>
83        <td>Yes</td>
84        <td>Orthogonalization</td>
85        <td>Yes</td>
86        <td>Average</td>
87        <td>-</td>
88    </tr>
89
90    <tr class="alt">
91        <td>LLT</td>
92        <td>Positive definite</td>
93        <td>Very fast</td>
94        <td>Depends on condition number</td>
95        <td>-</td>
96        <td>-</td>
97        <td>Yes</td>
98        <td>Excellent</td>
99        <td>Blocking</td>
100    </tr>
101
102    <tr>
103        <td>LDLT</td>
104        <td>Positive or negative semidefinite<sup><a href="#note1">1</a></sup></td>
105        <td>Very fast</td>
106        <td>Good</td>
107        <td>-</td>
108        <td>-</td>
109        <td>Yes</td>
110        <td>Excellent</td>
111        <td><em>Soon: blocking</em></td>
112    </tr>
113
114    <tr><th class="inter" colspan="9">\n Singular values and eigenvalues decompositions</th></tr>
115
116    <tr>
117        <td>JacobiSVD (two-sided)</td>
118        <td>-</td>
119        <td>Slow (but fast for small matrices)</td>
120        <td>Proven<sup><a href="#note3">3</a></sup></td>
121        <td>Yes</td>
122        <td>Singular values/vectors, least squares</td>
123        <td>Yes (and does least squares)</td>
124        <td>Excellent</td>
125        <td>R-SVD</td>
126    </tr>
127
128    <tr class="alt">
129        <td>SelfAdjointEigenSolver</td>
130        <td>Self-adjoint</td>
131        <td>Fast-average<sup><a href="#note2">2</a></sup></td>
132        <td>Good</td>
133        <td>Yes</td>
134        <td>Eigenvalues/vectors</td>
135        <td>-</td>
136        <td>Excellent</td>
137        <td><em>Closed forms for 2x2 and 3x3</em></td>
138    </tr>
139
140    <tr>
141        <td>ComplexEigenSolver</td>
142        <td>Square</td>
143        <td>Slow-very slow<sup><a href="#note2">2</a></sup></td>
144        <td>Depends on condition number</td>
145        <td>Yes</td>
146        <td>Eigenvalues/vectors</td>
147        <td>-</td>
148        <td>Average</td>
149        <td>-</td>
150    </tr>
151
152    <tr class="alt">
153        <td>EigenSolver</td>
154        <td>Square and real</td>
155        <td>Average-slow<sup><a href="#note2">2</a></sup></td>
156        <td>Depends on condition number</td>
157        <td>Yes</td>
158        <td>Eigenvalues/vectors</td>
159        <td>-</td>
160        <td>Average</td>
161        <td>-</td>
162    </tr>
163
164    <tr>
165        <td>GeneralizedSelfAdjointEigenSolver</td>
166        <td>Square</td>
167        <td>Fast-average<sup><a href="#note2">2</a></sup></td>
168        <td>Depends on condition number</td>
169        <td>-</td>
170        <td>Generalized eigenvalues/vectors</td>
171        <td>-</td>
172        <td>Good</td>
173        <td>-</td>
174    </tr>
175
176    <tr><th class="inter" colspan="9">\n Helper decompositions</th></tr>
177
178    <tr>
179        <td>RealSchur</td>
180        <td>Square and real</td>
181        <td>Average-slow<sup><a href="#note2">2</a></sup></td>
182        <td>Depends on condition number</td>
183        <td>Yes</td>
184        <td>-</td>
185        <td>-</td>
186        <td>Average</td>
187        <td>-</td>
188    </tr>
189
190    <tr class="alt">
191        <td>ComplexSchur</td>
192        <td>Square</td>
193        <td>Slow-very slow<sup><a href="#note2">2</a></sup></td>
194        <td>Depends on condition number</td>
195        <td>Yes</td>
196        <td>-</td>
197        <td>-</td>
198        <td>Average</td>
199        <td>-</td>
200    </tr>
201
202    <tr class="alt">
203        <td>Tridiagonalization</td>
204        <td>Self-adjoint</td>
205        <td>Fast</td>
206        <td>Good</td>
207        <td>-</td>
208        <td>-</td>
209        <td>-</td>
210        <td>Good</td>
211        <td><em>Soon: blocking</em></td>
212    </tr>
213
214    <tr>
215        <td>HessenbergDecomposition</td>
216        <td>Square</td>
217        <td>Average</td>
218        <td>Good</td>
219        <td>-</td>
220        <td>-</td>
221        <td>-</td>
222        <td>Good</td>
223        <td><em>Soon: blocking</em></td>
224    </tr>
225
226</table>
227
228\b Notes:
229<ul>
230<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>
231<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>
232<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.
233</ul>
234
235\section TopicLinAlgTerminology Terminology
236
237<dl>
238  <dt><b>Selfadjoint</b></dt>
239    <dd>For a real matrix, selfadjoint is a synonym for symmetric. For a complex matrix, selfadjoint is a synonym for \em hermitian.
240        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>
241  <dt><b>Positive/negative definite</b></dt>
242    <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$.
243        In the same vein, it is negative definite if \f$ v^* A v < 0 \f$ for any non zero vector \f$ v \f$ </dd>
244  <dt><b>Positive/negative semidefinite</b></dt>
245    <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$.
246        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>
247
248  <dt><b>Blocking</b></dt>
249    <dd>Means the algorithm can work per block, whence guaranteeing a good scaling of the performance for large matrices.</dd>
250  <dt><b>Implicit Multi Threading (MT)</b></dt>
251    <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>
252  <dt><b>Explicit Multi Threading (MT)</b></dt>
253    <dd>Means the algorithm is explicitly parallelized to take advantage of multicore processors via OpenMP.</dd>
254  <dt><b>Meta-unroller</b></dt>
255    <dd>Means the algorithm is automatically and explicitly unrolled for very small fixed size matrices.</dd>
256  <dt><b></b></dt>
257    <dd></dd>
258</dl>
259
260
261*/
262
263}
264