blob: 467899fb588937ebe9c9d24539146a16e9219e71 [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
Jitse Niesen49747fa2010-07-06 13:10:08 +010021The Eigen library is divided in a Core module and several additional modules. Each module has a corresponding header file which has to be included in order to use the module. The \c %Dense and \c Eigen header files are provided to conveniently gain access to several modules at once.
22
Gael Guennebaudca4bd582010-10-19 11:59:11 +020023<table class="manual">
24<tr><th>Module</th><th>Header file</th><th>Contents</th></tr>
Jitse Niesen49747fa2010-07-06 13:10:08 +010025<tr><td>\link Core_Module Core \endlink</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>
Gael Guennebaudca4bd582010-10-19 11:59:11 +020026<tr class="alt"><td>\link Geometry_Module Geometry \endlink</td><td>\code#include <Eigen/Geometry>\endcode</td><td>Transform, Translation, Scaling, Rotation2D and 3D rotations (Quaternion, AngleAxis)</td></tr>
Jitse Niesen49747fa2010-07-06 13:10:08 +010027<tr><td>\link LU_Module LU \endlink</td><td>\code#include <Eigen/LU>\endcode</td><td>Inverse, determinant, LU decompositions with solver (FullPivLU, PartialPivLU)</td></tr>
28<tr><td>\link Cholesky_Module Cholesky \endlink</td><td>\code#include <Eigen/Cholesky>\endcode</td><td>LLT and LDLT Cholesky factorization with solver</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +020029<tr class="alt"><td>\link Householder_Module Householder \endlink</td><td>\code#include <Eigen/Householder>\endcode</td><td>Householder transformations; this module is used by several linear algebra modules</td></tr>
Benoit Jacob26129222010-10-15 09:44:43 -040030<tr><td>\link SVD_Module SVD \endlink</td><td>\code#include <Eigen/SVD>\endcode</td><td>SVD decomposition with least-squares solver (JacobiSVD)</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +020031<tr class="alt"><td>\link QR_Module QR \endlink</td><td>\code#include <Eigen/QR>\endcode</td><td>QR decomposition with solver (HouseholderQR, ColPivHouseholderQR, FullPivHouseholderQR)</td></tr>
Jitse Niesen49747fa2010-07-06 13:10:08 +010032<tr><td>\link Eigenvalues_Module Eigenvalues \endlink</td><td>\code#include <Eigen/Eigenvalues>\endcode</td><td>Eigenvalue, eigenvector decompositions (EigenSolver, SelfAdjointEigenSolver, ComplexEigenSolver)</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +020033<tr class="alt"><td>\link Sparse_Module Sparse \endlink</td><td>\code#include <Eigen/Sparse>\endcode</td><td>%Sparse matrix storage and related basic linear algebra (SparseMatrix, DynamicSparseMatrix, SparseVector)</td></tr>
Benoit Jacob26129222010-10-15 09:44:43 -040034<tr><td></td><td>\code#include <Eigen/Dense>\endcode</td><td>Includes Core, Geometry, LU, Cholesky, SVD, QR, and Eigenvalues header files</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +020035<tr class="alt"><td></td><td>\code#include <Eigen/Eigen>\endcode</td><td>Includes %Dense and %Sparse header files (the whole Eigen library)</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +020036</table>
37
38<a href="#" class="top">top</a>
39\section QuickRef_Types Array, matrix and vector types
40
Gael Guennebauddbefd7a2010-06-28 13:30:10 +020041
Gael Guennebaud5c866f22010-06-26 16:59:18 +020042\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:
Hauke Heibel7b74d372010-06-30 18:27:27 +020043<div class="desired_tutorial_width">
Gael Guennebaud5c866f22010-06-26 16:59:18 +020044\code
45typedef Matrix<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyMatrixType;
46typedef Array<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyArrayType;
47\endcode
Hauke Heibel7b74d372010-06-30 18:27:27 +020048</div>
Gael Guennebaud5c866f22010-06-26 16:59:18 +020049
50\li \c Scalar is the scalar type of the coefficients (e.g., \c float, \c double, \c bool, \c int, etc.).
51\li \c RowsAtCompileTime and \c ColsAtCompileTime are the number of rows and columns of the matrix as known at compile-time or \c Dynamic.
52\li \c Options can be \c ColMajor or \c RowMajor, default is \c ColMajor. (see class Matrix for more options)
53
54All 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:
Hauke Heibel7b74d372010-06-30 18:27:27 +020055<div class="desired_tutorial_width">
Gael Guennebaud5c866f22010-06-26 16:59:18 +020056\code
57Matrix<double, 6, Dynamic> // Dynamic number of columns (heap allocation)
58Matrix<double, Dynamic, 2> // Dynamic number of rows (heap allocation)
59Matrix<double, Dynamic, Dynamic, RowMajor> // Fully dynamic, row major (heap allocation)
60Matrix<double, 13, 3> // Fully fixed (static allocation)
61\endcode
Hauke Heibel7b74d372010-06-30 18:27:27 +020062</div>
Gael Guennebaud5c866f22010-06-26 16:59:18 +020063
64In most cases, you can simply use one of the convenience typedefs for \ref matrixtypedefs "matrices" and \ref arraytypedefs "arrays". Some examples:
65<table class="tutorial_code">
66<tr><td>\code
67Matrix<float,Dynamic,Dynamic> <=> MatrixXf
68Matrix<double,Dynamic,1> <=> VectorXd
69Matrix<int,1,Dynamic> <=> RowVectorXi
70Matrix<float,3,3> <=> Matrix3f
71Matrix<float,4,1> <=> Vector4f
72\endcode</td><td>\code
73Array<float,Dynamic,Dynamic> <=> ArrayXXf
74Array<double,Dynamic,1> <=> ArrayXd
75Array<int,1,Dynamic> <=> RowArrayXi
76Array<float,3,3> <=> Array33f
77Array<float,4,1> <=> Array4f
78\endcode</td></tr>
79</table>
80
81Conversion between the matrix and array worlds:
Hauke Heibel7b74d372010-06-30 18:27:27 +020082<div class="desired_tutorial_width">
Gael Guennebaud5c866f22010-06-26 16:59:18 +020083\code
84Array44f a1, a1;
85Matrix4f m1, m2;
Gael Guennebaudf3c64c72010-06-26 22:19:03 +020086m1 = a1 * a2; // coeffwise product, implicit conversion from array to matrix.
87a1 = m1 * m2; // matrix product, implicit conversion from matrix to array.
88a2 = a1 + m1.array(); // mixing array and matrix is forbidden
89m2 = a1.matrix() + m1; // and explicit conversion is required.
90ArrayWrapper<Matrix4f> m1a(m1); // m1a is an alias for m1.array(), they share the same coefficients
Gael Guennebaud5c866f22010-06-26 16:59:18 +020091MatrixWrapper<Array44f> a1m(a1);
92\endcode
Hauke Heibel7b74d372010-06-30 18:27:27 +020093</div>
Gael Guennebaud5c866f22010-06-26 16:59:18 +020094
95In the rest of this document we will use the following symbols to emphasize the features which are specifics to a given kind of object:
96\li <a name="matrixonly"><a/>\matrixworld linear algebra matrix and vector only
97\li <a name="arrayonly"><a/>\arrayworld array objects only
98
Gael Guennebaud5c866f22010-06-26 16:59:18 +020099\subsection QuickRef_Basics Basic matrix manipulation
100
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200101<table class="manual">
102<tr><th></th><th>1D objects</th><th>2D objects</th><th>Notes</th></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200103<tr><td>Constructors</td>
104<td>\code
105Vector4d v4;
106Vector2f v1(x, y);
107Array3i v2(x, y, z);
108Vector4d v3(x, y, z, w);
109
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200110VectorXf v5; // empty object
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200111ArrayXf v6(size);
112\endcode</td><td>\code
113Matrix4f m1;
114
115
116
117
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200118MatrixXf m5; // empty object
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200119MatrixXf m6(nb_rows, nb_columns);
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200120\endcode</td><td class="note">
121By default, the coefficients \n are left uninitialized</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200122<tr class="alt"><td>Comma initializer</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200123<td>\code
124Vector3f v1; v1 << x, y, z;
125ArrayXf v2(4); v2 << 1, 2, 3, 4;
126
127\endcode</td><td>\code
128Matrix3f m1; m1 << 1, 2, 3,
129 4, 5, 6,
130 7, 8, 9;
131\endcode</td><td></td></tr>
132
133<tr><td>Comma initializer (bis)</td>
134<td colspan="2">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200135\include Tutorial_commainit_02.cpp
136</td>
137<td>
138output:
139\verbinclude Tutorial_commainit_02.out
Hauke Heibele1348b92010-06-30 18:04:36 +0200140</td>
141</tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200142
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200143<tr class="alt"><td>Runtime info</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200144<td>\code
145vector.size();
146
147vector.innerStride();
148vector.data();
149\endcode</td><td>\code
150matrix.rows(); matrix.cols();
151matrix.innerSize(); matrix.outerSize();
152matrix.innerStride(); matrix.outerStride();
153matrix.data();
Hauke Heibelb1741c12010-06-30 15:52:00 +0200154\endcode</td><td class="note">Inner/Outer* are storage order dependent</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200155<tr><td>Compile-time info</td>
156<td colspan="2">\code
157ObjectType::Scalar ObjectType::RowsAtCompileTime
158ObjectType::RealScalar ObjectType::ColsAtCompileTime
159ObjectType::Index ObjectType::SizeAtCompileTime
160\endcode</td><td></td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200161<tr class="alt"><td>Resizing</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200162<td>\code
163vector.resize(size);
164
165
166vector.resizeLike(other_vector);
167vector.conservativeResize(size);
168\endcode</td><td>\code
169matrix.resize(nb_rows, nb_cols);
170matrix.resize(Eigen::NoChange, nb_cols);
171matrix.resize(nb_rows, Eigen::NoChange);
172matrix.resizeLike(other_matrix);
173matrix.conservativeResize(nb_rows, nb_cols);
Hauke Heibelb1741c12010-06-30 15:52:00 +0200174\endcode</td><td class="note">no-op if the new sizes match,<br/>otherwise data are lost<br/><br/>resizing with data preservation</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200175
176<tr><td>Coeff access with \n range checking</td>
177<td>\code
178vector(i) vector.x()
179vector[i] vector.y()
180 vector.z()
181 vector.w()
182\endcode</td><td>\code
183matrix(i,j)
Jitse Niesen76fbe942010-08-10 11:37:23 +0100184\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 +0200185
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200186<tr class="alt"><td>Coeff access without \n range checking</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200187<td>\code
188vector.coeff(i)
189vector.coeffRef(i)
190\endcode</td><td>\code
191matrix.coeff(i,j)
192matrix.coeffRef(i,j)
193\endcode</td><td></td></tr>
194
195<tr><td>Assignment/copy</td>
196<td colspan="2">\code
197object = expression;
198object_of_float = expression_of_double.cast<float>();
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200199\endcode</td><td class="note">the destination is automatically resized (if possible)</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200200
201</table>
202
203\subsection QuickRef_PredefMat Predefined Matrices
204
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200205<table class="manual">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200206<tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200207 <th>Fixed-size matrix or vector</th>
208 <th>Dynamic-size matrix</th>
209 <th>Dynamic-size vector</th>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200210</tr>
211<tr style="border-bottom-style: none;">
212 <td>
213\code
214typedef {Matrix3f|Array33f} FixedXD;
215FixedXD x;
216
217x = FixedXD::Zero();
218x = FixedXD::Ones();
219x = FixedXD::Constant(value);
220x = FixedXD::Random();
Jitse Niesen3abbdfd2010-07-20 21:55:22 +0100221x = FixedXD::LinSpaced(low, high, size);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200222
223x.setZero();
224x.setOnes();
225x.setConstant(value);
226x.setRandom();
Jitse Niesen3abbdfd2010-07-20 21:55:22 +0100227x.setLinSpaced(low, high, size);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200228\endcode
229 </td>
230 <td>
231\code
232typedef {MatrixXf|ArrayXXf} Dynamic2D;
233Dynamic2D x;
234
235x = Dynamic2D::Zero(rows, cols);
236x = Dynamic2D::Ones(rows, cols);
237x = Dynamic2D::Constant(rows, cols, value);
238x = Dynamic2D::Random(rows, cols);
Jitse Niesen3abbdfd2010-07-20 21:55:22 +0100239N/A
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200240
241x.setZero(rows, cols);
242x.setOnes(rows, cols);
243x.setConstant(rows, cols, value);
244x.setRandom(rows, cols);
Jitse Niesen3abbdfd2010-07-20 21:55:22 +0100245N/A
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200246\endcode
247 </td>
248 <td>
249\code
250typedef {VectorXf|ArrayXf} Dynamic1D;
251Dynamic1D x;
252
253x = Dynamic1D::Zero(size);
254x = Dynamic1D::Ones(size);
255x = Dynamic1D::Constant(size, value);
256x = Dynamic1D::Random(size);
Jitse Niesen3abbdfd2010-07-20 21:55:22 +0100257x = Dynamic1D::LinSpaced(low, high, size);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200258
259x.setZero(size);
260x.setOnes(size);
261x.setConstant(size, value);
262x.setRandom(size);
Jitse Niesen3abbdfd2010-07-20 21:55:22 +0100263x.setLinSpaced(low, high, size);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200264\endcode
265 </td>
266</tr>
267
268<tr><td colspan="3">Identity and \link MatrixBase::Unit basis vectors \endlink \matrixworld</td></tr>
269<tr style="border-bottom-style: none;">
270 <td>
271\code
272x = FixedXD::Identity();
273x.setIdentity();
274
275Vector3f::UnitX() // 1 0 0
276Vector3f::UnitY() // 0 1 0
277Vector3f::UnitZ() // 0 0 1
278\endcode
279 </td>
280 <td>
281\code
282x = Dynamic2D::Identity(rows, cols);
283x.setIdentity(rows, cols);
284
285
286
287N/A
288\endcode
289 </td>
290 <td>\code
291N/A
292
293
294VectorXf::Unit(size,i)
295VectorXf::Unit(4,1) == Vector4f(0,1,0,0)
296 == Vector4f::UnitY()
297\endcode
298 </td>
299</tr>
300</table>
301
302
303
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200304\subsection QuickRef_Map Mapping external arrays
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200305
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200306<table class="manual">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200307<tr>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200308<td>Contiguous \n memory</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200309<td>\code
310float data[] = {1,2,3,4};
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200311Map<Vector3f> v1(data); // uses v1 as a Vector3f object
312Map<ArrayXf> v2(data,3); // uses v2 as a ArrayXf object
313Map<Array22f> m1(data); // uses m1 as a Array22f object
314Map<MatrixXf> m2(data,2,2); // uses m2 as a MatrixXf object
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200315\endcode</td>
316</tr>
317<tr>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200318<td>Typical usage \n of strides</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200319<td>\code
320float data[] = {1,2,3,4,5,6,7,8,9};
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200321Map<VectorXf,0,InnerStride<2> > v1(data,3); // = [1,3,5]
322Map<VectorXf,0,InnerStride<> > v2(data,3,InnerStride<>(3)); // = [1,4,7]
323Map<MatrixXf,0,OuterStride<3> > m2(data,2,3); // both lines |1,4,7|
324Map<MatrixXf,0,OuterStride<> > m1(data,2,3,OuterStride<>(3)); // are equal to: |2,5,8|
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200325\endcode</td>
326</tr>
327</table>
328
329
330<a href="#" class="top">top</a>
331\section QuickRef_ArithmeticOperators Arithmetic Operators
332
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200333<table class="manual">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200334<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200335add \n subtract</td><td>\code
336mat3 = mat1 + mat2; mat3 += mat1;
337mat3 = mat1 - mat2; mat3 -= mat1;\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200338</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200339<tr class="alt"><td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200340scalar product</td><td>\code
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200341mat3 = mat1 * s1; mat3 *= s1; mat3 = s1 * mat1;
342mat3 = mat1 / s1; mat3 /= s1;\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200343</td></tr>
344<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200345matrix/vector \n products \matrixworld</td><td>\code
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200346col2 = mat1 * col1;
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200347row2 = row1 * mat1; row1 *= mat1;
348mat3 = mat1 * mat2; mat3 *= mat1; \endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200349</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200350<tr class="alt"><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200351transposition \n adjoint \matrixworld</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200352mat1 = mat2.transpose(); mat1.transposeInPlace();
353mat1 = mat2.adjoint(); mat1.adjointInPlace();
354\endcode
355</td></tr>
356<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200357\link MatrixBase::dot() dot \endlink product \n inner product \matrixworld</td><td>\code
358scalar = vec1.dot(vec2);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200359scalar = col1.adjoint() * col2;
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200360scalar = (col1.adjoint() * col2).value();\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200361</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200362<tr class="alt"><td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200363outer product \matrixworld</td><td>\code
364mat = col1 * col2.transpose();\endcode
365</td></tr>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200366
367<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200368\link MatrixBase::norm() norm \endlink \n \link MatrixBase::normalized() normalization \endlink \matrixworld</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200369scalar = vec1.norm(); scalar = vec1.squaredNorm()
370vec2 = vec1.normalized(); vec1.normalize(); // inplace \endcode
371</td></tr>
372
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200373<tr class="alt"><td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200374\link MatrixBase::cross() cross product \endlink \matrixworld</td><td>\code
375#include <Eigen/Geometry>
376vec3 = vec1.cross(vec2);\endcode</td></tr>
377</table>
378
379<a href="#" class="top">top</a>
380\section QuickRef_Coeffwise Coefficient-wise \& Array operators
381Coefficient-wise operators for matrices and vectors:
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200382<table class="manual">
383<tr><th>Matrix API \matrixworld</th><th>Via Array conversions</th></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200384<tr><td>\code
385mat1.cwiseMin(mat2)
386mat1.cwiseMax(mat2)
Gael Guennebauddbefd7a2010-06-28 13:30:10 +0200387mat1.cwiseAbs2()
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200388mat1.cwiseAbs()
389mat1.cwiseSqrt()
390mat1.cwiseProduct(mat2)
391mat1.cwiseQuotient(mat2)\endcode
392</td><td>\code
393mat1.array().min(mat2.array())
394mat1.array().max(mat2.array())
395mat1.array().abs2()
396mat1.array().abs()
397mat1.array().sqrt()
398mat1.array() * mat2.array()
399mat1.array() / mat2.array()
400\endcode</td></tr>
401</table>
402
403Array operators:\arrayworld
404
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200405<table class="manual">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200406<tr><td>Arithmetic operators</td><td>\code
407array1 * array2 array1 / array2 array1 *= array2 array1 /= array2
408array1 + scalar array1 - scalar array1 += scalar array1 -= scalar
409\endcode</td></tr>
410<tr><td>Comparisons</td><td>\code
411array1 < array2 array1 > array2 array1 < scalar array1 > scalar
412array1 <= array2 array1 >= array2 array1 <= scalar array1 >= scalar
413array1 == array2 array1 != array2 array1 == scalar array1 != scalar
414\endcode</td></tr>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200415<tr><td>Trigo, power, and \n misc functions \n and the STL variants</td><td>\code
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200416array1.min(array2) std::min(array1,array2)
417array1.max(array2) std::max(array1,array2)
418array1.abs2()
419array1.abs() std::abs(array1)
420array1.sqrt() std::sqrt(array1)
421array1.log() std::log(array1)
422array1.exp() std::exp(array1)
423array1.pow(exponent) std::pow(array1,exponent)
424array1.square()
425array1.cube()
426array1.inverse()
427array1.sin() std::sin(array1)
428array1.cos() std::cos(array1)
429array1.tan() std::tan(array1)
430\endcode
431</td></tr>
432</table>
433
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200434<a href="#" class="top">top</a>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200435\section QuickRef_Reductions Reductions
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200436
437Eigen provides several reduction methods such as:
438\link DenseBase::minCoeff() minCoeff() \endlink, \link DenseBase::maxCoeff() maxCoeff() \endlink,
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200439\link DenseBase::sum() sum() \endlink, \link DenseBase::prod() prod() \endlink,
440\link MatrixBase::trace() trace() \endlink \matrixworld,
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200441\link MatrixBase::norm() norm() \endlink \matrixworld, \link MatrixBase::squaredNorm() squaredNorm() \endlink \matrixworld,
442\link DenseBase::all() all() \endlink \redstar,and \link DenseBase::any() any() \endlink \redstar.
443All reduction operations can be done matrix-wise,
444\link DenseBase::colwise() column-wise \endlink \redstar or
445\link DenseBase::rowwise() row-wise \endlink \redstar. Usage example:
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200446<table class="manual">
447<tr><td rowspan="3" style="border-right-style:dashed;vertical-align:middle">\code
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200448 5 3 1
449mat = 2 7 8
450 9 4 6 \endcode
451</td> <td>\code mat.minCoeff(); \endcode</td><td>\code 1 \endcode</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200452<tr class="alt"><td>\code mat.colwise().minCoeff(); \endcode</td><td>\code 2 3 1 \endcode</td></tr>
453<tr style="vertical-align:middle"><td>\code mat.rowwise().minCoeff(); \endcode</td><td>\code
Gael Guennebaud5c866f22010-06-26 16:59:18 +02004541
4552
4564
457\endcode</td></tr>
458</table>
459
Jitse Niesen1420f8b2010-07-25 20:29:07 +0100460Special versions of \link DenseBase::minCoeff(Index*,Index*) minCoeff \endlink and \link DenseBase::maxCoeff(Index*,Index*) maxCoeff \endlink:
Hauke Heibel7b74d372010-06-30 18:27:27 +0200461<div class="desired_tutorial_width">
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200462\code
463int i, j;
464s = vector.minCoeff(&i); // s == vector[i]
465s = matrix.maxCoeff(&i, &j); // s == matrix(i,j)
466\endcode
Hauke Heibel7b74d372010-06-30 18:27:27 +0200467</div>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200468Typical use cases of all() and any():
Hauke Heibel7b74d372010-06-30 18:27:27 +0200469<div class="desired_tutorial_width">
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200470\code
471if((array1 > 0).all()) ... // if all coefficients of array1 are greater than 0 ...
472if((array1 < array2).any()) ... // if there exist a pair i,j such that array1(i,j) < array2(i,j) ...
473\endcode
Hauke Heibel7b74d372010-06-30 18:27:27 +0200474</div>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200475
476
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200477<a href="#" class="top">top</a>\section QuickRef_Blocks Sub-matrices
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200478
Jitse Niesen1420f8b2010-07-25 20:29:07 +0100479Read-write access to a \link DenseBase::col(Index) column \endlink
480or a \link DenseBase::row(Index) row \endlink of a matrix (or array):
Hauke Heibel7b74d372010-06-30 18:27:27 +0200481<div class="desired_tutorial_width">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200482\code
483mat1.row(i) = mat2.col(j);
484mat1.col(j1).swap(mat1.col(j2));
485\endcode
Hauke Heibel7b74d372010-06-30 18:27:27 +0200486</div>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200487
488Read-write access to sub-vectors:
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200489<table class="manual">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200490<tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200491<th>Default versions</th>
492<th>Optimized versions when the size \n is known at compile time</th></tr>
493<th></th>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200494
495<tr><td>\code vec1.head(n)\endcode</td><td>\code vec1.head<n>()\endcode</td><td>the first \c n coeffs </td></tr>
496<tr><td>\code vec1.tail(n)\endcode</td><td>\code vec1.tail<n>()\endcode</td><td>the last \c n coeffs </td></tr>
497<tr><td>\code vec1.segment(pos,n)\endcode</td><td>\code vec1.segment<n>(pos)\endcode</td>
Jitse Niesen49747fa2010-07-06 13:10:08 +0100498 <td>the \c n coeffs in \n the range [\c pos : \c pos + \c n [</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200499<tr class="alt"><td colspan="3">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200500
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200501Read-write access to sub-matrices:</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200502<tr>
503 <td>\code mat1.block(i,j,rows,cols)\endcode
Jitse Niesen1420f8b2010-07-25 20:29:07 +0100504 \link DenseBase::block(Index,Index,Index,Index) (more) \endlink</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200505 <td>\code mat1.block<rows,cols>(i,j)\endcode
Jitse Niesen1420f8b2010-07-25 20:29:07 +0100506 \link DenseBase::block(Index,Index) (more) \endlink</td>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200507 <td>the \c rows x \c cols sub-matrix \n starting from position (\c i,\c j)</td></tr>
508<tr><td>\code
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200509 mat1.topLeftCorner(rows,cols)
510 mat1.topRightCorner(rows,cols)
511 mat1.bottomLeftCorner(rows,cols)
512 mat1.bottomRightCorner(rows,cols)\endcode
513 <td>\code
514 mat1.topLeftCorner<rows,cols>()
515 mat1.topRightCorner<rows,cols>()
516 mat1.bottomLeftCorner<rows,cols>()
517 mat1.bottomRightCorner<rows,cols>()\endcode
518 <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 +0200519 <tr><td>\code
520 mat1.topRows(rows)
521 mat1.bottomRows(rows)
522 mat1.leftCols(cols)
523 mat1.rightCols(cols)\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200524 <td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200525 mat1.topRows<rows>()
526 mat1.bottomRows<rows>()
527 mat1.leftCols<cols>()
528 mat1.rightCols<cols>()\endcode
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200529 <td>specialized versions of block() \n when the block fit two corners</td></tr>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200530</table>
531
532
533
534
535<a href="#" class="top">top</a>\section QuickRef_DiagTriSymm Diagonal, Triangular, and Self-adjoint matrices
536(matrix world \matrixworld)
537
538\subsection QuickRef_Diagonal Diagonal matrices
539
540<table class="tutorial_code">
541<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200542view a vector \link MatrixBase::asDiagonal() as a diagonal matrix \endlink \n </td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200543mat1 = vec1.asDiagonal();\endcode
544</td></tr>
545<tr><td>
546Declare a diagonal matrix</td><td>\code
547DiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
548diag1.diagonal() = vector;\endcode
549</td></tr>
Jitse Niesen1420f8b2010-07-25 20:29:07 +0100550<tr><td>Access the \link MatrixBase::diagonal() diagonal \endlink and \link MatrixBase::diagonal(Index) super/sub diagonals \endlink of a matrix as a vector (read/write)</td>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200551 <td>\code
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200552vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal
553vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal
554vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal
555vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal
556vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200557\endcode</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200558</tr>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200559
560<tr><td>Optimized products and inverse</td>
561 <td>\code
562mat3 = scalar * diag1 * mat1;
563mat3 += scalar * mat1 * vec1.asDiagonal();
564mat3 = vec1.asDiagonal().inverse() * mat1
565mat3 = mat1 * diag1.inverse()
566\endcode</td>
567</tr>
568
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200569</table>
570
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200571\subsection QuickRef_TriangularView Triangular views
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200572
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200573TriangularView 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 +0200574
575<table class="tutorial_code">
576<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200577Reference to a triangular with optional \n
578unit or null diagonal (read/write):
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200579</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200580m.triangularView<Xxx>()
581\endcode \n
582\c Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200583</td></tr>
584<tr><td>
585Writing to a specific triangular part:\n (only the referenced triangular part is evaluated)
586</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200587m1.triangularView<Eigen::Lower>() = m2 + m3 \endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200588</td></tr>
589<tr><td>
590Conversion to a dense matrix setting the opposite triangular part to zero:
591</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200592m2 = m1.triangularView<Eigen::UnitUpper>()\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200593</td></tr>
594<tr><td>
595Products:
596</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200597m3 += s1 * m1.adjoint().triangularView<Eigen::UnitUpper>() * m2
598m3 -= s1 * m2.conjugate() * m1.adjoint().triangularView<Eigen::Lower>() \endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200599</td></tr>
600<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200601Solving linear equations:\n
602\f$ M_2 := L_1^{-1} M_2 \f$ \n
603\f$ M_3 := {L_1^*}^{-1} M_3 \f$ \n
Thomas Capricelli8e21cef2010-07-22 13:15:15 +0200604\f$ M_4 := M_4 U_1^{-1} \f$
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200605</td><td>\n \code
606L1.triangularView<Eigen::UnitLower>().solveInPlace(M2)
607L1.triangularView<Eigen::Lower>().adjoint().solveInPlace(M3)
608U1.triangularView<Eigen::Upper>().solveInPlace<OnTheRight>(M4)\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200609</td></tr>
610</table>
611
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200612\subsection QuickRef_SelfadjointMatrix Symmetric/selfadjoint views
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200613
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200614Just as for triangular matrix, you can reference any triangular part of a square matrix to see it as a selfadjoint
615matrix and perform special and optimized operations. Again the opposite triangular part is never referenced and can be
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200616used to store other information.
617
618<table class="tutorial_code">
619<tr><td>
620Conversion to a dense matrix:
621</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200622m2 = m.selfadjointView<Eigen::Lower>();\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200623</td></tr>
624<tr><td>
625Product with another general matrix or vector:
626</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200627m3 = s1 * m1.conjugate().selfadjointView<Eigen::Upper>() * m3;
628m3 -= s1 * m3.adjoint() * m1.selfadjointView<Eigen::Lower>();\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200629</td></tr>
630<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200631Rank 1 and rank K update: \n
632\f$ upper(M_1) += s1 M_2^* M_2 \f$ \n
633\f$ lower(M_1) -= M_2 M_2^* \f$
634</td><td>\n \code
635M1.selfadjointView<Eigen::Upper>().rankUpdate(M2,s1);
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200636m1.selfadjointView<Eigen::Lower>().rankUpdate(m2.adjoint(),-1); \endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200637</td></tr>
638<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200639Rank 2 update: (\f$ M += s u v^* + s v u^* \f$)
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200640</td><td>\code
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200641M.selfadjointView<Eigen::Upper>().rankUpdate(u,v,s);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200642\endcode
643</td></tr>
644<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200645Solving linear equations:\n(\f$ M_2 := M_1^{-1} M_2 \f$)
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200646</td><td>\code
647// via a standard Cholesky factorization
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200648m2 = m1.selfadjointView<Eigen::Upper>().llt().solve(m2);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200649// via a Cholesky factorization with pivoting
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200650m2 = m1.selfadjointView<Eigen::Lower>().ldlt().solve(m2);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200651\endcode
652</td></tr>
653</table>
654
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200655*/
656
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200657/*
658<table class="tutorial_code">
659<tr><td>
660\link MatrixBase::asDiagonal() make a diagonal matrix \endlink \n from a vector </td><td>\code
661mat1 = vec1.asDiagonal();\endcode
662</td></tr>
663<tr><td>
664Declare a diagonal matrix</td><td>\code
665DiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
666diag1.diagonal() = vector;\endcode
667</td></tr>
668<tr><td>Access \link MatrixBase::diagonal() the diagonal and super/sub diagonals of a matrix \endlink as a vector (read/write)</td>
669 <td>\code
670vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal
671vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal
672vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal
673vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal
674vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal
675\endcode</td>
676</tr>
677
678<tr><td>View on a triangular part of a matrix (read/write)</td>
679 <td>\code
680mat2 = mat1.triangularView<Xxx>();
681// Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower
682mat1.triangularView<Upper>() = mat2 + mat3; // only the upper part is evaluated and referenced
683\endcode</td></tr>
684
685<tr><td>View a triangular part as a symmetric/self-adjoint matrix (read/write)</td>
686 <td>\code
687mat2 = mat1.selfadjointView<Xxx>(); // Xxx = Upper or Lower
688mat1.selfadjointView<Upper>() = mat2 + mat2.adjoint(); // evaluated and write to the upper triangular part only
689\endcode</td></tr>
690
691</table>
692
693Optimized products:
694\code
695mat3 += scalar * vec1.asDiagonal() * mat1
696mat3 += scalar * mat1 * vec1.asDiagonal()
697mat3.noalias() += scalar * mat1.triangularView<Xxx>() * mat2
698mat3.noalias() += scalar * mat2 * mat1.triangularView<Xxx>()
699mat3.noalias() += scalar * mat1.selfadjointView<Upper or Lower>() * mat2
700mat3.noalias() += scalar * mat2 * mat1.selfadjointView<Upper or Lower>()
701mat1.selfadjointView<Upper or Lower>().rankUpdate(mat2);
702mat1.selfadjointView<Upper or Lower>().rankUpdate(mat2.adjoint(), scalar);
703\endcode
704
705Inverse products: (all are optimized)
706\code
707mat3 = vec1.asDiagonal().inverse() * mat1
708mat3 = mat1 * diag1.inverse()
709mat1.triangularView<Xxx>().solveInPlace(mat2)
710mat1.triangularView<Xxx>().solveInPlace<OnTheRight>(mat2)
711mat2 = mat1.selfadjointView<Upper or Lower>().llt().solve(mat2)
712\endcode
713
714*/
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200715}