xref: /aosp_15_r20/external/eigen/doc/QuickReference.dox (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1*bf2c3715SXin Linamespace Eigen {
2*bf2c3715SXin Li
3*bf2c3715SXin Li/** \eigenManualPage QuickRefPage Quick reference guide
4*bf2c3715SXin Li
5*bf2c3715SXin Li\eigenAutoToc
6*bf2c3715SXin Li
7*bf2c3715SXin Li<hr>
8*bf2c3715SXin Li
9*bf2c3715SXin Li<a href="#" class="top">top</a>
10*bf2c3715SXin Li\section QuickRef_Headers Modules and Header files
11*bf2c3715SXin Li
12*bf2c3715SXin LiThe 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.
13*bf2c3715SXin Li
14*bf2c3715SXin Li<table class="manual">
15*bf2c3715SXin Li<tr><th>Module</th><th>Header file</th><th>Contents</th></tr>
16*bf2c3715SXin Li<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>
17*bf2c3715SXin Li<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>
18*bf2c3715SXin Li<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>
19*bf2c3715SXin Li<tr class="alt"><td>\link Cholesky_Module Cholesky \endlink</td><td>\code#include <Eigen/Cholesky>\endcode</td><td>LLT and LDLT Cholesky factorization with solver</td></tr>
20*bf2c3715SXin Li<tr            ><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>
21*bf2c3715SXin Li<tr class="alt"><td>\link SVD_Module SVD \endlink</td><td>\code#include <Eigen/SVD>\endcode</td><td>SVD decompositions with least-squares solver (JacobiSVD, BDCSVD)</td></tr>
22*bf2c3715SXin Li<tr            ><td>\link QR_Module QR \endlink</td><td>\code#include <Eigen/QR>\endcode</td><td>QR decomposition with solver (HouseholderQR, ColPivHouseholderQR, FullPivHouseholderQR)</td></tr>
23*bf2c3715SXin Li<tr class="alt"><td>\link Eigenvalues_Module Eigenvalues \endlink</td><td>\code#include <Eigen/Eigenvalues>\endcode</td><td>Eigenvalue, eigenvector decompositions (EigenSolver, SelfAdjointEigenSolver, ComplexEigenSolver)</td></tr>
24*bf2c3715SXin Li<tr            ><td>\link Sparse_Module Sparse \endlink</td><td>\code#include <Eigen/Sparse>\endcode</td><td>%Sparse matrix storage and related basic linear algebra (SparseMatrix, SparseVector) \n (see \ref SparseQuickRefPage for details on sparse modules)</td></tr>
25*bf2c3715SXin Li<tr class="alt"><td></td><td>\code#include <Eigen/Dense>\endcode</td><td>Includes Core, Geometry, LU, Cholesky, SVD, QR, and Eigenvalues header files</td></tr>
26*bf2c3715SXin Li<tr            ><td></td><td>\code#include <Eigen/Eigen>\endcode</td><td>Includes %Dense and %Sparse header files (the whole Eigen library)</td></tr>
27*bf2c3715SXin Li</table>
28*bf2c3715SXin Li
29*bf2c3715SXin Li<a href="#" class="top">top</a>
30*bf2c3715SXin Li\section QuickRef_Types Array, matrix and vector types
31*bf2c3715SXin Li
32*bf2c3715SXin Li
33*bf2c3715SXin Li\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:
34*bf2c3715SXin Li\code
35*bf2c3715SXin Litypedef Matrix<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyMatrixType;
36*bf2c3715SXin Litypedef Array<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyArrayType;
37*bf2c3715SXin Li\endcode
38*bf2c3715SXin Li
39*bf2c3715SXin Li\li \c Scalar is the scalar type of the coefficients (e.g., \c float, \c double, \c bool, \c int, etc.).
40*bf2c3715SXin Li\li \c RowsAtCompileTime and \c ColsAtCompileTime are the number of rows and columns of the matrix as known at compile-time or \c Dynamic.
41*bf2c3715SXin Li\li \c Options can be \c ColMajor or \c RowMajor, default is \c ColMajor. (see class Matrix for more options)
42*bf2c3715SXin Li
43*bf2c3715SXin LiAll 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:
44*bf2c3715SXin Li\code
45*bf2c3715SXin LiMatrix<double, 6, Dynamic>                  // Dynamic number of columns (heap allocation)
46*bf2c3715SXin LiMatrix<double, Dynamic, 2>                  // Dynamic number of rows (heap allocation)
47*bf2c3715SXin LiMatrix<double, Dynamic, Dynamic, RowMajor>  // Fully dynamic, row major (heap allocation)
48*bf2c3715SXin LiMatrix<double, 13, 3>                       // Fully fixed (usually allocated on stack)
49*bf2c3715SXin Li\endcode
50*bf2c3715SXin Li
51*bf2c3715SXin LiIn most cases, you can simply use one of the convenience typedefs for \ref matrixtypedefs "matrices" and \ref arraytypedefs "arrays". Some examples:
52*bf2c3715SXin Li<table class="example">
53*bf2c3715SXin Li<tr><th>Matrices</th><th>Arrays</th></tr>
54*bf2c3715SXin Li<tr><td>\code
55*bf2c3715SXin LiMatrix<float,Dynamic,Dynamic>   <=>   MatrixXf
56*bf2c3715SXin LiMatrix<double,Dynamic,1>        <=>   VectorXd
57*bf2c3715SXin LiMatrix<int,1,Dynamic>           <=>   RowVectorXi
58*bf2c3715SXin LiMatrix<float,3,3>               <=>   Matrix3f
59*bf2c3715SXin LiMatrix<float,4,1>               <=>   Vector4f
60*bf2c3715SXin Li\endcode</td><td>\code
61*bf2c3715SXin LiArray<float,Dynamic,Dynamic>    <=>   ArrayXXf
62*bf2c3715SXin LiArray<double,Dynamic,1>         <=>   ArrayXd
63*bf2c3715SXin LiArray<int,1,Dynamic>            <=>   RowArrayXi
64*bf2c3715SXin LiArray<float,3,3>                <=>   Array33f
65*bf2c3715SXin LiArray<float,4,1>                <=>   Array4f
66*bf2c3715SXin Li\endcode</td></tr>
67*bf2c3715SXin Li</table>
68*bf2c3715SXin Li
69*bf2c3715SXin LiConversion between the matrix and array worlds:
70*bf2c3715SXin Li\code
71*bf2c3715SXin LiArray44f a1, a2;
72*bf2c3715SXin LiMatrix4f m1, m2;
73*bf2c3715SXin Lim1 = a1 * a2;                     // coeffwise product, implicit conversion from array to matrix.
74*bf2c3715SXin Lia1 = m1 * m2;                     // matrix product, implicit conversion from matrix to array.
75*bf2c3715SXin Lia2 = a1 + m1.array();             // mixing array and matrix is forbidden
76*bf2c3715SXin Lim2 = a1.matrix() + m1;            // and explicit conversion is required.
77*bf2c3715SXin LiArrayWrapper<Matrix4f> m1a(m1);   // m1a is an alias for m1.array(), they share the same coefficients
78*bf2c3715SXin LiMatrixWrapper<Array44f> a1m(a1);
79*bf2c3715SXin Li\endcode
80*bf2c3715SXin Li
81*bf2c3715SXin LiIn the rest of this document we will use the following symbols to emphasize the features which are specifics to a given kind of object:
82*bf2c3715SXin Li\li <a name="matrixonly"></a>\matrixworld linear algebra matrix and vector only
83*bf2c3715SXin Li\li <a name="arrayonly"></a>\arrayworld array objects only
84*bf2c3715SXin Li
85*bf2c3715SXin Li\subsection QuickRef_Basics Basic matrix manipulation
86*bf2c3715SXin Li
87*bf2c3715SXin Li<table class="manual">
88*bf2c3715SXin Li<tr><th></th><th>1D objects</th><th>2D objects</th><th>Notes</th></tr>
89*bf2c3715SXin Li<tr><td>Constructors</td>
90*bf2c3715SXin Li<td>\code
91*bf2c3715SXin LiVector4d  v4;
92*bf2c3715SXin LiVector2f  v1(x, y);
93*bf2c3715SXin LiArray3i   v2(x, y, z);
94*bf2c3715SXin LiVector4d  v3(x, y, z, w);
95*bf2c3715SXin Li
96*bf2c3715SXin LiVectorXf  v5; // empty object
97*bf2c3715SXin LiArrayXf   v6(size);
98*bf2c3715SXin Li\endcode</td><td>\code
99*bf2c3715SXin LiMatrix4f  m1;
100*bf2c3715SXin Li
101*bf2c3715SXin Li
102*bf2c3715SXin Li
103*bf2c3715SXin Li
104*bf2c3715SXin LiMatrixXf  m5; // empty object
105*bf2c3715SXin LiMatrixXf  m6(nb_rows, nb_columns);
106*bf2c3715SXin Li\endcode</td><td class="note">
107*bf2c3715SXin LiBy default, the coefficients \n are left uninitialized</td></tr>
108*bf2c3715SXin Li<tr class="alt"><td>Comma initializer</td>
109*bf2c3715SXin Li<td>\code
110*bf2c3715SXin LiVector3f  v1;     v1 << x, y, z;
111*bf2c3715SXin LiArrayXf   v2(4);  v2 << 1, 2, 3, 4;
112*bf2c3715SXin Li
113*bf2c3715SXin Li\endcode</td><td>\code
114*bf2c3715SXin LiMatrix3f  m1;   m1 << 1, 2, 3,
115*bf2c3715SXin Li                      4, 5, 6,
116*bf2c3715SXin Li                      7, 8, 9;
117*bf2c3715SXin Li\endcode</td><td></td></tr>
118*bf2c3715SXin Li
119*bf2c3715SXin Li<tr><td>Comma initializer (bis)</td>
120*bf2c3715SXin Li<td colspan="2">
121*bf2c3715SXin Li\include Tutorial_commainit_02.cpp
122*bf2c3715SXin Li</td>
123*bf2c3715SXin Li<td>
124*bf2c3715SXin Lioutput:
125*bf2c3715SXin Li\verbinclude Tutorial_commainit_02.out
126*bf2c3715SXin Li</td>
127*bf2c3715SXin Li</tr>
128*bf2c3715SXin Li
129*bf2c3715SXin Li<tr class="alt"><td>Runtime info</td>
130*bf2c3715SXin Li<td>\code
131*bf2c3715SXin Livector.size();
132*bf2c3715SXin Li
133*bf2c3715SXin Livector.innerStride();
134*bf2c3715SXin Livector.data();
135*bf2c3715SXin Li\endcode</td><td>\code
136*bf2c3715SXin Limatrix.rows();          matrix.cols();
137*bf2c3715SXin Limatrix.innerSize();     matrix.outerSize();
138*bf2c3715SXin Limatrix.innerStride();   matrix.outerStride();
139*bf2c3715SXin Limatrix.data();
140*bf2c3715SXin Li\endcode</td><td class="note">Inner/Outer* are storage order dependent</td></tr>
141*bf2c3715SXin Li<tr><td>Compile-time info</td>
142*bf2c3715SXin Li<td colspan="2">\code
143*bf2c3715SXin LiObjectType::Scalar              ObjectType::RowsAtCompileTime
144*bf2c3715SXin LiObjectType::RealScalar          ObjectType::ColsAtCompileTime
145*bf2c3715SXin LiObjectType::Index               ObjectType::SizeAtCompileTime
146*bf2c3715SXin Li\endcode</td><td></td></tr>
147*bf2c3715SXin Li<tr class="alt"><td>Resizing</td>
148*bf2c3715SXin Li<td>\code
149*bf2c3715SXin Livector.resize(size);
150*bf2c3715SXin Li
151*bf2c3715SXin Li
152*bf2c3715SXin Livector.resizeLike(other_vector);
153*bf2c3715SXin Livector.conservativeResize(size);
154*bf2c3715SXin Li\endcode</td><td>\code
155*bf2c3715SXin Limatrix.resize(nb_rows, nb_cols);
156*bf2c3715SXin Limatrix.resize(Eigen::NoChange, nb_cols);
157*bf2c3715SXin Limatrix.resize(nb_rows, Eigen::NoChange);
158*bf2c3715SXin Limatrix.resizeLike(other_matrix);
159*bf2c3715SXin Limatrix.conservativeResize(nb_rows, nb_cols);
160*bf2c3715SXin Li\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>
161*bf2c3715SXin Li
162*bf2c3715SXin Li<tr><td>Coeff access with \n range checking</td>
163*bf2c3715SXin Li<td>\code
164*bf2c3715SXin Livector(i)     vector.x()
165*bf2c3715SXin Livector[i]     vector.y()
166*bf2c3715SXin Li              vector.z()
167*bf2c3715SXin Li              vector.w()
168*bf2c3715SXin Li\endcode</td><td>\code
169*bf2c3715SXin Limatrix(i,j)
170*bf2c3715SXin Li\endcode</td><td class="note">Range checking is disabled if \n NDEBUG or EIGEN_NO_DEBUG is defined</td></tr>
171*bf2c3715SXin Li
172*bf2c3715SXin Li<tr class="alt"><td>Coeff access without \n range checking</td>
173*bf2c3715SXin Li<td>\code
174*bf2c3715SXin Livector.coeff(i)
175*bf2c3715SXin Livector.coeffRef(i)
176*bf2c3715SXin Li\endcode</td><td>\code
177*bf2c3715SXin Limatrix.coeff(i,j)
178*bf2c3715SXin Limatrix.coeffRef(i,j)
179*bf2c3715SXin Li\endcode</td><td></td></tr>
180*bf2c3715SXin Li
181*bf2c3715SXin Li<tr><td>Assignment/copy</td>
182*bf2c3715SXin Li<td colspan="2">\code
183*bf2c3715SXin Liobject = expression;
184*bf2c3715SXin Liobject_of_float = expression_of_double.cast<float>();
185*bf2c3715SXin Li\endcode</td><td class="note">the destination is automatically resized (if possible)</td></tr>
186*bf2c3715SXin Li
187*bf2c3715SXin Li</table>
188*bf2c3715SXin Li
189*bf2c3715SXin Li\subsection QuickRef_PredefMat Predefined Matrices
190*bf2c3715SXin Li
191*bf2c3715SXin Li<table class="manual">
192*bf2c3715SXin Li<tr>
193*bf2c3715SXin Li  <th>Fixed-size matrix or vector</th>
194*bf2c3715SXin Li  <th>Dynamic-size matrix</th>
195*bf2c3715SXin Li  <th>Dynamic-size vector</th>
196*bf2c3715SXin Li</tr>
197*bf2c3715SXin Li<tr style="border-bottom-style: none;">
198*bf2c3715SXin Li  <td>
199*bf2c3715SXin Li\code
200*bf2c3715SXin Litypedef {Matrix3f|Array33f} FixedXD;
201*bf2c3715SXin LiFixedXD x;
202*bf2c3715SXin Li
203*bf2c3715SXin Lix = FixedXD::Zero();
204*bf2c3715SXin Lix = FixedXD::Ones();
205*bf2c3715SXin Lix = FixedXD::Constant(value);
206*bf2c3715SXin Lix = FixedXD::Random();
207*bf2c3715SXin Lix = FixedXD::LinSpaced(size, low, high);
208*bf2c3715SXin Li
209*bf2c3715SXin Lix.setZero();
210*bf2c3715SXin Lix.setOnes();
211*bf2c3715SXin Lix.setConstant(value);
212*bf2c3715SXin Lix.setRandom();
213*bf2c3715SXin Lix.setLinSpaced(size, low, high);
214*bf2c3715SXin Li\endcode
215*bf2c3715SXin Li  </td>
216*bf2c3715SXin Li  <td>
217*bf2c3715SXin Li\code
218*bf2c3715SXin Litypedef {MatrixXf|ArrayXXf} Dynamic2D;
219*bf2c3715SXin LiDynamic2D x;
220*bf2c3715SXin Li
221*bf2c3715SXin Lix = Dynamic2D::Zero(rows, cols);
222*bf2c3715SXin Lix = Dynamic2D::Ones(rows, cols);
223*bf2c3715SXin Lix = Dynamic2D::Constant(rows, cols, value);
224*bf2c3715SXin Lix = Dynamic2D::Random(rows, cols);
225*bf2c3715SXin LiN/A
226*bf2c3715SXin Li
227*bf2c3715SXin Lix.setZero(rows, cols);
228*bf2c3715SXin Lix.setOnes(rows, cols);
229*bf2c3715SXin Lix.setConstant(rows, cols, value);
230*bf2c3715SXin Lix.setRandom(rows, cols);
231*bf2c3715SXin LiN/A
232*bf2c3715SXin Li\endcode
233*bf2c3715SXin Li  </td>
234*bf2c3715SXin Li  <td>
235*bf2c3715SXin Li\code
236*bf2c3715SXin Litypedef {VectorXf|ArrayXf} Dynamic1D;
237*bf2c3715SXin LiDynamic1D x;
238*bf2c3715SXin Li
239*bf2c3715SXin Lix = Dynamic1D::Zero(size);
240*bf2c3715SXin Lix = Dynamic1D::Ones(size);
241*bf2c3715SXin Lix = Dynamic1D::Constant(size, value);
242*bf2c3715SXin Lix = Dynamic1D::Random(size);
243*bf2c3715SXin Lix = Dynamic1D::LinSpaced(size, low, high);
244*bf2c3715SXin Li
245*bf2c3715SXin Lix.setZero(size);
246*bf2c3715SXin Lix.setOnes(size);
247*bf2c3715SXin Lix.setConstant(size, value);
248*bf2c3715SXin Lix.setRandom(size);
249*bf2c3715SXin Lix.setLinSpaced(size, low, high);
250*bf2c3715SXin Li\endcode
251*bf2c3715SXin Li  </td>
252*bf2c3715SXin Li</tr>
253*bf2c3715SXin Li
254*bf2c3715SXin Li<tr><td colspan="3">Identity and \link MatrixBase::Unit basis vectors \endlink \matrixworld</td></tr>
255*bf2c3715SXin Li<tr style="border-bottom-style: none;">
256*bf2c3715SXin Li  <td>
257*bf2c3715SXin Li\code
258*bf2c3715SXin Lix = FixedXD::Identity();
259*bf2c3715SXin Lix.setIdentity();
260*bf2c3715SXin Li
261*bf2c3715SXin LiVector3f::UnitX() // 1 0 0
262*bf2c3715SXin LiVector3f::UnitY() // 0 1 0
263*bf2c3715SXin LiVector3f::UnitZ() // 0 0 1
264*bf2c3715SXin LiVector4f::Unit(i)
265*bf2c3715SXin Lix.setUnit(i);
266*bf2c3715SXin Li\endcode
267*bf2c3715SXin Li  </td>
268*bf2c3715SXin Li  <td>
269*bf2c3715SXin Li\code
270*bf2c3715SXin Lix = Dynamic2D::Identity(rows, cols);
271*bf2c3715SXin Lix.setIdentity(rows, cols);
272*bf2c3715SXin Li
273*bf2c3715SXin Li
274*bf2c3715SXin Li
275*bf2c3715SXin LiN/A
276*bf2c3715SXin Li\endcode
277*bf2c3715SXin Li  </td>
278*bf2c3715SXin Li  <td>\code
279*bf2c3715SXin LiN/A
280*bf2c3715SXin Li
281*bf2c3715SXin Li
282*bf2c3715SXin LiVectorXf::Unit(size,i)
283*bf2c3715SXin Lix.setUnit(size,i);
284*bf2c3715SXin LiVectorXf::Unit(4,1) == Vector4f(0,1,0,0)
285*bf2c3715SXin Li                    == Vector4f::UnitY()
286*bf2c3715SXin Li\endcode
287*bf2c3715SXin Li  </td>
288*bf2c3715SXin Li</tr>
289*bf2c3715SXin Li</table>
290*bf2c3715SXin Li
291*bf2c3715SXin LiNote that it is allowed to call any of the \c set* functions to a dynamic-sized vector or matrix without passing new sizes.
292*bf2c3715SXin LiFor instance:
293*bf2c3715SXin Li\code
294*bf2c3715SXin LiMatrixXi M(3,3);
295*bf2c3715SXin LiM.setIdentity();
296*bf2c3715SXin Li\endcode
297*bf2c3715SXin Li
298*bf2c3715SXin Li\subsection QuickRef_Map Mapping external arrays
299*bf2c3715SXin Li
300*bf2c3715SXin Li<table class="manual">
301*bf2c3715SXin Li<tr>
302*bf2c3715SXin Li<td>Contiguous \n memory</td>
303*bf2c3715SXin Li<td>\code
304*bf2c3715SXin Lifloat data[] = {1,2,3,4};
305*bf2c3715SXin LiMap<Vector3f> v1(data);       // uses v1 as a Vector3f object
306*bf2c3715SXin LiMap<ArrayXf>  v2(data,3);     // uses v2 as a ArrayXf object
307*bf2c3715SXin LiMap<Array22f> m1(data);       // uses m1 as a Array22f object
308*bf2c3715SXin LiMap<MatrixXf> m2(data,2,2);   // uses m2 as a MatrixXf object
309*bf2c3715SXin Li\endcode</td>
310*bf2c3715SXin Li</tr>
311*bf2c3715SXin Li<tr>
312*bf2c3715SXin Li<td>Typical usage \n of strides</td>
313*bf2c3715SXin Li<td>\code
314*bf2c3715SXin Lifloat data[] = {1,2,3,4,5,6,7,8,9};
315*bf2c3715SXin LiMap<VectorXf,0,InnerStride<2> >  v1(data,3);                      // = [1,3,5]
316*bf2c3715SXin LiMap<VectorXf,0,InnerStride<> >   v2(data,3,InnerStride<>(3));     // = [1,4,7]
317*bf2c3715SXin LiMap<MatrixXf,0,OuterStride<3> >  m2(data,2,3);                    // both lines     |1,4,7|
318*bf2c3715SXin LiMap<MatrixXf,0,OuterStride<> >   m1(data,2,3,OuterStride<>(3));   // are equal to:  |2,5,8|
319*bf2c3715SXin Li\endcode</td>
320*bf2c3715SXin Li</tr>
321*bf2c3715SXin Li</table>
322*bf2c3715SXin Li
323*bf2c3715SXin Li
324*bf2c3715SXin Li<a href="#" class="top">top</a>
325*bf2c3715SXin Li\section QuickRef_ArithmeticOperators Arithmetic Operators
326*bf2c3715SXin Li
327*bf2c3715SXin Li<table class="manual">
328*bf2c3715SXin Li<tr><td>
329*bf2c3715SXin Liadd \n subtract</td><td>\code
330*bf2c3715SXin Limat3 = mat1 + mat2;           mat3 += mat1;
331*bf2c3715SXin Limat3 = mat1 - mat2;           mat3 -= mat1;\endcode
332*bf2c3715SXin Li</td></tr>
333*bf2c3715SXin Li<tr class="alt"><td>
334*bf2c3715SXin Liscalar product</td><td>\code
335*bf2c3715SXin Limat3 = mat1 * s1;             mat3 *= s1;           mat3 = s1 * mat1;
336*bf2c3715SXin Limat3 = mat1 / s1;             mat3 /= s1;\endcode
337*bf2c3715SXin Li</td></tr>
338*bf2c3715SXin Li<tr><td>
339*bf2c3715SXin Limatrix/vector \n products \matrixworld</td><td>\code
340*bf2c3715SXin Licol2 = mat1 * col1;
341*bf2c3715SXin Lirow2 = row1 * mat1;           row1 *= mat1;
342*bf2c3715SXin Limat3 = mat1 * mat2;           mat3 *= mat1; \endcode
343*bf2c3715SXin Li</td></tr>
344*bf2c3715SXin Li<tr class="alt"><td>
345*bf2c3715SXin Litransposition \n adjoint \matrixworld</td><td>\code
346*bf2c3715SXin Limat1 = mat2.transpose();      mat1.transposeInPlace();
347*bf2c3715SXin Limat1 = mat2.adjoint();        mat1.adjointInPlace();
348*bf2c3715SXin Li\endcode
349*bf2c3715SXin Li</td></tr>
350*bf2c3715SXin Li<tr><td>
351*bf2c3715SXin Li\link MatrixBase::dot dot \endlink product \n inner product \matrixworld</td><td>\code
352*bf2c3715SXin Liscalar = vec1.dot(vec2);
353*bf2c3715SXin Liscalar = col1.adjoint() * col2;
354*bf2c3715SXin Liscalar = (col1.adjoint() * col2).value();\endcode
355*bf2c3715SXin Li</td></tr>
356*bf2c3715SXin Li<tr class="alt"><td>
357*bf2c3715SXin Liouter product \matrixworld</td><td>\code
358*bf2c3715SXin Limat = col1 * col2.transpose();\endcode
359*bf2c3715SXin Li</td></tr>
360*bf2c3715SXin Li
361*bf2c3715SXin Li<tr><td>
362*bf2c3715SXin Li\link MatrixBase::norm() norm \endlink \n \link MatrixBase::normalized() normalization \endlink \matrixworld</td><td>\code
363*bf2c3715SXin Liscalar = vec1.norm();         scalar = vec1.squaredNorm()
364*bf2c3715SXin Livec2 = vec1.normalized();     vec1.normalize(); // inplace \endcode
365*bf2c3715SXin Li</td></tr>
366*bf2c3715SXin Li
367*bf2c3715SXin Li<tr class="alt"><td>
368*bf2c3715SXin Li\link MatrixBase::cross() cross product \endlink \matrixworld</td><td>\code
369*bf2c3715SXin Li#include <Eigen/Geometry>
370*bf2c3715SXin Livec3 = vec1.cross(vec2);\endcode</td></tr>
371*bf2c3715SXin Li</table>
372*bf2c3715SXin Li
373*bf2c3715SXin Li<a href="#" class="top">top</a>
374*bf2c3715SXin Li\section QuickRef_Coeffwise Coefficient-wise \& Array operators
375*bf2c3715SXin Li
376*bf2c3715SXin LiIn addition to the aforementioned operators, Eigen supports numerous coefficient-wise operator and functions.
377*bf2c3715SXin LiMost of them unambiguously makes sense in array-world\arrayworld. The following operators are readily available for arrays,
378*bf2c3715SXin Lior available through .array() for vectors and matrices:
379*bf2c3715SXin Li
380*bf2c3715SXin Li<table class="manual">
381*bf2c3715SXin Li<tr><td>Arithmetic operators</td><td>\code
382*bf2c3715SXin Liarray1 * array2     array1 / array2     array1 *= array2    array1 /= array2
383*bf2c3715SXin Liarray1 + scalar     array1 - scalar     array1 += scalar    array1 -= scalar
384*bf2c3715SXin Li\endcode</td></tr>
385*bf2c3715SXin Li<tr><td>Comparisons</td><td>\code
386*bf2c3715SXin Liarray1 < array2     array1 > array2     array1 < scalar     array1 > scalar
387*bf2c3715SXin Liarray1 <= array2    array1 >= array2    array1 <= scalar    array1 >= scalar
388*bf2c3715SXin Liarray1 == array2    array1 != array2    array1 == scalar    array1 != scalar
389*bf2c3715SXin Liarray1.min(array2)  array1.max(array2)  array1.min(scalar)  array1.max(scalar)
390*bf2c3715SXin Li\endcode</td></tr>
391*bf2c3715SXin Li<tr><td>Trigo, power, and \n misc functions \n and the STL-like variants</td><td>\code
392*bf2c3715SXin Liarray1.abs2()
393*bf2c3715SXin Liarray1.abs()                  abs(array1)
394*bf2c3715SXin Liarray1.sqrt()                 sqrt(array1)
395*bf2c3715SXin Liarray1.log()                  log(array1)
396*bf2c3715SXin Liarray1.log10()                log10(array1)
397*bf2c3715SXin Liarray1.exp()                  exp(array1)
398*bf2c3715SXin Liarray1.pow(array2)            pow(array1,array2)
399*bf2c3715SXin Liarray1.pow(scalar)            pow(array1,scalar)
400*bf2c3715SXin Li                              pow(scalar,array2)
401*bf2c3715SXin Liarray1.square()
402*bf2c3715SXin Liarray1.cube()
403*bf2c3715SXin Liarray1.inverse()
404*bf2c3715SXin Li
405*bf2c3715SXin Liarray1.sin()                  sin(array1)
406*bf2c3715SXin Liarray1.cos()                  cos(array1)
407*bf2c3715SXin Liarray1.tan()                  tan(array1)
408*bf2c3715SXin Liarray1.asin()                 asin(array1)
409*bf2c3715SXin Liarray1.acos()                 acos(array1)
410*bf2c3715SXin Liarray1.atan()                 atan(array1)
411*bf2c3715SXin Liarray1.sinh()                 sinh(array1)
412*bf2c3715SXin Liarray1.cosh()                 cosh(array1)
413*bf2c3715SXin Liarray1.tanh()                 tanh(array1)
414*bf2c3715SXin Liarray1.arg()                  arg(array1)
415*bf2c3715SXin Li
416*bf2c3715SXin Liarray1.floor()                floor(array1)
417*bf2c3715SXin Liarray1.ceil()                 ceil(array1)
418*bf2c3715SXin Liarray1.round()                round(aray1)
419*bf2c3715SXin Li
420*bf2c3715SXin Liarray1.isFinite()             isfinite(array1)
421*bf2c3715SXin Liarray1.isInf()                isinf(array1)
422*bf2c3715SXin Liarray1.isNaN()                isnan(array1)
423*bf2c3715SXin Li\endcode
424*bf2c3715SXin Li</td></tr>
425*bf2c3715SXin Li</table>
426*bf2c3715SXin Li
427*bf2c3715SXin Li
428*bf2c3715SXin LiThe following coefficient-wise operators are available for all kind of expressions (matrices, vectors, and arrays), and for both real or complex scalar types:
429*bf2c3715SXin Li
430*bf2c3715SXin Li<table class="manual">
431*bf2c3715SXin Li<tr><th>Eigen's API</th><th>STL-like APIs\arrayworld </th><th>Comments</th></tr>
432*bf2c3715SXin Li<tr><td>\code
433*bf2c3715SXin Limat1.real()
434*bf2c3715SXin Limat1.imag()
435*bf2c3715SXin Limat1.conjugate()
436*bf2c3715SXin Li\endcode
437*bf2c3715SXin Li</td><td>\code
438*bf2c3715SXin Lireal(array1)
439*bf2c3715SXin Liimag(array1)
440*bf2c3715SXin Liconj(array1)
441*bf2c3715SXin Li\endcode
442*bf2c3715SXin Li</td><td>
443*bf2c3715SXin Li\code
444*bf2c3715SXin Li // read-write, no-op for real expressions
445*bf2c3715SXin Li // read-only for real, read-write for complexes
446*bf2c3715SXin Li // no-op for real expressions
447*bf2c3715SXin Li\endcode
448*bf2c3715SXin Li</td></tr>
449*bf2c3715SXin Li</table>
450*bf2c3715SXin Li
451*bf2c3715SXin LiSome coefficient-wise operators are readily available for for matrices and vectors through the following cwise* methods:
452*bf2c3715SXin Li<table class="manual">
453*bf2c3715SXin Li<tr><th>Matrix API \matrixworld</th><th>Via Array conversions</th></tr>
454*bf2c3715SXin Li<tr><td>\code
455*bf2c3715SXin Limat1.cwiseMin(mat2)         mat1.cwiseMin(scalar)
456*bf2c3715SXin Limat1.cwiseMax(mat2)         mat1.cwiseMax(scalar)
457*bf2c3715SXin Limat1.cwiseAbs2()
458*bf2c3715SXin Limat1.cwiseAbs()
459*bf2c3715SXin Limat1.cwiseSqrt()
460*bf2c3715SXin Limat1.cwiseInverse()
461*bf2c3715SXin Limat1.cwiseProduct(mat2)
462*bf2c3715SXin Limat1.cwiseQuotient(mat2)
463*bf2c3715SXin Limat1.cwiseEqual(mat2)       mat1.cwiseEqual(scalar)
464*bf2c3715SXin Limat1.cwiseNotEqual(mat2)
465*bf2c3715SXin Li\endcode
466*bf2c3715SXin Li</td><td>\code
467*bf2c3715SXin Limat1.array().min(mat2.array())    mat1.array().min(scalar)
468*bf2c3715SXin Limat1.array().max(mat2.array())    mat1.array().max(scalar)
469*bf2c3715SXin Limat1.array().abs2()
470*bf2c3715SXin Limat1.array().abs()
471*bf2c3715SXin Limat1.array().sqrt()
472*bf2c3715SXin Limat1.array().inverse()
473*bf2c3715SXin Limat1.array() * mat2.array()
474*bf2c3715SXin Limat1.array() / mat2.array()
475*bf2c3715SXin Limat1.array() == mat2.array()      mat1.array() == scalar
476*bf2c3715SXin Limat1.array() != mat2.array()
477*bf2c3715SXin Li\endcode</td></tr>
478*bf2c3715SXin Li</table>
479*bf2c3715SXin LiThe main difference between the two API is that the one based on cwise* methods returns an expression in the matrix world,
480*bf2c3715SXin Liwhile the second one (based on .array()) returns an array expression.
481*bf2c3715SXin LiRecall that .array() has no cost, it only changes the available API and interpretation of the data.
482*bf2c3715SXin Li
483*bf2c3715SXin LiIt is also very simple to apply any user defined function \c foo using DenseBase::unaryExpr together with <a href="http://en.cppreference.com/w/cpp/utility/functional/ptr_fun">std::ptr_fun</a> (c++03, deprecated or removed in newer C++ versions), <a href="http://en.cppreference.com/w/cpp/utility/functional/ref">std::ref</a> (c++11), or <a href="http://en.cppreference.com/w/cpp/language/lambda">lambdas</a> (c++11):
484*bf2c3715SXin Li\code
485*bf2c3715SXin Limat1.unaryExpr(std::ptr_fun(foo));
486*bf2c3715SXin Limat1.unaryExpr(std::ref(foo));
487*bf2c3715SXin Limat1.unaryExpr([](double x) { return foo(x); });
488*bf2c3715SXin Li\endcode
489*bf2c3715SXin Li
490*bf2c3715SXin LiPlease note that it's not possible to pass a raw function pointer to \c unaryExpr, so please warp it as shown above.
491*bf2c3715SXin Li
492*bf2c3715SXin Li<a href="#" class="top">top</a>
493*bf2c3715SXin Li\section QuickRef_Reductions Reductions
494*bf2c3715SXin Li
495*bf2c3715SXin LiEigen provides several reduction methods such as:
496*bf2c3715SXin Li\link DenseBase::minCoeff() minCoeff() \endlink, \link DenseBase::maxCoeff() maxCoeff() \endlink,
497*bf2c3715SXin Li\link DenseBase::sum() sum() \endlink, \link DenseBase::prod() prod() \endlink,
498*bf2c3715SXin Li\link MatrixBase::trace() trace() \endlink \matrixworld,
499*bf2c3715SXin Li\link MatrixBase::norm() norm() \endlink \matrixworld, \link MatrixBase::squaredNorm() squaredNorm() \endlink \matrixworld,
500*bf2c3715SXin Li\link DenseBase::all() all() \endlink, and \link DenseBase::any() any() \endlink.
501*bf2c3715SXin LiAll reduction operations can be done matrix-wise,
502*bf2c3715SXin Li\link DenseBase::colwise() column-wise \endlink or
503*bf2c3715SXin Li\link DenseBase::rowwise() row-wise \endlink. Usage example:
504*bf2c3715SXin Li<table class="manual">
505*bf2c3715SXin Li<tr><td rowspan="3" style="border-right-style:dashed;vertical-align:middle">\code
506*bf2c3715SXin Li      5 3 1
507*bf2c3715SXin Limat = 2 7 8
508*bf2c3715SXin Li      9 4 6 \endcode
509*bf2c3715SXin Li</td> <td>\code mat.minCoeff(); \endcode</td><td>\code 1 \endcode</td></tr>
510*bf2c3715SXin Li<tr class="alt"><td>\code mat.colwise().minCoeff(); \endcode</td><td>\code 2 3 1 \endcode</td></tr>
511*bf2c3715SXin Li<tr style="vertical-align:middle"><td>\code mat.rowwise().minCoeff(); \endcode</td><td>\code
512*bf2c3715SXin Li1
513*bf2c3715SXin Li2
514*bf2c3715SXin Li4
515*bf2c3715SXin Li\endcode</td></tr>
516*bf2c3715SXin Li</table>
517*bf2c3715SXin Li
518*bf2c3715SXin LiSpecial versions of \link DenseBase::minCoeff(IndexType*,IndexType*) const minCoeff \endlink and \link DenseBase::maxCoeff(IndexType*,IndexType*) const maxCoeff \endlink:
519*bf2c3715SXin Li\code
520*bf2c3715SXin Liint i, j;
521*bf2c3715SXin Lis = vector.minCoeff(&i);        // s == vector[i]
522*bf2c3715SXin Lis = matrix.maxCoeff(&i, &j);    // s == matrix(i,j)
523*bf2c3715SXin Li\endcode
524*bf2c3715SXin LiTypical use cases of all() and any():
525*bf2c3715SXin Li\code
526*bf2c3715SXin Liif((array1 > 0).all()) ...      // if all coefficients of array1 are greater than 0 ...
527*bf2c3715SXin Liif((array1 < array2).any()) ... // if there exist a pair i,j such that array1(i,j) < array2(i,j) ...
528*bf2c3715SXin Li\endcode
529*bf2c3715SXin Li
530*bf2c3715SXin Li
531*bf2c3715SXin Li<a href="#" class="top">top</a>\section QuickRef_Blocks Sub-matrices
532*bf2c3715SXin Li
533*bf2c3715SXin Li<div class="warningbox">
534*bf2c3715SXin Li<strong>PLEASE HELP US IMPROVING THIS SECTION.</strong>
535*bf2c3715SXin Li%Eigen 3.4 supports a much improved API for sub-matrices, including,
536*bf2c3715SXin Lislicing and indexing from arrays: \ref TutorialSlicingIndexing
537*bf2c3715SXin Li</div>
538*bf2c3715SXin Li
539*bf2c3715SXin LiRead-write access to a \link DenseBase::col(Index) column \endlink
540*bf2c3715SXin Lior a \link DenseBase::row(Index) row \endlink of a matrix (or array):
541*bf2c3715SXin Li\code
542*bf2c3715SXin Limat1.row(i) = mat2.col(j);
543*bf2c3715SXin Limat1.col(j1).swap(mat1.col(j2));
544*bf2c3715SXin Li\endcode
545*bf2c3715SXin Li
546*bf2c3715SXin LiRead-write access to sub-vectors:
547*bf2c3715SXin Li<table class="manual">
548*bf2c3715SXin Li<tr>
549*bf2c3715SXin Li<th>Default versions</th>
550*bf2c3715SXin Li<th>Optimized versions when the size \n is known at compile time</th></tr>
551*bf2c3715SXin Li<th></th>
552*bf2c3715SXin Li
553*bf2c3715SXin Li<tr><td>\code vec1.head(n)\endcode</td><td>\code vec1.head<n>()\endcode</td><td>the first \c n coeffs </td></tr>
554*bf2c3715SXin Li<tr><td>\code vec1.tail(n)\endcode</td><td>\code vec1.tail<n>()\endcode</td><td>the last \c n coeffs </td></tr>
555*bf2c3715SXin Li<tr><td>\code vec1.segment(pos,n)\endcode</td><td>\code vec1.segment<n>(pos)\endcode</td>
556*bf2c3715SXin Li    <td>the \c n coeffs in the \n range [\c pos : \c pos + \c n - 1]</td></tr>
557*bf2c3715SXin Li<tr class="alt"><td colspan="3">
558*bf2c3715SXin Li
559*bf2c3715SXin LiRead-write access to sub-matrices:</td></tr>
560*bf2c3715SXin Li<tr>
561*bf2c3715SXin Li  <td>\code mat1.block(i,j,rows,cols)\endcode
562*bf2c3715SXin Li      \link DenseBase::block(Index,Index,Index,Index) (more) \endlink</td>
563*bf2c3715SXin Li  <td>\code mat1.block<rows,cols>(i,j)\endcode
564*bf2c3715SXin Li      \link DenseBase::block(Index,Index) (more) \endlink</td>
565*bf2c3715SXin Li  <td>the \c rows x \c cols sub-matrix \n starting from position (\c i,\c j)</td></tr>
566*bf2c3715SXin Li<tr><td>\code
567*bf2c3715SXin Li mat1.topLeftCorner(rows,cols)
568*bf2c3715SXin Li mat1.topRightCorner(rows,cols)
569*bf2c3715SXin Li mat1.bottomLeftCorner(rows,cols)
570*bf2c3715SXin Li mat1.bottomRightCorner(rows,cols)\endcode
571*bf2c3715SXin Li <td>\code
572*bf2c3715SXin Li mat1.topLeftCorner<rows,cols>()
573*bf2c3715SXin Li mat1.topRightCorner<rows,cols>()
574*bf2c3715SXin Li mat1.bottomLeftCorner<rows,cols>()
575*bf2c3715SXin Li mat1.bottomRightCorner<rows,cols>()\endcode
576*bf2c3715SXin Li <td>the \c rows x \c cols sub-matrix \n taken in one of the four corners</td></tr>
577*bf2c3715SXin Li <tr><td>\code
578*bf2c3715SXin Li mat1.topRows(rows)
579*bf2c3715SXin Li mat1.bottomRows(rows)
580*bf2c3715SXin Li mat1.leftCols(cols)
581*bf2c3715SXin Li mat1.rightCols(cols)\endcode
582*bf2c3715SXin Li <td>\code
583*bf2c3715SXin Li mat1.topRows<rows>()
584*bf2c3715SXin Li mat1.bottomRows<rows>()
585*bf2c3715SXin Li mat1.leftCols<cols>()
586*bf2c3715SXin Li mat1.rightCols<cols>()\endcode
587*bf2c3715SXin Li <td>specialized versions of block() \n when the block fit two corners</td></tr>
588*bf2c3715SXin Li</table>
589*bf2c3715SXin Li
590*bf2c3715SXin Li
591*bf2c3715SXin Li
592*bf2c3715SXin Li<a href="#" class="top">top</a>\section QuickRef_Misc Miscellaneous operations
593*bf2c3715SXin Li
594*bf2c3715SXin Li<div class="warningbox">
595*bf2c3715SXin Li<strong>PLEASE HELP US IMPROVING THIS SECTION.</strong>
596*bf2c3715SXin Li%Eigen 3.4 supports a new API for reshaping: \ref TutorialReshape
597*bf2c3715SXin Li</div>
598*bf2c3715SXin Li
599*bf2c3715SXin Li\subsection QuickRef_Reverse Reverse
600*bf2c3715SXin LiVectors, rows, and/or columns of a matrix can be reversed (see DenseBase::reverse(), DenseBase::reverseInPlace(), VectorwiseOp::reverse()).
601*bf2c3715SXin Li\code
602*bf2c3715SXin Livec.reverse()           mat.colwise().reverse()   mat.rowwise().reverse()
603*bf2c3715SXin Livec.reverseInPlace()
604*bf2c3715SXin Li\endcode
605*bf2c3715SXin Li
606*bf2c3715SXin Li\subsection QuickRef_Replicate Replicate
607*bf2c3715SXin LiVectors, matrices, rows, and/or columns can be replicated in any direction (see DenseBase::replicate(), VectorwiseOp::replicate())
608*bf2c3715SXin Li\code
609*bf2c3715SXin Livec.replicate(times)                                          vec.replicate<Times>
610*bf2c3715SXin Limat.replicate(vertical_times, horizontal_times)               mat.replicate<VerticalTimes, HorizontalTimes>()
611*bf2c3715SXin Limat.colwise().replicate(vertical_times, horizontal_times)     mat.colwise().replicate<VerticalTimes, HorizontalTimes>()
612*bf2c3715SXin Limat.rowwise().replicate(vertical_times, horizontal_times)     mat.rowwise().replicate<VerticalTimes, HorizontalTimes>()
613*bf2c3715SXin Li\endcode
614*bf2c3715SXin Li
615*bf2c3715SXin Li
616*bf2c3715SXin Li<a href="#" class="top">top</a>\section QuickRef_DiagTriSymm Diagonal, Triangular, and Self-adjoint matrices
617*bf2c3715SXin Li(matrix world \matrixworld)
618*bf2c3715SXin Li
619*bf2c3715SXin Li\subsection QuickRef_Diagonal Diagonal matrices
620*bf2c3715SXin Li
621*bf2c3715SXin Li<table class="example">
622*bf2c3715SXin Li<tr><th>Operation</th><th>Code</th></tr>
623*bf2c3715SXin Li<tr><td>
624*bf2c3715SXin Liview a vector \link MatrixBase::asDiagonal() as a diagonal matrix \endlink \n </td><td>\code
625*bf2c3715SXin Limat1 = vec1.asDiagonal();\endcode
626*bf2c3715SXin Li</td></tr>
627*bf2c3715SXin Li<tr><td>
628*bf2c3715SXin LiDeclare a diagonal matrix</td><td>\code
629*bf2c3715SXin LiDiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
630*bf2c3715SXin Lidiag1.diagonal() = vector;\endcode
631*bf2c3715SXin Li</td></tr>
632*bf2c3715SXin Li<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>
633*bf2c3715SXin Li <td>\code
634*bf2c3715SXin Livec1 = mat1.diagonal();        mat1.diagonal() = vec1;      // main diagonal
635*bf2c3715SXin Livec1 = mat1.diagonal(+n);      mat1.diagonal(+n) = vec1;    // n-th super diagonal
636*bf2c3715SXin Livec1 = mat1.diagonal(-n);      mat1.diagonal(-n) = vec1;    // n-th sub diagonal
637*bf2c3715SXin Livec1 = mat1.diagonal<1>();     mat1.diagonal<1>() = vec1;   // first super diagonal
638*bf2c3715SXin Livec1 = mat1.diagonal<-2>();    mat1.diagonal<-2>() = vec1;  // second sub diagonal
639*bf2c3715SXin Li\endcode</td>
640*bf2c3715SXin Li</tr>
641*bf2c3715SXin Li
642*bf2c3715SXin Li<tr><td>Optimized products and inverse</td>
643*bf2c3715SXin Li <td>\code
644*bf2c3715SXin Limat3  = scalar * diag1 * mat1;
645*bf2c3715SXin Limat3 += scalar * mat1 * vec1.asDiagonal();
646*bf2c3715SXin Limat3 = vec1.asDiagonal().inverse() * mat1
647*bf2c3715SXin Limat3 = mat1 * diag1.inverse()
648*bf2c3715SXin Li\endcode</td>
649*bf2c3715SXin Li</tr>
650*bf2c3715SXin Li
651*bf2c3715SXin Li</table>
652*bf2c3715SXin Li
653*bf2c3715SXin Li\subsection QuickRef_TriangularView Triangular views
654*bf2c3715SXin Li
655*bf2c3715SXin LiTriangularView 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.
656*bf2c3715SXin Li
657*bf2c3715SXin Li\note The .triangularView() template member function requires the \c template keyword if it is used on an
658*bf2c3715SXin Liobject of a type that depends on a template parameter; see \ref TopicTemplateKeyword for details.
659*bf2c3715SXin Li
660*bf2c3715SXin Li<table class="example">
661*bf2c3715SXin Li<tr><th>Operation</th><th>Code</th></tr>
662*bf2c3715SXin Li<tr><td>
663*bf2c3715SXin LiReference to a triangular with optional \n
664*bf2c3715SXin Liunit or null diagonal (read/write):
665*bf2c3715SXin Li</td><td>\code
666*bf2c3715SXin Lim.triangularView<Xxx>()
667*bf2c3715SXin Li\endcode \n
668*bf2c3715SXin Li\c Xxx = ::Upper, ::Lower, ::StrictlyUpper, ::StrictlyLower, ::UnitUpper, ::UnitLower
669*bf2c3715SXin Li</td></tr>
670*bf2c3715SXin Li<tr><td>
671*bf2c3715SXin LiWriting to a specific triangular part:\n (only the referenced triangular part is evaluated)
672*bf2c3715SXin Li</td><td>\code
673*bf2c3715SXin Lim1.triangularView<Eigen::Lower>() = m2 + m3 \endcode
674*bf2c3715SXin Li</td></tr>
675*bf2c3715SXin Li<tr><td>
676*bf2c3715SXin LiConversion to a dense matrix setting the opposite triangular part to zero:
677*bf2c3715SXin Li</td><td>\code
678*bf2c3715SXin Lim2 = m1.triangularView<Eigen::UnitUpper>()\endcode
679*bf2c3715SXin Li</td></tr>
680*bf2c3715SXin Li<tr><td>
681*bf2c3715SXin LiProducts:
682*bf2c3715SXin Li</td><td>\code
683*bf2c3715SXin Lim3 += s1 * m1.adjoint().triangularView<Eigen::UnitUpper>() * m2
684*bf2c3715SXin Lim3 -= s1 * m2.conjugate() * m1.adjoint().triangularView<Eigen::Lower>() \endcode
685*bf2c3715SXin Li</td></tr>
686*bf2c3715SXin Li<tr><td>
687*bf2c3715SXin LiSolving linear equations:\n
688*bf2c3715SXin Li\f$ M_2 := L_1^{-1} M_2 \f$ \n
689*bf2c3715SXin Li\f$ M_3 := {L_1^*}^{-1} M_3 \f$ \n
690*bf2c3715SXin Li\f$ M_4 := M_4 U_1^{-1} \f$
691*bf2c3715SXin Li</td><td>\n \code
692*bf2c3715SXin LiL1.triangularView<Eigen::UnitLower>().solveInPlace(M2)
693*bf2c3715SXin LiL1.triangularView<Eigen::Lower>().adjoint().solveInPlace(M3)
694*bf2c3715SXin LiU1.triangularView<Eigen::Upper>().solveInPlace<OnTheRight>(M4)\endcode
695*bf2c3715SXin Li</td></tr>
696*bf2c3715SXin Li</table>
697*bf2c3715SXin Li
698*bf2c3715SXin Li\subsection QuickRef_SelfadjointMatrix Symmetric/selfadjoint views
699*bf2c3715SXin Li
700*bf2c3715SXin LiJust as for triangular matrix, you can reference any triangular part of a square matrix to see it as a selfadjoint
701*bf2c3715SXin Limatrix and perform special and optimized operations. Again the opposite triangular part is never referenced and can be
702*bf2c3715SXin Liused to store other information.
703*bf2c3715SXin Li
704*bf2c3715SXin Li\note The .selfadjointView() template member function requires the \c template keyword if it is used on an
705*bf2c3715SXin Liobject of a type that depends on a template parameter; see \ref TopicTemplateKeyword for details.
706*bf2c3715SXin Li
707*bf2c3715SXin Li<table class="example">
708*bf2c3715SXin Li<tr><th>Operation</th><th>Code</th></tr>
709*bf2c3715SXin Li<tr><td>
710*bf2c3715SXin LiConversion to a dense matrix:
711*bf2c3715SXin Li</td><td>\code
712*bf2c3715SXin Lim2 = m.selfadjointView<Eigen::Lower>();\endcode
713*bf2c3715SXin Li</td></tr>
714*bf2c3715SXin Li<tr><td>
715*bf2c3715SXin LiProduct with another general matrix or vector:
716*bf2c3715SXin Li</td><td>\code
717*bf2c3715SXin Lim3  = s1 * m1.conjugate().selfadjointView<Eigen::Upper>() * m3;
718*bf2c3715SXin Lim3 -= s1 * m3.adjoint() * m1.selfadjointView<Eigen::Lower>();\endcode
719*bf2c3715SXin Li</td></tr>
720*bf2c3715SXin Li<tr><td>
721*bf2c3715SXin LiRank 1 and rank K update: \n
722*bf2c3715SXin Li\f$ upper(M_1) \mathrel{{+}{=}} s_1 M_2 M_2^* \f$ \n
723*bf2c3715SXin Li\f$ lower(M_1) \mathbin{{-}{=}} M_2^* M_2 \f$
724*bf2c3715SXin Li</td><td>\n \code
725*bf2c3715SXin LiM1.selfadjointView<Eigen::Upper>().rankUpdate(M2,s1);
726*bf2c3715SXin LiM1.selfadjointView<Eigen::Lower>().rankUpdate(M2.adjoint(),-1); \endcode
727*bf2c3715SXin Li</td></tr>
728*bf2c3715SXin Li<tr><td>
729*bf2c3715SXin LiRank 2 update: (\f$ M \mathrel{{+}{=}} s u v^* + s v u^* \f$)
730*bf2c3715SXin Li</td><td>\code
731*bf2c3715SXin LiM.selfadjointView<Eigen::Upper>().rankUpdate(u,v,s);
732*bf2c3715SXin Li\endcode
733*bf2c3715SXin Li</td></tr>
734*bf2c3715SXin Li<tr><td>
735*bf2c3715SXin LiSolving linear equations:\n(\f$ M_2 := M_1^{-1} M_2 \f$)
736*bf2c3715SXin Li</td><td>\code
737*bf2c3715SXin Li// via a standard Cholesky factorization
738*bf2c3715SXin Lim2 = m1.selfadjointView<Eigen::Upper>().llt().solve(m2);
739*bf2c3715SXin Li// via a Cholesky factorization with pivoting
740*bf2c3715SXin Lim2 = m1.selfadjointView<Eigen::Lower>().ldlt().solve(m2);
741*bf2c3715SXin Li\endcode
742*bf2c3715SXin Li</td></tr>
743*bf2c3715SXin Li</table>
744*bf2c3715SXin Li
745*bf2c3715SXin Li*/
746*bf2c3715SXin Li
747*bf2c3715SXin Li/*
748*bf2c3715SXin Li<table class="tutorial_code">
749*bf2c3715SXin Li<tr><td>
750*bf2c3715SXin Li\link MatrixBase::asDiagonal() make a diagonal matrix \endlink \n from a vector </td><td>\code
751*bf2c3715SXin Limat1 = vec1.asDiagonal();\endcode
752*bf2c3715SXin Li</td></tr>
753*bf2c3715SXin Li<tr><td>
754*bf2c3715SXin LiDeclare a diagonal matrix</td><td>\code
755*bf2c3715SXin LiDiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
756*bf2c3715SXin Lidiag1.diagonal() = vector;\endcode
757*bf2c3715SXin Li</td></tr>
758*bf2c3715SXin Li<tr><td>Access \link MatrixBase::diagonal() the diagonal and super/sub diagonals of a matrix \endlink as a vector (read/write)</td>
759*bf2c3715SXin Li <td>\code
760*bf2c3715SXin Livec1 = mat1.diagonal();            mat1.diagonal() = vec1;      // main diagonal
761*bf2c3715SXin Livec1 = mat1.diagonal(+n);          mat1.diagonal(+n) = vec1;    // n-th super diagonal
762*bf2c3715SXin Livec1 = mat1.diagonal(-n);          mat1.diagonal(-n) = vec1;    // n-th sub diagonal
763*bf2c3715SXin Livec1 = mat1.diagonal<1>();         mat1.diagonal<1>() = vec1;   // first super diagonal
764*bf2c3715SXin Livec1 = mat1.diagonal<-2>();        mat1.diagonal<-2>() = vec1;  // second sub diagonal
765*bf2c3715SXin Li\endcode</td>
766*bf2c3715SXin Li</tr>
767*bf2c3715SXin Li
768*bf2c3715SXin Li<tr><td>View on a triangular part of a matrix (read/write)</td>
769*bf2c3715SXin Li <td>\code
770*bf2c3715SXin Limat2 = mat1.triangularView<Xxx>();
771*bf2c3715SXin Li// Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower
772*bf2c3715SXin Limat1.triangularView<Upper>() = mat2 + mat3; // only the upper part is evaluated and referenced
773*bf2c3715SXin Li\endcode</td></tr>
774*bf2c3715SXin Li
775*bf2c3715SXin Li<tr><td>View a triangular part as a symmetric/self-adjoint matrix (read/write)</td>
776*bf2c3715SXin Li <td>\code
777*bf2c3715SXin Limat2 = mat1.selfadjointView<Xxx>();     // Xxx = Upper or Lower
778*bf2c3715SXin Limat1.selfadjointView<Upper>() = mat2 + mat2.adjoint();  // evaluated and write to the upper triangular part only
779*bf2c3715SXin Li\endcode</td></tr>
780*bf2c3715SXin Li
781*bf2c3715SXin Li</table>
782*bf2c3715SXin Li
783*bf2c3715SXin LiOptimized products:
784*bf2c3715SXin Li\code
785*bf2c3715SXin Limat3 += scalar * vec1.asDiagonal() * mat1
786*bf2c3715SXin Limat3 += scalar * mat1 * vec1.asDiagonal()
787*bf2c3715SXin Limat3.noalias() += scalar * mat1.triangularView<Xxx>() * mat2
788*bf2c3715SXin Limat3.noalias() += scalar * mat2 * mat1.triangularView<Xxx>()
789*bf2c3715SXin Limat3.noalias() += scalar * mat1.selfadjointView<Upper or Lower>() * mat2
790*bf2c3715SXin Limat3.noalias() += scalar * mat2 * mat1.selfadjointView<Upper or Lower>()
791*bf2c3715SXin Limat1.selfadjointView<Upper or Lower>().rankUpdate(mat2);
792*bf2c3715SXin Limat1.selfadjointView<Upper or Lower>().rankUpdate(mat2.adjoint(), scalar);
793*bf2c3715SXin Li\endcode
794*bf2c3715SXin Li
795*bf2c3715SXin LiInverse products: (all are optimized)
796*bf2c3715SXin Li\code
797*bf2c3715SXin Limat3 = vec1.asDiagonal().inverse() * mat1
798*bf2c3715SXin Limat3 = mat1 * diag1.inverse()
799*bf2c3715SXin Limat1.triangularView<Xxx>().solveInPlace(mat2)
800*bf2c3715SXin Limat1.triangularView<Xxx>().solveInPlace<OnTheRight>(mat2)
801*bf2c3715SXin Limat2 = mat1.selfadjointView<Upper or Lower>().llt().solve(mat2)
802*bf2c3715SXin Li\endcode
803*bf2c3715SXin Li
804*bf2c3715SXin Li*/
805*bf2c3715SXin Li}
806