blob: 59dc74fa20a4ac2feb975943cb0f22a20449c0fc [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:
43\code
44typedef Matrix<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyMatrixType;
45typedef Array<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyArrayType;
46\endcode
47
48\li \c Scalar is the scalar type of the coefficients (e.g., \c float, \c double, \c bool, \c int, etc.).
49\li \c RowsAtCompileTime and \c ColsAtCompileTime are the number of rows and columns of the matrix as known at compile-time or \c Dynamic.
50\li \c Options can be \c ColMajor or \c RowMajor, default is \c ColMajor. (see class Matrix for more options)
51
52All 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:
Gael Guennebaud5c866f22010-06-26 16:59:18 +020053\code
54Matrix<double, 6, Dynamic> // Dynamic number of columns (heap allocation)
55Matrix<double, Dynamic, 2> // Dynamic number of rows (heap allocation)
56Matrix<double, Dynamic, Dynamic, RowMajor> // Fully dynamic, row major (heap allocation)
57Matrix<double, 13, 3> // Fully fixed (static allocation)
58\endcode
59
60In most cases, you can simply use one of the convenience typedefs for \ref matrixtypedefs "matrices" and \ref arraytypedefs "arrays". Some examples:
Benoit Jacob9fa54d42010-10-19 08:42:49 -040061<table class="example">
62<tr><th>Matrices</th><th>Arrays</th></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +020063<tr><td>\code
64Matrix<float,Dynamic,Dynamic> <=> MatrixXf
65Matrix<double,Dynamic,1> <=> VectorXd
66Matrix<int,1,Dynamic> <=> RowVectorXi
67Matrix<float,3,3> <=> Matrix3f
68Matrix<float,4,1> <=> Vector4f
69\endcode</td><td>\code
70Array<float,Dynamic,Dynamic> <=> ArrayXXf
71Array<double,Dynamic,1> <=> ArrayXd
72Array<int,1,Dynamic> <=> RowArrayXi
73Array<float,3,3> <=> Array33f
74Array<float,4,1> <=> Array4f
75\endcode</td></tr>
76</table>
77
78Conversion between the matrix and array worlds:
79\code
80Array44f a1, a1;
81Matrix4f m1, m2;
Gael Guennebaudf3c64c72010-06-26 22:19:03 +020082m1 = a1 * a2; // coeffwise product, implicit conversion from array to matrix.
83a1 = m1 * m2; // matrix product, implicit conversion from matrix to array.
84a2 = a1 + m1.array(); // mixing array and matrix is forbidden
85m2 = a1.matrix() + m1; // and explicit conversion is required.
86ArrayWrapper<Matrix4f> m1a(m1); // m1a is an alias for m1.array(), they share the same coefficients
Gael Guennebaud5c866f22010-06-26 16:59:18 +020087MatrixWrapper<Array44f> a1m(a1);
88\endcode
89
90In the rest of this document we will use the following symbols to emphasize the features which are specifics to a given kind of object:
91\li <a name="matrixonly"><a/>\matrixworld linear algebra matrix and vector only
92\li <a name="arrayonly"><a/>\arrayworld array objects only
93
Gael Guennebaud5c866f22010-06-26 16:59:18 +020094\subsection QuickRef_Basics Basic matrix manipulation
95
Gael Guennebaudca4bd582010-10-19 11:59:11 +020096<table class="manual">
97<tr><th></th><th>1D objects</th><th>2D objects</th><th>Notes</th></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +020098<tr><td>Constructors</td>
99<td>\code
100Vector4d v4;
101Vector2f v1(x, y);
102Array3i v2(x, y, z);
103Vector4d v3(x, y, z, w);
104
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200105VectorXf v5; // empty object
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200106ArrayXf v6(size);
107\endcode</td><td>\code
108Matrix4f m1;
109
110
111
112
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200113MatrixXf m5; // empty object
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200114MatrixXf m6(nb_rows, nb_columns);
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200115\endcode</td><td class="note">
116By default, the coefficients \n are left uninitialized</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200117<tr class="alt"><td>Comma initializer</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200118<td>\code
119Vector3f v1; v1 << x, y, z;
120ArrayXf v2(4); v2 << 1, 2, 3, 4;
121
122\endcode</td><td>\code
123Matrix3f m1; m1 << 1, 2, 3,
124 4, 5, 6,
125 7, 8, 9;
126\endcode</td><td></td></tr>
127
128<tr><td>Comma initializer (bis)</td>
129<td colspan="2">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200130\include Tutorial_commainit_02.cpp
131</td>
132<td>
133output:
134\verbinclude Tutorial_commainit_02.out
Hauke Heibele1348b92010-06-30 18:04:36 +0200135</td>
136</tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200137
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200138<tr class="alt"><td>Runtime info</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200139<td>\code
140vector.size();
141
142vector.innerStride();
143vector.data();
144\endcode</td><td>\code
145matrix.rows(); matrix.cols();
146matrix.innerSize(); matrix.outerSize();
147matrix.innerStride(); matrix.outerStride();
148matrix.data();
Hauke Heibelb1741c12010-06-30 15:52:00 +0200149\endcode</td><td class="note">Inner/Outer* are storage order dependent</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200150<tr><td>Compile-time info</td>
151<td colspan="2">\code
152ObjectType::Scalar ObjectType::RowsAtCompileTime
153ObjectType::RealScalar ObjectType::ColsAtCompileTime
154ObjectType::Index ObjectType::SizeAtCompileTime
155\endcode</td><td></td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200156<tr class="alt"><td>Resizing</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200157<td>\code
158vector.resize(size);
159
160
161vector.resizeLike(other_vector);
162vector.conservativeResize(size);
163\endcode</td><td>\code
164matrix.resize(nb_rows, nb_cols);
165matrix.resize(Eigen::NoChange, nb_cols);
166matrix.resize(nb_rows, Eigen::NoChange);
167matrix.resizeLike(other_matrix);
168matrix.conservativeResize(nb_rows, nb_cols);
Hauke Heibelb1741c12010-06-30 15:52:00 +0200169\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 +0200170
171<tr><td>Coeff access with \n range checking</td>
172<td>\code
173vector(i) vector.x()
174vector[i] vector.y()
175 vector.z()
176 vector.w()
177\endcode</td><td>\code
178matrix(i,j)
Jitse Niesen8db9acb2010-12-27 15:07:11 +0000179\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 +0200180
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200181<tr class="alt"><td>Coeff access without \n range checking</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200182<td>\code
183vector.coeff(i)
184vector.coeffRef(i)
185\endcode</td><td>\code
186matrix.coeff(i,j)
187matrix.coeffRef(i,j)
188\endcode</td><td></td></tr>
189
190<tr><td>Assignment/copy</td>
191<td colspan="2">\code
192object = expression;
193object_of_float = expression_of_double.cast<float>();
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200194\endcode</td><td class="note">the destination is automatically resized (if possible)</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200195
196</table>
197
198\subsection QuickRef_PredefMat Predefined Matrices
199
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200200<table class="manual">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200201<tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200202 <th>Fixed-size matrix or vector</th>
203 <th>Dynamic-size matrix</th>
204 <th>Dynamic-size vector</th>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200205</tr>
206<tr style="border-bottom-style: none;">
207 <td>
208\code
209typedef {Matrix3f|Array33f} FixedXD;
210FixedXD x;
211
212x = FixedXD::Zero();
213x = FixedXD::Ones();
214x = FixedXD::Constant(value);
215x = FixedXD::Random();
Trevor Ironse112ad82010-12-29 10:08:41 -0700216x = FixedXD::LinSpaced(size, low, high);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200217
218x.setZero();
219x.setOnes();
220x.setConstant(value);
221x.setRandom();
Trevor Ironse112ad82010-12-29 10:08:41 -0700222x.setLinSpaced(size, low, high);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200223\endcode
224 </td>
225 <td>
226\code
227typedef {MatrixXf|ArrayXXf} Dynamic2D;
228Dynamic2D x;
229
230x = Dynamic2D::Zero(rows, cols);
231x = Dynamic2D::Ones(rows, cols);
232x = Dynamic2D::Constant(rows, cols, value);
233x = Dynamic2D::Random(rows, cols);
Jitse Niesen3abbdfd2010-07-20 21:55:22 +0100234N/A
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200235
236x.setZero(rows, cols);
237x.setOnes(rows, cols);
238x.setConstant(rows, cols, value);
239x.setRandom(rows, cols);
Jitse Niesen3abbdfd2010-07-20 21:55:22 +0100240N/A
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200241\endcode
242 </td>
243 <td>
244\code
245typedef {VectorXf|ArrayXf} Dynamic1D;
246Dynamic1D x;
247
248x = Dynamic1D::Zero(size);
249x = Dynamic1D::Ones(size);
250x = Dynamic1D::Constant(size, value);
251x = Dynamic1D::Random(size);
Trevor Ironse112ad82010-12-29 10:08:41 -0700252x = Dynamic1D::LinSpaced(size, low, high);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200253
254x.setZero(size);
255x.setOnes(size);
256x.setConstant(size, value);
257x.setRandom(size);
Trevor Ironse112ad82010-12-29 10:08:41 -0700258x.setLinSpaced(size, low, high);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200259\endcode
260 </td>
261</tr>
262
263<tr><td colspan="3">Identity and \link MatrixBase::Unit basis vectors \endlink \matrixworld</td></tr>
264<tr style="border-bottom-style: none;">
265 <td>
266\code
267x = FixedXD::Identity();
268x.setIdentity();
269
270Vector3f::UnitX() // 1 0 0
271Vector3f::UnitY() // 0 1 0
272Vector3f::UnitZ() // 0 0 1
273\endcode
274 </td>
275 <td>
276\code
277x = Dynamic2D::Identity(rows, cols);
278x.setIdentity(rows, cols);
279
280
281
282N/A
283\endcode
284 </td>
285 <td>\code
286N/A
287
288
289VectorXf::Unit(size,i)
290VectorXf::Unit(4,1) == Vector4f(0,1,0,0)
291 == Vector4f::UnitY()
292\endcode
293 </td>
294</tr>
295</table>
296
297
298
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200299\subsection QuickRef_Map Mapping external arrays
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200300
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200301<table class="manual">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200302<tr>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200303<td>Contiguous \n memory</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200304<td>\code
305float data[] = {1,2,3,4};
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200306Map<Vector3f> v1(data); // uses v1 as a Vector3f object
307Map<ArrayXf> v2(data,3); // uses v2 as a ArrayXf object
308Map<Array22f> m1(data); // uses m1 as a Array22f object
309Map<MatrixXf> m2(data,2,2); // uses m2 as a MatrixXf object
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200310\endcode</td>
311</tr>
312<tr>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200313<td>Typical usage \n of strides</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200314<td>\code
315float data[] = {1,2,3,4,5,6,7,8,9};
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200316Map<VectorXf,0,InnerStride<2> > v1(data,3); // = [1,3,5]
317Map<VectorXf,0,InnerStride<> > v2(data,3,InnerStride<>(3)); // = [1,4,7]
318Map<MatrixXf,0,OuterStride<3> > m2(data,2,3); // both lines |1,4,7|
319Map<MatrixXf,0,OuterStride<> > m1(data,2,3,OuterStride<>(3)); // are equal to: |2,5,8|
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200320\endcode</td>
321</tr>
322</table>
323
324
325<a href="#" class="top">top</a>
326\section QuickRef_ArithmeticOperators Arithmetic Operators
327
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200328<table class="manual">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200329<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200330add \n subtract</td><td>\code
331mat3 = mat1 + mat2; mat3 += mat1;
332mat3 = mat1 - mat2; mat3 -= mat1;\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200333</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200334<tr class="alt"><td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200335scalar product</td><td>\code
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200336mat3 = mat1 * s1; mat3 *= s1; mat3 = s1 * mat1;
337mat3 = mat1 / s1; mat3 /= s1;\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200338</td></tr>
339<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200340matrix/vector \n products \matrixworld</td><td>\code
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200341col2 = mat1 * col1;
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200342row2 = row1 * mat1; row1 *= mat1;
343mat3 = mat1 * mat2; mat3 *= mat1; \endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200344</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200345<tr class="alt"><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200346transposition \n adjoint \matrixworld</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200347mat1 = mat2.transpose(); mat1.transposeInPlace();
348mat1 = mat2.adjoint(); mat1.adjointInPlace();
349\endcode
350</td></tr>
351<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200352\link MatrixBase::dot() dot \endlink product \n inner product \matrixworld</td><td>\code
353scalar = vec1.dot(vec2);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200354scalar = col1.adjoint() * col2;
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200355scalar = (col1.adjoint() * col2).value();\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200356</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200357<tr class="alt"><td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200358outer product \matrixworld</td><td>\code
359mat = col1 * col2.transpose();\endcode
360</td></tr>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200361
362<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200363\link MatrixBase::norm() norm \endlink \n \link MatrixBase::normalized() normalization \endlink \matrixworld</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200364scalar = vec1.norm(); scalar = vec1.squaredNorm()
365vec2 = vec1.normalized(); vec1.normalize(); // inplace \endcode
366</td></tr>
367
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200368<tr class="alt"><td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200369\link MatrixBase::cross() cross product \endlink \matrixworld</td><td>\code
370#include <Eigen/Geometry>
371vec3 = vec1.cross(vec2);\endcode</td></tr>
372</table>
373
374<a href="#" class="top">top</a>
375\section QuickRef_Coeffwise Coefficient-wise \& Array operators
376Coefficient-wise operators for matrices and vectors:
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200377<table class="manual">
378<tr><th>Matrix API \matrixworld</th><th>Via Array conversions</th></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200379<tr><td>\code
380mat1.cwiseMin(mat2)
381mat1.cwiseMax(mat2)
Gael Guennebauddbefd7a2010-06-28 13:30:10 +0200382mat1.cwiseAbs2()
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200383mat1.cwiseAbs()
384mat1.cwiseSqrt()
385mat1.cwiseProduct(mat2)
386mat1.cwiseQuotient(mat2)\endcode
387</td><td>\code
388mat1.array().min(mat2.array())
389mat1.array().max(mat2.array())
390mat1.array().abs2()
391mat1.array().abs()
392mat1.array().sqrt()
393mat1.array() * mat2.array()
394mat1.array() / mat2.array()
395\endcode</td></tr>
396</table>
397
Gael Guennebaudeda59ff2011-02-17 18:07:21 +0100398It is also very simple to apply any user defined function \c foo using DenseBase::unaryExpr together with std::ptr_fun:
399\code mat1.unaryExpr(std::ptr_fun(foo))\endcode
400
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200401Array operators:\arrayworld
402
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200403<table class="manual">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200404<tr><td>Arithmetic operators</td><td>\code
405array1 * array2 array1 / array2 array1 *= array2 array1 /= array2
406array1 + scalar array1 - scalar array1 += scalar array1 -= scalar
407\endcode</td></tr>
408<tr><td>Comparisons</td><td>\code
409array1 < array2 array1 > array2 array1 < scalar array1 > scalar
410array1 <= array2 array1 >= array2 array1 <= scalar array1 >= scalar
411array1 == array2 array1 != array2 array1 == scalar array1 != scalar
412\endcode</td></tr>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200413<tr><td>Trigo, power, and \n misc functions \n and the STL variants</td><td>\code
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200414array1.min(array2) std::min(array1,array2)
415array1.max(array2) std::max(array1,array2)
416array1.abs2()
417array1.abs() std::abs(array1)
418array1.sqrt() std::sqrt(array1)
419array1.log() std::log(array1)
420array1.exp() std::exp(array1)
421array1.pow(exponent) std::pow(array1,exponent)
422array1.square()
423array1.cube()
424array1.inverse()
425array1.sin() std::sin(array1)
426array1.cos() std::cos(array1)
427array1.tan() std::tan(array1)
Gael Guennebaudeda59ff2011-02-17 18:07:21 +0100428array1.asin() std::asin(array1)
429array1.acos() std::acos(array1)
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200430\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:
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200461\code
462int i, j;
463s = vector.minCoeff(&i); // s == vector[i]
464s = matrix.maxCoeff(&i, &j); // s == matrix(i,j)
465\endcode
466Typical use cases of all() and any():
467\code
468if((array1 > 0).all()) ... // if all coefficients of array1 are greater than 0 ...
469if((array1 < array2).any()) ... // if there exist a pair i,j such that array1(i,j) < array2(i,j) ...
470\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200471
472
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200473<a href="#" class="top">top</a>\section QuickRef_Blocks Sub-matrices
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200474
Jitse Niesen1420f8b2010-07-25 20:29:07 +0100475Read-write access to a \link DenseBase::col(Index) column \endlink
476or a \link DenseBase::row(Index) row \endlink of a matrix (or array):
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200477\code
478mat1.row(i) = mat2.col(j);
479mat1.col(j1).swap(mat1.col(j2));
480\endcode
481
482Read-write access to sub-vectors:
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200483<table class="manual">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200484<tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200485<th>Default versions</th>
486<th>Optimized versions when the size \n is known at compile time</th></tr>
487<th></th>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200488
489<tr><td>\code vec1.head(n)\endcode</td><td>\code vec1.head<n>()\endcode</td><td>the first \c n coeffs </td></tr>
490<tr><td>\code vec1.tail(n)\endcode</td><td>\code vec1.tail<n>()\endcode</td><td>the last \c n coeffs </td></tr>
491<tr><td>\code vec1.segment(pos,n)\endcode</td><td>\code vec1.segment<n>(pos)\endcode</td>
Jitse Niesen49747fa2010-07-06 13:10:08 +0100492 <td>the \c n coeffs in \n the range [\c pos : \c pos + \c n [</td></tr>
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200493<tr class="alt"><td colspan="3">
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200494
Gael Guennebaudca4bd582010-10-19 11:59:11 +0200495Read-write access to sub-matrices:</td></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200496<tr>
497 <td>\code mat1.block(i,j,rows,cols)\endcode
Jitse Niesen1420f8b2010-07-25 20:29:07 +0100498 \link DenseBase::block(Index,Index,Index,Index) (more) \endlink</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200499 <td>\code mat1.block<rows,cols>(i,j)\endcode
Jitse Niesen1420f8b2010-07-25 20:29:07 +0100500 \link DenseBase::block(Index,Index) (more) \endlink</td>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200501 <td>the \c rows x \c cols sub-matrix \n starting from position (\c i,\c j)</td></tr>
502<tr><td>\code
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200503 mat1.topLeftCorner(rows,cols)
504 mat1.topRightCorner(rows,cols)
505 mat1.bottomLeftCorner(rows,cols)
506 mat1.bottomRightCorner(rows,cols)\endcode
507 <td>\code
508 mat1.topLeftCorner<rows,cols>()
509 mat1.topRightCorner<rows,cols>()
510 mat1.bottomLeftCorner<rows,cols>()
511 mat1.bottomRightCorner<rows,cols>()\endcode
512 <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 +0200513 <tr><td>\code
514 mat1.topRows(rows)
515 mat1.bottomRows(rows)
516 mat1.leftCols(cols)
517 mat1.rightCols(cols)\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200518 <td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200519 mat1.topRows<rows>()
520 mat1.bottomRows<rows>()
521 mat1.leftCols<cols>()
522 mat1.rightCols<cols>()\endcode
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200523 <td>specialized versions of block() \n when the block fit two corners</td></tr>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200524</table>
525
526
527
528
529<a href="#" class="top">top</a>\section QuickRef_DiagTriSymm Diagonal, Triangular, and Self-adjoint matrices
530(matrix world \matrixworld)
531
532\subsection QuickRef_Diagonal Diagonal matrices
533
Benoit Jacob9fa54d42010-10-19 08:42:49 -0400534<table class="example">
535<tr><th>Operation</th><th>Code</th></tr>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200536<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200537view a vector \link MatrixBase::asDiagonal() as a diagonal matrix \endlink \n </td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200538mat1 = vec1.asDiagonal();\endcode
539</td></tr>
540<tr><td>
541Declare a diagonal matrix</td><td>\code
542DiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
543diag1.diagonal() = vector;\endcode
544</td></tr>
Jitse Niesen1420f8b2010-07-25 20:29:07 +0100545<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 +0200546 <td>\code
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200547vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal
548vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal
549vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal
550vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal
551vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200552\endcode</td>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200553</tr>
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200554
555<tr><td>Optimized products and inverse</td>
556 <td>\code
557mat3 = scalar * diag1 * mat1;
558mat3 += scalar * mat1 * vec1.asDiagonal();
559mat3 = vec1.asDiagonal().inverse() * mat1
560mat3 = mat1 * diag1.inverse()
561\endcode</td>
562</tr>
563
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200564</table>
565
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200566\subsection QuickRef_TriangularView Triangular views
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200567
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200568TriangularView 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 +0200569
Benoit Jacob9fa54d42010-10-19 08:42:49 -0400570<table class="example">
571<tr><th>Operation</th><th>Code</th></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200572<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200573Reference to a triangular with optional \n
574unit or null diagonal (read/write):
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200575</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200576m.triangularView<Xxx>()
577\endcode \n
578\c Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200579</td></tr>
580<tr><td>
581Writing to a specific triangular part:\n (only the referenced triangular part is evaluated)
582</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200583m1.triangularView<Eigen::Lower>() = m2 + m3 \endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200584</td></tr>
585<tr><td>
586Conversion to a dense matrix setting the opposite triangular part to zero:
587</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200588m2 = m1.triangularView<Eigen::UnitUpper>()\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200589</td></tr>
590<tr><td>
591Products:
592</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200593m3 += s1 * m1.adjoint().triangularView<Eigen::UnitUpper>() * m2
594m3 -= s1 * m2.conjugate() * m1.adjoint().triangularView<Eigen::Lower>() \endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200595</td></tr>
596<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200597Solving linear equations:\n
598\f$ M_2 := L_1^{-1} M_2 \f$ \n
599\f$ M_3 := {L_1^*}^{-1} M_3 \f$ \n
Thomas Capricelli8e21cef2010-07-22 13:15:15 +0200600\f$ M_4 := M_4 U_1^{-1} \f$
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200601</td><td>\n \code
602L1.triangularView<Eigen::UnitLower>().solveInPlace(M2)
603L1.triangularView<Eigen::Lower>().adjoint().solveInPlace(M3)
604U1.triangularView<Eigen::Upper>().solveInPlace<OnTheRight>(M4)\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200605</td></tr>
606</table>
607
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200608\subsection QuickRef_SelfadjointMatrix Symmetric/selfadjoint views
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200609
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200610Just as for triangular matrix, you can reference any triangular part of a square matrix to see it as a selfadjoint
611matrix and perform special and optimized operations. Again the opposite triangular part is never referenced and can be
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200612used to store other information.
613
Benoit Jacob9fa54d42010-10-19 08:42:49 -0400614<table class="example">
615<tr><th>Operation</th><th>Code</th></tr>
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200616<tr><td>
617Conversion to a dense matrix:
618</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200619m2 = m.selfadjointView<Eigen::Lower>();\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200620</td></tr>
621<tr><td>
622Product with another general matrix or vector:
623</td><td>\code
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200624m3 = s1 * m1.conjugate().selfadjointView<Eigen::Upper>() * m3;
625m3 -= s1 * m3.adjoint() * m1.selfadjointView<Eigen::Lower>();\endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200626</td></tr>
627<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200628Rank 1 and rank K update: \n
629\f$ upper(M_1) += s1 M_2^* M_2 \f$ \n
630\f$ lower(M_1) -= M_2 M_2^* \f$
631</td><td>\n \code
632M1.selfadjointView<Eigen::Upper>().rankUpdate(M2,s1);
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200633m1.selfadjointView<Eigen::Lower>().rankUpdate(m2.adjoint(),-1); \endcode
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200634</td></tr>
635<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200636Rank 2 update: (\f$ M += s u v^* + s v u^* \f$)
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200637</td><td>\code
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200638M.selfadjointView<Eigen::Upper>().rankUpdate(u,v,s);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200639\endcode
640</td></tr>
641<tr><td>
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200642Solving linear equations:\n(\f$ M_2 := M_1^{-1} M_2 \f$)
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200643</td><td>\code
644// via a standard Cholesky factorization
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200645m2 = m1.selfadjointView<Eigen::Upper>().llt().solve(m2);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200646// via a Cholesky factorization with pivoting
Gael Guennebaudf3c64c72010-06-26 22:19:03 +0200647m2 = m1.selfadjointView<Eigen::Lower>().ldlt().solve(m2);
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200648\endcode
649</td></tr>
650</table>
651
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200652*/
653
Gael Guennebaud1c783e22010-06-26 18:49:50 +0200654/*
655<table class="tutorial_code">
656<tr><td>
657\link MatrixBase::asDiagonal() make a diagonal matrix \endlink \n from a vector </td><td>\code
658mat1 = vec1.asDiagonal();\endcode
659</td></tr>
660<tr><td>
661Declare a diagonal matrix</td><td>\code
662DiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
663diag1.diagonal() = vector;\endcode
664</td></tr>
665<tr><td>Access \link MatrixBase::diagonal() the diagonal and super/sub diagonals of a matrix \endlink as a vector (read/write)</td>
666 <td>\code
667vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal
668vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal
669vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal
670vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal
671vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal
672\endcode</td>
673</tr>
674
675<tr><td>View on a triangular part of a matrix (read/write)</td>
676 <td>\code
677mat2 = mat1.triangularView<Xxx>();
678// Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower
679mat1.triangularView<Upper>() = mat2 + mat3; // only the upper part is evaluated and referenced
680\endcode</td></tr>
681
682<tr><td>View a triangular part as a symmetric/self-adjoint matrix (read/write)</td>
683 <td>\code
684mat2 = mat1.selfadjointView<Xxx>(); // Xxx = Upper or Lower
685mat1.selfadjointView<Upper>() = mat2 + mat2.adjoint(); // evaluated and write to the upper triangular part only
686\endcode</td></tr>
687
688</table>
689
690Optimized products:
691\code
692mat3 += scalar * vec1.asDiagonal() * mat1
693mat3 += scalar * mat1 * vec1.asDiagonal()
694mat3.noalias() += scalar * mat1.triangularView<Xxx>() * mat2
695mat3.noalias() += scalar * mat2 * mat1.triangularView<Xxx>()
696mat3.noalias() += scalar * mat1.selfadjointView<Upper or Lower>() * mat2
697mat3.noalias() += scalar * mat2 * mat1.selfadjointView<Upper or Lower>()
698mat1.selfadjointView<Upper or Lower>().rankUpdate(mat2);
699mat1.selfadjointView<Upper or Lower>().rankUpdate(mat2.adjoint(), scalar);
700\endcode
701
702Inverse products: (all are optimized)
703\code
704mat3 = vec1.asDiagonal().inverse() * mat1
705mat3 = mat1 * diag1.inverse()
706mat1.triangularView<Xxx>().solveInPlace(mat2)
707mat1.triangularView<Xxx>().solveInPlace<OnTheRight>(mat2)
708mat2 = mat1.selfadjointView<Upper or Lower>().llt().solve(mat2)
709\endcode
710
711*/
Gael Guennebaud5c866f22010-06-26 16:59:18 +0200712}