blob: e24c4f00f9d846141e997f45f3ad3f8a689ca936 [file] [log] [blame]
Gael Guennebaud5c866f22010-06-26 16:59:18 +02001namespace Eigen {
2
3/** \page QuickRefPage Quick reference guide
4
5\b Table \b of \b contents
6 - \ref QuickRef_Headers
7 - \ref QuickRef_Types
8 - \ref QuickRef_Map
9 - \ref QuickRef_ArithmeticOperators
10 - \ref QuickRef_Coeffwise
Gael Guennebaud1c783e22010-06-26 18:49:50 +020011 - \ref QuickRef_Reductions
12 - \ref QuickRef_Blocks
13 - \ref QuickRef_DiagTriSymm
Gael Guennebaud5c866f22010-06-26 16:59:18 +020014\n
15
16<hr>
17
18<a href="#" class="top">top</a>
19\section QuickRef_Headers Modules and Header files
20
21<table class="tutorial_code">
22<tr><td>Module</td><td>Header file</td><td>Contents</td></tr>
23<tr><td>Core</td><td>\code#include <Eigen/Core>\endcode</td><td>Matrix and Array classes, basic linear algebra (including triangular and selfadjoint products), array manipulation</td></tr>
24<tr><td>Geometry</td><td>\code#include <Eigen/Geometry>\endcode</td><td>Transformation, Translation, Scaling, 2D and 3D rotations (Quaternion, AngleAxis)</td></tr>
25<tr><td>LU</td><td>\code#include <Eigen/LU>\endcode</td><td>Inverse, determinant, LU decompositions (FullPivLU, PartialPivLU) with solver</td></tr>
26<tr><td>Cholesky</td><td>\code#include <Eigen/Cholesky>\endcode</td><td>LLT and LDLT Cholesky factorization with solver</td></tr>
27<tr><td>SVD</td><td>\code#include <Eigen/SVD>\endcode</td><td>SVD decomposition with solver (HouseholderSVD, JacobiSVD)</td></tr>
28<tr><td>QR</td><td>\code#include <Eigen/QR>\endcode</td><td>QR decomposition with solver (HouseholderQR, ColPivHouseholerQR, FullPivHouseholderQR)</td></tr>
29<tr><td>Eigenvalues</td><td>\code#include <Eigen/Eigenvalues>\endcode</td><td>Eigenvalue, eigenvector decompositions for selfadjoint and non selfadjoint real or complex matrices.</td></tr>
30<tr><td>Sparse</td><td>\code#include <Eigen/Sparse>\endcode</td><td>Sparse matrix storage and related basic linear algebra.</td></tr>
31<tr><td></td><td>\code#include <Eigen/Dense>\endcode</td><td>Includes Core, Geometry, LU, Cholesky, SVD, QR, and Eigenvalues</td></tr>
32<tr><td></td><td>\code#include <Eigen/Eigen>\endcode</td><td>Includes Dense and Sparse</td></tr>
33</table>
34
35<a href="#" class="top">top</a>
36\section QuickRef_Types Array, matrix and vector types
37
Gael Guennebauddbefd7a2010-06-28 13:30:10 +020038
Gael Guennebaud5c866f22010-06-26 16:59:18 +020039\b Recall: Eigen provides two kinds of dense objects: mathematical matrices and vectors which are both represented by the template class Matrix, and general 1D and 2D arrays represented by the template class Array:
40\code
41typedef Matrix<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyMatrixType;
42typedef Array<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyArrayType;
43\endcode
44
45\li \c Scalar is the scalar type of the coefficients (e.g., \c float, \c double, \c bool, \c int, etc.).
46\li \c RowsAtCompileTime and \c ColsAtCompileTime are the number of rows and columns of the matrix as known at compile-time or \c Dynamic.
47\li \c Options can be \c ColMajor or \c RowMajor, default is \c ColMajor. (see class Matrix for more options)
48
49All combinations are allowed: you can have a matrix with a fixed number of rows and a dynamic number of columns, etc. The following are all valid:
50
51\code
52Matrix<double, 6, Dynamic> // Dynamic number of columns (heap allocation)
53Matrix<double, Dynamic, 2> // Dynamic number of rows (heap allocation)
54Matrix<double, Dynamic, Dynamic, RowMajor> // Fully dynamic, row major (heap allocation)
55Matrix<double, 13, 3> // Fully fixed (static allocation)
56\endcode
57
58In most cases, you can simply use one of the convenience typedefs for \ref matrixtypedefs "matrices" and \ref arraytypedefs "arrays". Some examples:
59<table class="tutorial_code">
60<tr><td>\code
61Matrix<float,Dynamic,Dynamic> <=> MatrixXf
62Matrix<double,Dynamic,1> <=> VectorXd
63Matrix<int,1,Dynamic> <=> RowVectorXi
64Matrix<float,3,3> <=> Matrix3f
65Matrix<float,4,1> <=> Vector4f
66\endcode</td><td>\code
67Array<float,Dynamic,Dynamic> <=> ArrayXXf
68Array<double,Dynamic,1> <=> ArrayXd
69Array<int,1,Dynamic> <=> RowArrayXi
70Array<float,3,3> <=> Array33f
71Array<float,4,1> <=> Array4f
72\endcode</td></tr>
73</table>
74
75Conversion between the matrix and array worlds:
76\code
77Array44f a1, a1;
78Matrix4f m1, m2;
Gael Guennebaudf3c64c72010-06-26 22:19:03 +020079m1 = a1 * a2; // coeffwise product, implicit conversion from array to matrix.
80a1 = m1 * m2; // matrix product, implicit conversion from matrix to array.
81a2 = a1 + m1.array(); // mixing array and matrix is forbidden
82m2 = a1.matrix() + m1; // and explicit conversion is required.
83ArrayWrapper<Matrix4f> m1a(m1); // m1a is an alias for m1.array(), they share the same coefficients
Gael Guennebaud5c866f22010-06-26 16:59:18 +020084MatrixWrapper<Array44f> a1m(a1);
85\endcode
86
87In the rest of this document we will use the following symbols to emphasize the features which are specifics to a given kind of object:
88\li <a name="matrixonly"><a/>\matrixworld linear algebra matrix and vector only
89\li <a name="arrayonly"><a/>\arrayworld array objects only
90
Gael Guennebaud5c866f22010-06-26 16:59:18 +020091\subsection QuickRef_Basics Basic matrix manipulation
92
93<table class="tutorial_code">
94<tr><td></td><td>1D objects</td><td>2D objects</td><td>Notes</td></tr>
95<tr><td>Constructors</td>
96<td>\code
97Vector4d v4;
98Vector2f v1(x, y);
99Array3i v2(x, y, z);
100Vector4d v3(x, y, z, w);
101
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200102VectorXf v5; // empty object
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200103ArrayXf v6(size);
104\endcode</td><td>\code
105Matrix4f m1;
106
107
108
109
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200110MatrixXf m5; // empty object
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200111MatrixXf m6(nb_rows, nb_columns);
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200112\endcode</td><td class="note">
113By default, the coefficients \n are left uninitialized</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200114<tr><td>Comma initializer</td>
115<td>\code
116Vector3f v1; v1 << x, y, z;
117ArrayXf v2(4); v2 << 1, 2, 3, 4;
118
119\endcode</td><td>\code
120Matrix3f m1; m1 << 1, 2, 3,
121 4, 5, 6,
122 7, 8, 9;
123\endcode</td><td></td></tr>
124
125<tr><td>Comma initializer (bis)</td>
126<td colspan="2">
127<table class="noborder"><tr><td>
128\include Tutorial_commainit_02.cpp
129</td>
130<td>
131output:
132\verbinclude Tutorial_commainit_02.out
133</td></tr></table>
134<td></td></tr>
135
136<tr><td>Runtime info</td>
137<td>\code
138vector.size();
139
140vector.innerStride();
141vector.data();
142\endcode</td><td>\code
143matrix.rows(); matrix.cols();
144matrix.innerSize(); matrix.outerSize();
145matrix.innerStride(); matrix.outerStride();
146matrix.data();
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200147\endcode</td><td class="note">\n Inner/Outer* are storage order dependent</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200148<tr><td>Compile-time info</td>
149<td colspan="2">\code
150ObjectType::Scalar ObjectType::RowsAtCompileTime
151ObjectType::RealScalar ObjectType::ColsAtCompileTime
152ObjectType::Index ObjectType::SizeAtCompileTime
153\endcode</td><td></td></tr>
154<tr><td>Resizing</td>
155<td>\code
156vector.resize(size);
157
158
159vector.resizeLike(other_vector);
160vector.conservativeResize(size);
161\endcode</td><td>\code
162matrix.resize(nb_rows, nb_cols);
163matrix.resize(Eigen::NoChange, nb_cols);
164matrix.resize(nb_rows, Eigen::NoChange);
165matrix.resizeLike(other_matrix);
166matrix.conservativeResize(nb_rows, nb_cols);
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200167\endcode</td><td class="note">no-op if the new sizes match,\n otherwise data are lost \n \n resizing with data preservation</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200168
169<tr><td>Coeff access with \n range checking</td>
170<td>\code
171vector(i) vector.x()
172vector[i] vector.y()
173 vector.z()
174 vector.w()
175\endcode</td><td>\code
176matrix(i,j)
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200177\endcode</td><td class="note">Range checking is disabled if \n NDEBUG or EIGEN_NO_DEBUG is defined</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200178
179<tr><td>Coeff access without \n range checking</td>
180<td>\code
181vector.coeff(i)
182vector.coeffRef(i)
183\endcode</td><td>\code
184matrix.coeff(i,j)
185matrix.coeffRef(i,j)
186\endcode</td><td></td></tr>
187
188<tr><td>Assignment/copy</td>
189<td colspan="2">\code
190object = expression;
191object_of_float = expression_of_double.cast<float>();
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200192\endcode</td><td class="note">the destination is automatically resized (if possible)</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200193
194</table>
195
196\subsection QuickRef_PredefMat Predefined Matrices
197
198<table class="tutorial_code">
199<tr>
200 <td>Fixed-size matrix or vector</td>
201 <td>Dynamic-size matrix</td>
202 <td>Dynamic-size vector</td>
203</tr>
204<tr style="border-bottom-style: none;">
205 <td>
206\code
207typedef {Matrix3f|Array33f} FixedXD;
208FixedXD x;
209
210x = FixedXD::Zero();
211x = FixedXD::Ones();
212x = FixedXD::Constant(value);
213x = FixedXD::Random();
214
215x.setZero();
216x.setOnes();
217x.setConstant(value);
218x.setRandom();
219\endcode
220 </td>
221 <td>
222\code
223typedef {MatrixXf|ArrayXXf} Dynamic2D;
224Dynamic2D x;
225
226x = Dynamic2D::Zero(rows, cols);
227x = Dynamic2D::Ones(rows, cols);
228x = Dynamic2D::Constant(rows, cols, value);
229x = Dynamic2D::Random(rows, cols);
230
231x.setZero(rows, cols);
232x.setOnes(rows, cols);
233x.setConstant(rows, cols, value);
234x.setRandom(rows, cols);
235\endcode
236 </td>
237 <td>
238\code
239typedef {VectorXf|ArrayXf} Dynamic1D;
240Dynamic1D x;
241
242x = Dynamic1D::Zero(size);
243x = Dynamic1D::Ones(size);
244x = Dynamic1D::Constant(size, value);
245x = Dynamic1D::Random(size);
246
247x.setZero(size);
248x.setOnes(size);
249x.setConstant(size, value);
250x.setRandom(size);
251\endcode
252 </td>
253</tr>
254
255<tr><td colspan="3">Identity and \link MatrixBase::Unit basis vectors \endlink \matrixworld</td></tr>
256<tr style="border-bottom-style: none;">
257 <td>
258\code
259x = FixedXD::Identity();
260x.setIdentity();
261
262Vector3f::UnitX() // 1 0 0
263Vector3f::UnitY() // 0 1 0
264Vector3f::UnitZ() // 0 0 1
265\endcode
266 </td>
267 <td>
268\code
269x = Dynamic2D::Identity(rows, cols);
270x.setIdentity(rows, cols);
271
272
273
274N/A
275\endcode
276 </td>
277 <td>\code
278N/A
279
280
281VectorXf::Unit(size,i)
282VectorXf::Unit(4,1) == Vector4f(0,1,0,0)
283 == Vector4f::UnitY()
284\endcode
285 </td>
286</tr>
287</table>
288
289
290
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200291\subsection QuickRef_Map Mapping external arrays
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200292
293<table class="tutorial_code">
294<tr>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200295<td>Contiguous \n memory</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200296<td>\code
297float data[] = {1,2,3,4};
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200298Map<Vector3f> v1(data); // uses v1 as a Vector3f object
299Map<ArrayXf> v2(data,3); // uses v2 as a ArrayXf object
300Map<Array22f> m1(data); // uses m1 as a Array22f object
301Map<MatrixXf> m2(data,2,2); // uses m2 as a MatrixXf object
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200302\endcode</td>
303</tr>
304<tr>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200305<td>Typical usage \n of strides</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200306<td>\code
307float data[] = {1,2,3,4,5,6,7,8,9};
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200308Map<VectorXf,0,InnerStride<2> > v1(data,3); // = [1,3,5]
309Map<VectorXf,0,InnerStride<> > v2(data,3,InnerStride<>(3)); // = [1,4,7]
310Map<MatrixXf,0,OuterStride<3> > m2(data,2,3); // both lines |1,4,7|
311Map<MatrixXf,0,OuterStride<> > m1(data,2,3,OuterStride<>(3)); // are equal to: |2,5,8|
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200312\endcode</td>
313</tr>
314</table>
315
316
317<a href="#" class="top">top</a>
318\section QuickRef_ArithmeticOperators Arithmetic Operators
319
320<table class="tutorial_code">
321<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200322add \n subtract</td><td>\code
323mat3 = mat1 + mat2; mat3 += mat1;
324mat3 = mat1 - mat2; mat3 -= mat1;\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200325</td></tr>
326<tr><td>
327scalar product</td><td>\code
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200328mat3 = mat1 * s1; mat3 *= s1; mat3 = s1 * mat1;
329mat3 = mat1 / s1; mat3 /= s1;\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200330</td></tr>
331<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200332matrix/vector \n products \matrixworld</td><td>\code
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200333col2 = mat1 * col1;
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200334row2 = row1 * mat1; row1 *= mat1;
335mat3 = mat1 * mat2; mat3 *= mat1; \endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200336</td></tr>
337<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200338transposition \n adjoint \matrixworld</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200339mat1 = mat2.transpose(); mat1.transposeInPlace();
340mat1 = mat2.adjoint(); mat1.adjointInPlace();
341\endcode
342</td></tr>
343<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200344\link MatrixBase::dot() dot \endlink product \n inner product \matrixworld</td><td>\code
345scalar = vec1.dot(vec2);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200346scalar = col1.adjoint() * col2;
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200347scalar = (col1.adjoint() * col2).value();\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200348</td></tr>
349<tr><td>
350outer product \matrixworld</td><td>\code
351mat = col1 * col2.transpose();\endcode
352</td></tr>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200353
354<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200355\link MatrixBase::norm() norm \endlink \n \link MatrixBase::normalized() normalization \endlink \matrixworld</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200356scalar = vec1.norm(); scalar = vec1.squaredNorm()
357vec2 = vec1.normalized(); vec1.normalize(); // inplace \endcode
358</td></tr>
359
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200360<tr><td>
361\link MatrixBase::cross() cross product \endlink \matrixworld</td><td>\code
362#include <Eigen/Geometry>
363vec3 = vec1.cross(vec2);\endcode</td></tr>
364</table>
365
366<a href="#" class="top">top</a>
367\section QuickRef_Coeffwise Coefficient-wise \& Array operators
368Coefficient-wise operators for matrices and vectors:
369<table class="tutorial_code">
370<tr><td>Matrix API \matrixworld</td><td>Via Array conversions</td></tr>
371<tr><td>\code
372mat1.cwiseMin(mat2)
373mat1.cwiseMax(mat2)
Gael Guennebauddbefd7a2010-06-28 13:30:10 +0200374mat1.cwiseAbs2()
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200375mat1.cwiseAbs()
376mat1.cwiseSqrt()
377mat1.cwiseProduct(mat2)
378mat1.cwiseQuotient(mat2)\endcode
379</td><td>\code
380mat1.array().min(mat2.array())
381mat1.array().max(mat2.array())
382mat1.array().abs2()
383mat1.array().abs()
384mat1.array().sqrt()
385mat1.array() * mat2.array()
386mat1.array() / mat2.array()
387\endcode</td></tr>
388</table>
389
390Array operators:\arrayworld
391
392<table class="tutorial_code">
393<tr><td>Arithmetic operators</td><td>\code
394array1 * array2 array1 / array2 array1 *= array2 array1 /= array2
395array1 + scalar array1 - scalar array1 += scalar array1 -= scalar
396\endcode</td></tr>
397<tr><td>Comparisons</td><td>\code
398array1 < array2 array1 > array2 array1 < scalar array1 > scalar
399array1 <= array2 array1 >= array2 array1 <= scalar array1 >= scalar
400array1 == array2 array1 != array2 array1 == scalar array1 != scalar
401\endcode</td></tr>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200402<tr><td>Trigo, power, and \n misc functions \n and the STL variants</td><td>\code
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200403array1.min(array2) std::min(array1,array2)
404array1.max(array2) std::max(array1,array2)
405array1.abs2()
406array1.abs() std::abs(array1)
407array1.sqrt() std::sqrt(array1)
408array1.log() std::log(array1)
409array1.exp() std::exp(array1)
410array1.pow(exponent) std::pow(array1,exponent)
411array1.square()
412array1.cube()
413array1.inverse()
414array1.sin() std::sin(array1)
415array1.cos() std::cos(array1)
416array1.tan() std::tan(array1)
417\endcode
418</td></tr>
419</table>
420
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200421<a href="#" class="top">top</a>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200422\section QuickRef_Reductions Reductions
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200423
424Eigen provides several reduction methods such as:
425\link DenseBase::minCoeff() minCoeff() \endlink, \link DenseBase::maxCoeff() maxCoeff() \endlink,
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200426\link DenseBase::sum() sum() \endlink, \link DenseBase::prod() prod() \endlink,
427\link MatrixBase::trace() trace() \endlink \matrixworld,
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200428\link MatrixBase::norm() norm() \endlink \matrixworld, \link MatrixBase::squaredNorm() squaredNorm() \endlink \matrixworld,
429\link DenseBase::all() all() \endlink \redstar,and \link DenseBase::any() any() \endlink \redstar.
430All reduction operations can be done matrix-wise,
431\link DenseBase::colwise() column-wise \endlink \redstar or
432\link DenseBase::rowwise() row-wise \endlink \redstar. Usage example:
433<table class="tutorial_code">
434<tr><td rowspan="3" style="border-right-style:dashed">\code
435 5 3 1
436mat = 2 7 8
437 9 4 6 \endcode
438</td> <td>\code mat.minCoeff(); \endcode</td><td>\code 1 \endcode</td></tr>
439<tr><td>\code mat.colwise().minCoeff(); \endcode</td><td>\code 2 3 1 \endcode</td></tr>
440<tr><td>\code mat.rowwise().minCoeff(); \endcode</td><td>\code
4411
4422
4434
444\endcode</td></tr>
445</table>
446
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200447Special versions of \link DenseBase::minCoeff(int*,int*) minCoeff \endlink and \link DenseBase::maxCoeff(int*,int*) maxCoeff \endlink:
448\code
449int i, j;
450s = vector.minCoeff(&i); // s == vector[i]
451s = matrix.maxCoeff(&i, &j); // s == matrix(i,j)
452\endcode
453Typical use cases of all() and any():
454\code
455if((array1 > 0).all()) ... // if all coefficients of array1 are greater than 0 ...
456if((array1 < array2).any()) ... // if there exist a pair i,j such that array1(i,j) < array2(i,j) ...
457\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200458
459
460
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200461<a href="#" class="top">top</a>\section QuickRef_Blocks Sub-matrices
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200462
463Read-write access to a \link DenseBase::col(int) column \endlink
464or a \link DenseBase::row(int) row \endlink of a matrix (or array):
465\code
466mat1.row(i) = mat2.col(j);
467mat1.col(j1).swap(mat1.col(j2));
468\endcode
469
470Read-write access to sub-vectors:
471<table class="tutorial_code">
472<tr>
473<td>Default versions</td>
474<td>Optimized versions when the size \n is known at compile time</td></tr>
475<td></td>
476
477<tr><td>\code vec1.head(n)\endcode</td><td>\code vec1.head<n>()\endcode</td><td>the first \c n coeffs </td></tr>
478<tr><td>\code vec1.tail(n)\endcode</td><td>\code vec1.tail<n>()\endcode</td><td>the last \c n coeffs </td></tr>
479<tr><td>\code vec1.segment(pos,n)\endcode</td><td>\code vec1.segment<n>(pos)\endcode</td>
480 <td>the \c size coeffs in \n the range [\c pos : \c pos + \c n [</td></tr>
481<tr style="border-style: dashed none dashed none;"><td>
482
483Read-write access to sub-matrices:</td><td></td><td></td></tr>
484<tr>
485 <td>\code mat1.block(i,j,rows,cols)\endcode
486 \link DenseBase::block(int,int,int,int) (more) \endlink</td>
487 <td>\code mat1.block<rows,cols>(i,j)\endcode
488 \link DenseBase::block(int,int) (more) \endlink</td>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200489 <td>the \c rows x \c cols sub-matrix \n starting from position (\c i,\c j)</td></tr>
490<tr><td>\code
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200491 mat1.topLeftCorner(rows,cols)
492 mat1.topRightCorner(rows,cols)
493 mat1.bottomLeftCorner(rows,cols)
494 mat1.bottomRightCorner(rows,cols)\endcode
495 <td>\code
496 mat1.topLeftCorner<rows,cols>()
497 mat1.topRightCorner<rows,cols>()
498 mat1.bottomLeftCorner<rows,cols>()
499 mat1.bottomRightCorner<rows,cols>()\endcode
500 <td>the \c rows x \c cols sub-matrix \n taken in one of the four corners</td></tr>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200501 <tr><td>\code
502 mat1.topRows(rows)
503 mat1.bottomRows(rows)
504 mat1.leftCols(cols)
505 mat1.rightCols(cols)\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200506 <td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200507 mat1.topRows<rows>()
508 mat1.bottomRows<rows>()
509 mat1.leftCols<cols>()
510 mat1.rightCols<cols>()\endcode
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200511 <td>specialized versions of block() \n when the block fit two corners</td></tr>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200512</table>
513
514
515
516
517<a href="#" class="top">top</a>\section QuickRef_DiagTriSymm Diagonal, Triangular, and Self-adjoint matrices
518(matrix world \matrixworld)
519
520\subsection QuickRef_Diagonal Diagonal matrices
521
522<table class="tutorial_code">
523<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200524view a vector \link MatrixBase::asDiagonal() as a diagonal matrix \endlink \n </td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200525mat1 = vec1.asDiagonal();\endcode
526</td></tr>
527<tr><td>
528Declare a diagonal matrix</td><td>\code
529DiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
530diag1.diagonal() = vector;\endcode
531</td></tr>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200532<tr><td>Access the \link MatrixBase::diagonal() diagonal \endlink and \link MatrixBase::diagonal(int) super/sub diagonals \endlink of a matrix as a vector (read/write)</td>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200533 <td>\code
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200534vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal
535vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal
536vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal
537vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal
538vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200539\endcode</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200540</tr>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200541
542<tr><td>Optimized products and inverse</td>
543 <td>\code
544mat3 = scalar * diag1 * mat1;
545mat3 += scalar * mat1 * vec1.asDiagonal();
546mat3 = vec1.asDiagonal().inverse() * mat1
547mat3 = mat1 * diag1.inverse()
548\endcode</td>
549</tr>
550
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200551</table>
552
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200553\subsection QuickRef_TriangularView Triangular views
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200554
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200555TriangularView gives a view on a triangular part of a dense matrix and allows to perform optimized operations on it. The opposite triangular part is never referenced and can be used to store other information.
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200556
557<table class="tutorial_code">
558<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200559Reference to a triangular with optional \n
560unit or null diagonal (read/write):
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200561</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200562m.triangularView<Xxx>()
563\endcode \n
564\c Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200565</td></tr>
566<tr><td>
567Writing to a specific triangular part:\n (only the referenced triangular part is evaluated)
568</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200569m1.triangularView<Eigen::Lower>() = m2 + m3 \endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200570</td></tr>
571<tr><td>
572Conversion to a dense matrix setting the opposite triangular part to zero:
573</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200574m2 = m1.triangularView<Eigen::UnitUpper>()\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200575</td></tr>
576<tr><td>
577Products:
578</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200579m3 += s1 * m1.adjoint().triangularView<Eigen::UnitUpper>() * m2
580m3 -= s1 * m2.conjugate() * m1.adjoint().triangularView<Eigen::Lower>() \endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200581</td></tr>
582<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200583Solving linear equations:\n
584\f$ M_2 := L_1^{-1} M_2 \f$ \n
585\f$ M_3 := {L_1^*}^{-1} M_3 \f$ \n
586\f$ M_4 := M_3 U_1^{-1} \f$
587</td><td>\n \code
588L1.triangularView<Eigen::UnitLower>().solveInPlace(M2)
589L1.triangularView<Eigen::Lower>().adjoint().solveInPlace(M3)
590U1.triangularView<Eigen::Upper>().solveInPlace<OnTheRight>(M4)\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200591</td></tr>
592</table>
593
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200594\subsection QuickRef_SelfadjointMatrix Symmetric/selfadjoint views
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200595
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200596Just as for triangular matrix, you can reference any triangular part of a square matrix to see it as a selfadjoint
597matrix and perform special and optimized operations. Again the opposite triangular part is never referenced and can be
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200598used to store other information.
599
600<table class="tutorial_code">
601<tr><td>
602Conversion to a dense matrix:
603</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200604m2 = m.selfadjointView<Eigen::Lower>();\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200605</td></tr>
606<tr><td>
607Product with another general matrix or vector:
608</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200609m3 = s1 * m1.conjugate().selfadjointView<Eigen::Upper>() * m3;
610m3 -= s1 * m3.adjoint() * m1.selfadjointView<Eigen::Lower>();\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200611</td></tr>
612<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200613Rank 1 and rank K update: \n
614\f$ upper(M_1) += s1 M_2^* M_2 \f$ \n
615\f$ lower(M_1) -= M_2 M_2^* \f$
616</td><td>\n \code
617M1.selfadjointView<Eigen::Upper>().rankUpdate(M2,s1);
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200618m1.selfadjointView<Eigen::Lower>().rankUpdate(m2.adjoint(),-1); \endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200619</td></tr>
620<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200621Rank 2 update: (\f$ M += s u v^* + s v u^* \f$)
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200622</td><td>\code
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200623M.selfadjointView<Eigen::Upper>().rankUpdate(u,v,s);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200624\endcode
625</td></tr>
626<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200627Solving linear equations:\n(\f$ M_2 := M_1^{-1} M_2 \f$)
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200628</td><td>\code
629// via a standard Cholesky factorization
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200630m2 = m1.selfadjointView<Eigen::Upper>().llt().solve(m2);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200631// via a Cholesky factorization with pivoting
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200632m2 = m1.selfadjointView<Eigen::Lower>().ldlt().solve(m2);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200633\endcode
634</td></tr>
635</table>
636
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200637*/
638
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200639/*
640<table class="tutorial_code">
641<tr><td>
642\link MatrixBase::asDiagonal() make a diagonal matrix \endlink \n from a vector </td><td>\code
643mat1 = vec1.asDiagonal();\endcode
644</td></tr>
645<tr><td>
646Declare a diagonal matrix</td><td>\code
647DiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
648diag1.diagonal() = vector;\endcode
649</td></tr>
650<tr><td>Access \link MatrixBase::diagonal() the diagonal and super/sub diagonals of a matrix \endlink as a vector (read/write)</td>
651 <td>\code
652vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal
653vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal
654vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal
655vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal
656vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal
657\endcode</td>
658</tr>
659
660<tr><td>View on a triangular part of a matrix (read/write)</td>
661 <td>\code
662mat2 = mat1.triangularView<Xxx>();
663// Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower
664mat1.triangularView<Upper>() = mat2 + mat3; // only the upper part is evaluated and referenced
665\endcode</td></tr>
666
667<tr><td>View a triangular part as a symmetric/self-adjoint matrix (read/write)</td>
668 <td>\code
669mat2 = mat1.selfadjointView<Xxx>(); // Xxx = Upper or Lower
670mat1.selfadjointView<Upper>() = mat2 + mat2.adjoint(); // evaluated and write to the upper triangular part only
671\endcode</td></tr>
672
673</table>
674
675Optimized products:
676\code
677mat3 += scalar * vec1.asDiagonal() * mat1
678mat3 += scalar * mat1 * vec1.asDiagonal()
679mat3.noalias() += scalar * mat1.triangularView<Xxx>() * mat2
680mat3.noalias() += scalar * mat2 * mat1.triangularView<Xxx>()
681mat3.noalias() += scalar * mat1.selfadjointView<Upper or Lower>() * mat2
682mat3.noalias() += scalar * mat2 * mat1.selfadjointView<Upper or Lower>()
683mat1.selfadjointView<Upper or Lower>().rankUpdate(mat2);
684mat1.selfadjointView<Upper or Lower>().rankUpdate(mat2.adjoint(), scalar);
685\endcode
686
687Inverse products: (all are optimized)
688\code
689mat3 = vec1.asDiagonal().inverse() * mat1
690mat3 = mat1 * diag1.inverse()
691mat1.triangularView<Xxx>().solveInPlace(mat2)
692mat1.triangularView<Xxx>().solveInPlace<OnTheRight>(mat2)
693mat2 = mat1.selfadjointView<Upper or Lower>().llt().solve(mat2)
694\endcode
695
696*/
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200697}