3 /** \eigenManualPage QuickRefPage Quick reference guide
9 <a href="#" class="top">top</a>
10 \section QuickRef_Headers Modules and Header files
12 The 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.
14 <table class="manual">
15 <tr><th>Module</th><th>Header file</th><th>Contents</th></tr>
16 <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 <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 <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 <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 <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 <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 <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 <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 <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 <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 <tr ><td></td><td>\code#include <Eigen/Eigen>\endcode</td><td>Includes %Dense and %Sparse header files (the whole Eigen library)</td></tr>
29 <a href="#" class="top">top</a>
30 \section QuickRef_Types Array, matrix and vector types
33 \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:
35 typedef Matrix<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyMatrixType;
36 typedef Array<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyArrayType;
39 \li \c Scalar is the scalar type of the coefficients (e.g., \c float, \c double, \c bool, \c int, etc.).
40 \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 \li \c Options can be \c ColMajor or \c RowMajor, default is \c ColMajor. (see class Matrix for more options)
43 All 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:
45 Matrix<double, 6, Dynamic> // Dynamic number of columns (heap allocation)
46 Matrix<double, Dynamic, 2> // Dynamic number of rows (heap allocation)
47 Matrix<double, Dynamic, Dynamic, RowMajor> // Fully dynamic, row major (heap allocation)
48 Matrix<double, 13, 3> // Fully fixed (usually allocated on stack)
51 In most cases, you can simply use one of the convenience typedefs for \ref matrixtypedefs "matrices" and \ref arraytypedefs "arrays". Some examples:
52 <table class="example">
53 <tr><th>Matrices</th><th>Arrays</th></tr>
55 Matrix<float,Dynamic,Dynamic> <=> MatrixXf
56 Matrix<double,Dynamic,1> <=> VectorXd
57 Matrix<int,1,Dynamic> <=> RowVectorXi
58 Matrix<float,3,3> <=> Matrix3f
59 Matrix<float,4,1> <=> Vector4f
60 \endcode</td><td>\code
61 Array<float,Dynamic,Dynamic> <=> ArrayXXf
62 Array<double,Dynamic,1> <=> ArrayXd
63 Array<int,1,Dynamic> <=> RowArrayXi
64 Array<float,3,3> <=> Array33f
65 Array<float,4,1> <=> Array4f
69 Conversion between the matrix and array worlds:
73 m1 = a1 * a2; // coeffwise product, implicit conversion from array to matrix.
74 a1 = m1 * m2; // matrix product, implicit conversion from matrix to array.
75 a2 = a1 + m1.array(); // mixing array and matrix is forbidden
76 m2 = a1.matrix() + m1; // and explicit conversion is required.
77 ArrayWrapper<Matrix4f> m1a(m1); // m1a is an alias for m1.array(), they share the same coefficients
78 MatrixWrapper<Array44f> a1m(a1);
81 In 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 \li <a name="matrixonly"></a>\matrixworld linear algebra matrix and vector only
83 \li <a name="arrayonly"></a>\arrayworld array objects only
85 \subsection QuickRef_Basics Basic matrix manipulation
87 <table class="manual">
88 <tr><th></th><th>1D objects</th><th>2D objects</th><th>Notes</th></tr>
89 <tr><td>Constructors</td>
94 Vector4d v3(x, y, z, w);
96 VectorXf v5; // empty object
98 \endcode</td><td>\code
104 MatrixXf m5; // empty object
105 MatrixXf m6(nb_rows, nb_columns);
106 \endcode</td><td class="note">
107 By default, the coefficients \n are left uninitialized</td></tr>
108 <tr class="alt"><td>Comma initializer</td>
110 Vector3f v1; v1 << x, y, z;
111 ArrayXf v2(4); v2 << 1, 2, 3, 4;
113 \endcode</td><td>\code
114 Matrix3f m1; m1 << 1, 2, 3,
117 \endcode</td><td></td></tr>
119 <tr><td>Comma initializer (bis)</td>
121 \include Tutorial_commainit_02.cpp
125 \verbinclude Tutorial_commainit_02.out
129 <tr class="alt"><td>Runtime info</td>
133 vector.innerStride();
135 \endcode</td><td>\code
136 matrix.rows(); matrix.cols();
137 matrix.innerSize(); matrix.outerSize();
138 matrix.innerStride(); matrix.outerStride();
140 \endcode</td><td class="note">Inner/Outer* are storage order dependent</td></tr>
141 <tr><td>Compile-time info</td>
142 <td colspan="2">\code
143 ObjectType::Scalar ObjectType::RowsAtCompileTime
144 ObjectType::RealScalar ObjectType::ColsAtCompileTime
145 ObjectType::Index ObjectType::SizeAtCompileTime
146 \endcode</td><td></td></tr>
147 <tr class="alt"><td>Resizing</td>
152 vector.resizeLike(other_vector);
153 vector.conservativeResize(size);
154 \endcode</td><td>\code
155 matrix.resize(nb_rows, nb_cols);
156 matrix.resize(Eigen::NoChange, nb_cols);
157 matrix.resize(nb_rows, Eigen::NoChange);
158 matrix.resizeLike(other_matrix);
159 matrix.conservativeResize(nb_rows, nb_cols);
160 \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>
162 <tr><td>Coeff access with \n range checking</td>
168 \endcode</td><td>\code
170 \endcode</td><td class="note">Range checking is disabled if \n NDEBUG or EIGEN_NO_DEBUG is defined</td></tr>
172 <tr class="alt"><td>Coeff access without \n range checking</td>
176 \endcode</td><td>\code
179 \endcode</td><td></td></tr>
181 <tr><td>Assignment/copy</td>
182 <td colspan="2">\code
184 object_of_float = expression_of_double.cast<float>();
185 \endcode</td><td class="note">the destination is automatically resized (if possible)</td></tr>
189 \subsection QuickRef_PredefMat Predefined Matrices
191 <table class="manual">
193 <th>Fixed-size matrix or vector</th>
194 <th>Dynamic-size matrix</th>
195 <th>Dynamic-size vector</th>
197 <tr style="border-bottom-style: none;">
200 typedef {Matrix3f|Array33f} FixedXD;
205 x = FixedXD::Constant(value);
206 x = FixedXD::Random();
207 x = FixedXD::LinSpaced(size, low, high);
211 x.setConstant(value);
213 x.setLinSpaced(size, low, high);
218 typedef {MatrixXf|ArrayXXf} Dynamic2D;
221 x = Dynamic2D::Zero(rows, cols);
222 x = Dynamic2D::Ones(rows, cols);
223 x = Dynamic2D::Constant(rows, cols, value);
224 x = Dynamic2D::Random(rows, cols);
227 x.setZero(rows, cols);
228 x.setOnes(rows, cols);
229 x.setConstant(rows, cols, value);
230 x.setRandom(rows, cols);
236 typedef {VectorXf|ArrayXf} Dynamic1D;
239 x = Dynamic1D::Zero(size);
240 x = Dynamic1D::Ones(size);
241 x = Dynamic1D::Constant(size, value);
242 x = Dynamic1D::Random(size);
243 x = Dynamic1D::LinSpaced(size, low, high);
247 x.setConstant(size, value);
249 x.setLinSpaced(size, low, high);
254 <tr><td colspan="3">Identity and \link MatrixBase::Unit basis vectors \endlink \matrixworld</td></tr>
255 <tr style="border-bottom-style: none;">
258 x = FixedXD::Identity();
261 Vector3f::UnitX() // 1 0 0
262 Vector3f::UnitY() // 0 1 0
263 Vector3f::UnitZ() // 0 0 1
268 x = Dynamic2D::Identity(rows, cols);
269 x.setIdentity(rows, cols);
280 VectorXf::Unit(size,i)
281 VectorXf::Unit(4,1) == Vector4f(0,1,0,0)
290 \subsection QuickRef_Map Mapping external arrays
292 <table class="manual">
294 <td>Contiguous \n memory</td>
296 float data[] = {1,2,3,4};
297 Map<Vector3f> v1(data); // uses v1 as a Vector3f object
298 Map<ArrayXf> v2(data,3); // uses v2 as a ArrayXf object
299 Map<Array22f> m1(data); // uses m1 as a Array22f object
300 Map<MatrixXf> m2(data,2,2); // uses m2 as a MatrixXf object
304 <td>Typical usage \n of strides</td>
306 float data[] = {1,2,3,4,5,6,7,8,9};
307 Map<VectorXf,0,InnerStride<2> > v1(data,3); // = [1,3,5]
308 Map<VectorXf,0,InnerStride<> > v2(data,3,InnerStride<>(3)); // = [1,4,7]
309 Map<MatrixXf,0,OuterStride<3> > m2(data,2,3); // both lines |1,4,7|
310 Map<MatrixXf,0,OuterStride<> > m1(data,2,3,OuterStride<>(3)); // are equal to: |2,5,8|
316 <a href="#" class="top">top</a>
317 \section QuickRef_ArithmeticOperators Arithmetic Operators
319 <table class="manual">
321 add \n subtract</td><td>\code
322 mat3 = mat1 + mat2; mat3 += mat1;
323 mat3 = mat1 - mat2; mat3 -= mat1;\endcode
326 scalar product</td><td>\code
327 mat3 = mat1 * s1; mat3 *= s1; mat3 = s1 * mat1;
328 mat3 = mat1 / s1; mat3 /= s1;\endcode
331 matrix/vector \n products \matrixworld</td><td>\code
333 row2 = row1 * mat1; row1 *= mat1;
334 mat3 = mat1 * mat2; mat3 *= mat1; \endcode
337 transposition \n adjoint \matrixworld</td><td>\code
338 mat1 = mat2.transpose(); mat1.transposeInPlace();
339 mat1 = mat2.adjoint(); mat1.adjointInPlace();
343 \link MatrixBase::dot dot \endlink product \n inner product \matrixworld</td><td>\code
344 scalar = vec1.dot(vec2);
345 scalar = col1.adjoint() * col2;
346 scalar = (col1.adjoint() * col2).value();\endcode
349 outer product \matrixworld</td><td>\code
350 mat = col1 * col2.transpose();\endcode
354 \link MatrixBase::norm() norm \endlink \n \link MatrixBase::normalized() normalization \endlink \matrixworld</td><td>\code
355 scalar = vec1.norm(); scalar = vec1.squaredNorm()
356 vec2 = vec1.normalized(); vec1.normalize(); // inplace \endcode
360 \link MatrixBase::cross() cross product \endlink \matrixworld</td><td>\code
361 #include <Eigen/Geometry>
362 vec3 = vec1.cross(vec2);\endcode</td></tr>
365 <a href="#" class="top">top</a>
366 \section QuickRef_Coeffwise Coefficient-wise \& Array operators
368 In addition to the aforementioned operators, Eigen supports numerous coefficient-wise operator and functions.
369 Most of them unambiguously makes sense in array-world\arrayworld. The following operators are readily available for arrays,
370 or available through .array() for vectors and matrices:
372 <table class="manual">
373 <tr><td>Arithmetic operators</td><td>\code
374 array1 * array2 array1 / array2 array1 *= array2 array1 /= array2
375 array1 + scalar array1 - scalar array1 += scalar array1 -= scalar
377 <tr><td>Comparisons</td><td>\code
378 array1 < array2 array1 > array2 array1 < scalar array1 > scalar
379 array1 <= array2 array1 >= array2 array1 <= scalar array1 >= scalar
380 array1 == array2 array1 != array2 array1 == scalar array1 != scalar
381 array1.min(array2) array1.max(array2) array1.min(scalar) array1.max(scalar)
383 <tr><td>Trigo, power, and \n misc functions \n and the STL-like variants</td><td>\code
385 array1.abs() abs(array1)
386 array1.sqrt() sqrt(array1)
387 array1.log() log(array1)
388 array1.log10() log10(array1)
389 array1.exp() exp(array1)
390 array1.pow(array2) pow(array1,array2)
391 array1.pow(scalar) pow(array1,scalar)
397 array1.sin() sin(array1)
398 array1.cos() cos(array1)
399 array1.tan() tan(array1)
400 array1.asin() asin(array1)
401 array1.acos() acos(array1)
402 array1.atan() atan(array1)
403 array1.sinh() sinh(array1)
404 array1.cosh() cosh(array1)
405 array1.tanh() tanh(array1)
406 array1.arg() arg(array1)
408 array1.floor() floor(array1)
409 array1.ceil() ceil(array1)
410 array1.round() round(aray1)
412 array1.isFinite() isfinite(array1)
413 array1.isInf() isinf(array1)
414 array1.isNaN() isnan(array1)
420 The following coefficient-wise operators are available for all kind of expressions (matrices, vectors, and arrays), and for both real or complex scalar types:
422 <table class="manual">
423 <tr><th>Eigen's API</th><th>STL-like APIs\arrayworld </th><th>Comments</th></tr>
436 // read-write, no-op for real expressions
437 // read-only for real, read-write for complexes
438 // no-op for real expressions
443 Some coefficient-wise operators are readily available for for matrices and vectors through the following cwise* methods:
444 <table class="manual">
445 <tr><th>Matrix API \matrixworld</th><th>Via Array conversions</th></tr>
447 mat1.cwiseMin(mat2) mat1.cwiseMin(scalar)
448 mat1.cwiseMax(mat2) mat1.cwiseMax(scalar)
453 mat1.cwiseProduct(mat2)
454 mat1.cwiseQuotient(mat2)
455 mat1.cwiseEqual(mat2) mat1.cwiseEqual(scalar)
456 mat1.cwiseNotEqual(mat2)
459 mat1.array().min(mat2.array()) mat1.array().min(scalar)
460 mat1.array().max(mat2.array()) mat1.array().max(scalar)
464 mat1.array().inverse()
465 mat1.array() * mat2.array()
466 mat1.array() / mat2.array()
467 mat1.array() == mat2.array() mat1.array() == scalar
468 mat1.array() != mat2.array()
471 The main difference between the two API is that the one based on cwise* methods returns an expression in the matrix world,
472 while the second one (based on .array()) returns an array expression.
473 Recall that .array() has no cost, it only changes the available API and interpretation of the data.
475 It 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), <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):
477 mat1.unaryExpr(std::ptr_fun(foo));
478 mat1.unaryExpr(std::ref(foo));
479 mat1.unaryExpr([](double x) { return foo(x); });
483 <a href="#" class="top">top</a>
484 \section QuickRef_Reductions Reductions
486 Eigen provides several reduction methods such as:
487 \link DenseBase::minCoeff() minCoeff() \endlink, \link DenseBase::maxCoeff() maxCoeff() \endlink,
488 \link DenseBase::sum() sum() \endlink, \link DenseBase::prod() prod() \endlink,
489 \link MatrixBase::trace() trace() \endlink \matrixworld,
490 \link MatrixBase::norm() norm() \endlink \matrixworld, \link MatrixBase::squaredNorm() squaredNorm() \endlink \matrixworld,
491 \link DenseBase::all() all() \endlink, and \link DenseBase::any() any() \endlink.
492 All reduction operations can be done matrix-wise,
493 \link DenseBase::colwise() column-wise \endlink or
494 \link DenseBase::rowwise() row-wise \endlink. Usage example:
495 <table class="manual">
496 <tr><td rowspan="3" style="border-right-style:dashed;vertical-align:middle">\code
500 </td> <td>\code mat.minCoeff(); \endcode</td><td>\code 1 \endcode</td></tr>
501 <tr class="alt"><td>\code mat.colwise().minCoeff(); \endcode</td><td>\code 2 3 1 \endcode</td></tr>
502 <tr style="vertical-align:middle"><td>\code mat.rowwise().minCoeff(); \endcode</td><td>\code
509 Special versions of \link DenseBase::minCoeff(IndexType*,IndexType*) const minCoeff \endlink and \link DenseBase::maxCoeff(IndexType*,IndexType*) const maxCoeff \endlink:
512 s = vector.minCoeff(&i); // s == vector[i]
513 s = matrix.maxCoeff(&i, &j); // s == matrix(i,j)
515 Typical use cases of all() and any():
517 if((array1 > 0).all()) ... // if all coefficients of array1 are greater than 0 ...
518 if((array1 < array2).any()) ... // if there exist a pair i,j such that array1(i,j) < array2(i,j) ...
522 <a href="#" class="top">top</a>\section QuickRef_Blocks Sub-matrices
524 Read-write access to a \link DenseBase::col(Index) column \endlink
525 or a \link DenseBase::row(Index) row \endlink of a matrix (or array):
527 mat1.row(i) = mat2.col(j);
528 mat1.col(j1).swap(mat1.col(j2));
531 Read-write access to sub-vectors:
532 <table class="manual">
534 <th>Default versions</th>
535 <th>Optimized versions when the size \n is known at compile time</th></tr>
538 <tr><td>\code vec1.head(n)\endcode</td><td>\code vec1.head<n>()\endcode</td><td>the first \c n coeffs </td></tr>
539 <tr><td>\code vec1.tail(n)\endcode</td><td>\code vec1.tail<n>()\endcode</td><td>the last \c n coeffs </td></tr>
540 <tr><td>\code vec1.segment(pos,n)\endcode</td><td>\code vec1.segment<n>(pos)\endcode</td>
541 <td>the \c n coeffs in the \n range [\c pos : \c pos + \c n - 1]</td></tr>
542 <tr class="alt"><td colspan="3">
544 Read-write access to sub-matrices:</td></tr>
546 <td>\code mat1.block(i,j,rows,cols)\endcode
547 \link DenseBase::block(Index,Index,Index,Index) (more) \endlink</td>
548 <td>\code mat1.block<rows,cols>(i,j)\endcode
549 \link DenseBase::block(Index,Index) (more) \endlink</td>
550 <td>the \c rows x \c cols sub-matrix \n starting from position (\c i,\c j)</td></tr>
552 mat1.topLeftCorner(rows,cols)
553 mat1.topRightCorner(rows,cols)
554 mat1.bottomLeftCorner(rows,cols)
555 mat1.bottomRightCorner(rows,cols)\endcode
557 mat1.topLeftCorner<rows,cols>()
558 mat1.topRightCorner<rows,cols>()
559 mat1.bottomLeftCorner<rows,cols>()
560 mat1.bottomRightCorner<rows,cols>()\endcode
561 <td>the \c rows x \c cols sub-matrix \n taken in one of the four corners</td></tr>
564 mat1.bottomRows(rows)
566 mat1.rightCols(cols)\endcode
569 mat1.bottomRows<rows>()
570 mat1.leftCols<cols>()
571 mat1.rightCols<cols>()\endcode
572 <td>specialized versions of block() \n when the block fit two corners</td></tr>
577 <a href="#" class="top">top</a>\section QuickRef_Misc Miscellaneous operations
579 \subsection QuickRef_Reverse Reverse
580 Vectors, rows, and/or columns of a matrix can be reversed (see DenseBase::reverse(), DenseBase::reverseInPlace(), VectorwiseOp::reverse()).
582 vec.reverse() mat.colwise().reverse() mat.rowwise().reverse()
586 \subsection QuickRef_Replicate Replicate
587 Vectors, matrices, rows, and/or columns can be replicated in any direction (see DenseBase::replicate(), VectorwiseOp::replicate())
589 vec.replicate(times) vec.replicate<Times>
590 mat.replicate(vertical_times, horizontal_times) mat.replicate<VerticalTimes, HorizontalTimes>()
591 mat.colwise().replicate(vertical_times, horizontal_times) mat.colwise().replicate<VerticalTimes, HorizontalTimes>()
592 mat.rowwise().replicate(vertical_times, horizontal_times) mat.rowwise().replicate<VerticalTimes, HorizontalTimes>()
596 <a href="#" class="top">top</a>\section QuickRef_DiagTriSymm Diagonal, Triangular, and Self-adjoint matrices
597 (matrix world \matrixworld)
599 \subsection QuickRef_Diagonal Diagonal matrices
601 <table class="example">
602 <tr><th>Operation</th><th>Code</th></tr>
604 view a vector \link MatrixBase::asDiagonal() as a diagonal matrix \endlink \n </td><td>\code
605 mat1 = vec1.asDiagonal();\endcode
608 Declare a diagonal matrix</td><td>\code
609 DiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
610 diag1.diagonal() = vector;\endcode
612 <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>
614 vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal
615 vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal
616 vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal
617 vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal
618 vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal
622 <tr><td>Optimized products and inverse</td>
624 mat3 = scalar * diag1 * mat1;
625 mat3 += scalar * mat1 * vec1.asDiagonal();
626 mat3 = vec1.asDiagonal().inverse() * mat1
627 mat3 = mat1 * diag1.inverse()
633 \subsection QuickRef_TriangularView Triangular views
635 TriangularView 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.
637 \note The .triangularView() template member function requires the \c template keyword if it is used on an
638 object of a type that depends on a template parameter; see \ref TopicTemplateKeyword for details.
640 <table class="example">
641 <tr><th>Operation</th><th>Code</th></tr>
643 Reference to a triangular with optional \n
644 unit or null diagonal (read/write):
646 m.triangularView<Xxx>()
648 \c Xxx = ::Upper, ::Lower, ::StrictlyUpper, ::StrictlyLower, ::UnitUpper, ::UnitLower
651 Writing to a specific triangular part:\n (only the referenced triangular part is evaluated)
653 m1.triangularView<Eigen::Lower>() = m2 + m3 \endcode
656 Conversion to a dense matrix setting the opposite triangular part to zero:
658 m2 = m1.triangularView<Eigen::UnitUpper>()\endcode
663 m3 += s1 * m1.adjoint().triangularView<Eigen::UnitUpper>() * m2
664 m3 -= s1 * m2.conjugate() * m1.adjoint().triangularView<Eigen::Lower>() \endcode
667 Solving linear equations:\n
668 \f$ M_2 := L_1^{-1} M_2 \f$ \n
669 \f$ M_3 := {L_1^*}^{-1} M_3 \f$ \n
670 \f$ M_4 := M_4 U_1^{-1} \f$
672 L1.triangularView<Eigen::UnitLower>().solveInPlace(M2)
673 L1.triangularView<Eigen::Lower>().adjoint().solveInPlace(M3)
674 U1.triangularView<Eigen::Upper>().solveInPlace<OnTheRight>(M4)\endcode
678 \subsection QuickRef_SelfadjointMatrix Symmetric/selfadjoint views
680 Just as for triangular matrix, you can reference any triangular part of a square matrix to see it as a selfadjoint
681 matrix and perform special and optimized operations. Again the opposite triangular part is never referenced and can be
682 used to store other information.
684 \note The .selfadjointView() template member function requires the \c template keyword if it is used on an
685 object of a type that depends on a template parameter; see \ref TopicTemplateKeyword for details.
687 <table class="example">
688 <tr><th>Operation</th><th>Code</th></tr>
690 Conversion to a dense matrix:
692 m2 = m.selfadjointView<Eigen::Lower>();\endcode
695 Product with another general matrix or vector:
697 m3 = s1 * m1.conjugate().selfadjointView<Eigen::Upper>() * m3;
698 m3 -= s1 * m3.adjoint() * m1.selfadjointView<Eigen::Lower>();\endcode
701 Rank 1 and rank K update: \n
702 \f$ upper(M_1) \mathrel{{+}{=}} s_1 M_2 M_2^* \f$ \n
703 \f$ lower(M_1) \mathbin{{-}{=}} M_2^* M_2 \f$
705 M1.selfadjointView<Eigen::Upper>().rankUpdate(M2,s1);
706 M1.selfadjointView<Eigen::Lower>().rankUpdate(M2.adjoint(),-1); \endcode
709 Rank 2 update: (\f$ M \mathrel{{+}{=}} s u v^* + s v u^* \f$)
711 M.selfadjointView<Eigen::Upper>().rankUpdate(u,v,s);
715 Solving linear equations:\n(\f$ M_2 := M_1^{-1} M_2 \f$)
717 // via a standard Cholesky factorization
718 m2 = m1.selfadjointView<Eigen::Upper>().llt().solve(m2);
719 // via a Cholesky factorization with pivoting
720 m2 = m1.selfadjointView<Eigen::Lower>().ldlt().solve(m2);
728 <table class="tutorial_code">
730 \link MatrixBase::asDiagonal() make a diagonal matrix \endlink \n from a vector </td><td>\code
731 mat1 = vec1.asDiagonal();\endcode
734 Declare a diagonal matrix</td><td>\code
735 DiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
736 diag1.diagonal() = vector;\endcode
738 <tr><td>Access \link MatrixBase::diagonal() the diagonal and super/sub diagonals of a matrix \endlink as a vector (read/write)</td>
740 vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal
741 vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal
742 vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal
743 vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal
744 vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal
748 <tr><td>View on a triangular part of a matrix (read/write)</td>
750 mat2 = mat1.triangularView<Xxx>();
751 // Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower
752 mat1.triangularView<Upper>() = mat2 + mat3; // only the upper part is evaluated and referenced
755 <tr><td>View a triangular part as a symmetric/self-adjoint matrix (read/write)</td>
757 mat2 = mat1.selfadjointView<Xxx>(); // Xxx = Upper or Lower
758 mat1.selfadjointView<Upper>() = mat2 + mat2.adjoint(); // evaluated and write to the upper triangular part only
765 mat3 += scalar * vec1.asDiagonal() * mat1
766 mat3 += scalar * mat1 * vec1.asDiagonal()
767 mat3.noalias() += scalar * mat1.triangularView<Xxx>() * mat2
768 mat3.noalias() += scalar * mat2 * mat1.triangularView<Xxx>()
769 mat3.noalias() += scalar * mat1.selfadjointView<Upper or Lower>() * mat2
770 mat3.noalias() += scalar * mat2 * mat1.selfadjointView<Upper or Lower>()
771 mat1.selfadjointView<Upper or Lower>().rankUpdate(mat2);
772 mat1.selfadjointView<Upper or Lower>().rankUpdate(mat2.adjoint(), scalar);
775 Inverse products: (all are optimized)
777 mat3 = vec1.asDiagonal().inverse() * mat1
778 mat3 = mat1 * diag1.inverse()
779 mat1.triangularView<Xxx>().solveInPlace(mat2)
780 mat1.triangularView<Xxx>().solveInPlace<OnTheRight>(mat2)
781 mat2 = mat1.selfadjointView<Upper or Lower>().llt().solve(mat2)