eigen/doc/QuickReference.dox
2010-06-26 16:59:18 +02:00

649 lines
21 KiB
Plaintext

namespace Eigen {
/** \page QuickRefPage Quick reference guide
\b Table \b of \b contents
- \ref QuickRef_Headers
- \ref QuickRef_Types
- \ref QuickRef_Map
- \ref QuickRef_ArithmeticOperators
- \ref QuickRef_Coeffwise
\n
<hr>
<a href="#" class="top">top</a>
\section QuickRef_Headers Modules and Header files
<table class="tutorial_code">
<tr><td>Module</td><td>Header file</td><td>Contents</td></tr>
<tr><td>Core</td><td>\code#include <Eigen/Core>\endcode</td><td>Matrix and Array classes, basic linear algebra (including triangular and selfadjoint products), array manipulation</td></tr>
<tr><td>Geometry</td><td>\code#include <Eigen/Geometry>\endcode</td><td>Transformation, Translation, Scaling, 2D and 3D rotations (Quaternion, AngleAxis)</td></tr>
<tr><td>LU</td><td>\code#include <Eigen/LU>\endcode</td><td>Inverse, determinant, LU decompositions (FullPivLU, PartialPivLU) with solver</td></tr>
<tr><td>Cholesky</td><td>\code#include <Eigen/Cholesky>\endcode</td><td>LLT and LDLT Cholesky factorization with solver</td></tr>
<tr><td>SVD</td><td>\code#include <Eigen/SVD>\endcode</td><td>SVD decomposition with solver (HouseholderSVD, JacobiSVD)</td></tr>
<tr><td>QR</td><td>\code#include <Eigen/QR>\endcode</td><td>QR decomposition with solver (HouseholderQR, ColPivHouseholerQR, FullPivHouseholderQR)</td></tr>
<tr><td>Eigenvalues</td><td>\code#include <Eigen/Eigenvalues>\endcode</td><td>Eigenvalue, eigenvector decompositions for selfadjoint and non selfadjoint real or complex matrices.</td></tr>
<tr><td>Sparse</td><td>\code#include <Eigen/Sparse>\endcode</td><td>Sparse matrix storage and related basic linear algebra.</td></tr>
<tr><td></td><td>\code#include <Eigen/Dense>\endcode</td><td>Includes Core, Geometry, LU, Cholesky, SVD, QR, and Eigenvalues</td></tr>
<tr><td></td><td>\code#include <Eigen/Eigen>\endcode</td><td>Includes Dense and Sparse</td></tr>
</table>
<a href="#" class="top">top</a>
\section QuickRef_Types Array, matrix and vector types
\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:
\code
typedef Matrix<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyMatrixType;
typedef Array<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyArrayType;
\endcode
\li \c Scalar is the scalar type of the coefficients (e.g., \c float, \c double, \c bool, \c int, etc.).
\li \c RowsAtCompileTime and \c ColsAtCompileTime are the number of rows and columns of the matrix as known at compile-time or \c Dynamic.
\li \c Options can be \c ColMajor or \c RowMajor, default is \c ColMajor. (see class Matrix for more options)
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:
\code
Matrix<double, 6, Dynamic> // Dynamic number of columns (heap allocation)
Matrix<double, Dynamic, 2> // Dynamic number of rows (heap allocation)
Matrix<double, Dynamic, Dynamic, RowMajor> // Fully dynamic, row major (heap allocation)
Matrix<double, 13, 3> // Fully fixed (static allocation)
\endcode
In most cases, you can simply use one of the convenience typedefs for \ref matrixtypedefs "matrices" and \ref arraytypedefs "arrays". Some examples:
<table class="tutorial_code">
<tr><td>\code
Matrix<float,Dynamic,Dynamic> <=> MatrixXf
Matrix<double,Dynamic,1> <=> VectorXd
Matrix<int,1,Dynamic> <=> RowVectorXi
Matrix<float,3,3> <=> Matrix3f
Matrix<float,4,1> <=> Vector4f
\endcode</td><td>\code
Array<float,Dynamic,Dynamic> <=> ArrayXXf
Array<double,Dynamic,1> <=> ArrayXd
Array<int,1,Dynamic> <=> RowArrayXi
Array<float,3,3> <=> Array33f
Array<float,4,1> <=> Array4f
\endcode</td></tr>
</table>
Conversion between the matrix and array worlds:
\code
Array44f a1, a1;
Matrix4f m1, m2;
m1 = a1 * a2; // OK, coeffwise product
a1 = m1 * m2; // OK, matrix product
a2 = a1 + m1.array();
m2 = a1.matrix() + m1;
ArrayWrapper<Matrix4f> m1a(m1); // m1a is an alias for m1.array(), they share the same coefficients
MatrixWrapper<Array44f> a1m(a1);
\endcode
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:
\li <a name="matrixonly"><a/>\matrixworld linear algebra matrix and vector only
\li <a name="arrayonly"><a/>\arrayworld array objects only
\subsection QuickRef_Basics Basic matrix manipulation
<table class="tutorial_code">
<tr><td></td><td>1D objects</td><td>2D objects</td><td>Notes</td></tr>
<tr><td>Constructors</td>
<td>\code
Vector4d v4;
Vector2f v1(x, y);
Array3i v2(x, y, z);
Vector4d v3(x, y, z, w);
VectorXf v5;
ArrayXf v6(size);
\endcode</td><td>\code
Matrix4f m1;
MatrixXf m5;
MatrixXf m6(nb_rows, nb_columns);
\endcode</td><td>
The coeffs are left uninitialized \n \n
\n \n Empty object \n The coeffs are left uninitialized</td></tr>
<tr><td>Comma initializer</td>
<td>\code
Vector3f v1; v1 << x, y, z;
ArrayXf v2(4); v2 << 1, 2, 3, 4;
\endcode</td><td>\code
Matrix3f m1; m1 << 1, 2, 3,
4, 5, 6,
7, 8, 9;
\endcode</td><td></td></tr>
<tr><td>Comma initializer (bis)</td>
<td colspan="2">
<table class="noborder"><tr><td>
\include Tutorial_commainit_02.cpp
</td>
<td>
output:
\verbinclude Tutorial_commainit_02.out
</td></tr></table>
<td></td></tr>
<tr><td>Runtime info</td>
<td>\code
vector.size();
vector.innerStride();
vector.data();
\endcode</td><td>\code
matrix.rows(); matrix.cols();
matrix.innerSize(); matrix.outerSize();
matrix.innerStride(); matrix.outerStride();
matrix.data();
\endcode</td><td>\n Inner/Outer* are storage order dependent</td></tr>
<tr><td>Compile-time info</td>
<td colspan="2">\code
ObjectType::Scalar ObjectType::RowsAtCompileTime
ObjectType::RealScalar ObjectType::ColsAtCompileTime
ObjectType::Index ObjectType::SizeAtCompileTime
\endcode</td><td></td></tr>
<tr><td>Resizing</td>
<td>\code
vector.resize(size);
vector.resizeLike(other_vector);
vector.conservativeResize(size);
\endcode</td><td>\code
matrix.resize(nb_rows, nb_cols);
matrix.resize(Eigen::NoChange, nb_cols);
matrix.resize(nb_rows, Eigen::NoChange);
matrix.resizeLike(other_matrix);
matrix.conservativeResize(nb_rows, nb_cols);
\endcode</td><td></td>no-op if the new sizes match,\n otherwise data are lost \n \n resizing with data preservation</tr>
<tr><td>Coeff access with \n range checking</td>
<td>\code
vector(i) vector.x()
vector[i] vector.y()
vector.z()
vector.w()
\endcode</td><td>\code
matrix(i,j)
\endcode</td><td><em class="note">Range checking is disabled if \n NDEBUG or EIGEN_NO_DEBUG is defined</em></td></tr>
<tr><td>Coeff access without \n range checking</td>
<td>\code
vector.coeff(i)
vector.coeffRef(i)
\endcode</td><td>\code
matrix.coeff(i,j)
matrix.coeffRef(i,j)
\endcode</td><td></td></tr>
<tr><td>Assignment/copy</td>
<td colspan="2">\code
object = expression;
object_of_float = expression_of_double.cast<float>();
\endcode</td><td>the destination is automatically resized (if possible)</td></tr>
</table>
\subsection QuickRef_PredefMat Predefined Matrices
<table class="tutorial_code">
<tr>
<td>Fixed-size matrix or vector</td>
<td>Dynamic-size matrix</td>
<td>Dynamic-size vector</td>
</tr>
<tr style="border-bottom-style: none;">
<td>
\code
typedef {Matrix3f|Array33f} FixedXD;
FixedXD x;
x = FixedXD::Zero();
x = FixedXD::Ones();
x = FixedXD::Constant(value);
x = FixedXD::Random();
x.setZero();
x.setOnes();
x.setConstant(value);
x.setRandom();
\endcode
</td>
<td>
\code
typedef {MatrixXf|ArrayXXf} Dynamic2D;
Dynamic2D x;
x = Dynamic2D::Zero(rows, cols);
x = Dynamic2D::Ones(rows, cols);
x = Dynamic2D::Constant(rows, cols, value);
x = Dynamic2D::Random(rows, cols);
x.setZero(rows, cols);
x.setOnes(rows, cols);
x.setConstant(rows, cols, value);
x.setRandom(rows, cols);
\endcode
</td>
<td>
\code
typedef {VectorXf|ArrayXf} Dynamic1D;
Dynamic1D x;
x = Dynamic1D::Zero(size);
x = Dynamic1D::Ones(size);
x = Dynamic1D::Constant(size, value);
x = Dynamic1D::Random(size);
x.setZero(size);
x.setOnes(size);
x.setConstant(size, value);
x.setRandom(size);
\endcode
</td>
</tr>
<tr><td colspan="3">Identity and \link MatrixBase::Unit basis vectors \endlink \matrixworld</td></tr>
<tr style="border-bottom-style: none;">
<td>
\code
x = FixedXD::Identity();
x.setIdentity();
Vector3f::UnitX() // 1 0 0
Vector3f::UnitY() // 0 1 0
Vector3f::UnitZ() // 0 0 1
\endcode
</td>
<td>
\code
x = Dynamic2D::Identity(rows, cols);
x.setIdentity(rows, cols);
N/A
\endcode
</td>
<td>\code
N/A
VectorXf::Unit(size,i)
VectorXf::Unit(4,1) == Vector4f(0,1,0,0)
== Vector4f::UnitY()
\endcode
</td>
</tr>
</table>
\subsection QuickRef_Map Map
<table class="tutorial_code">
<tr>
<td>Contiguous memory</td>
<td>\code
float data[] = {1,2,3,4};
Map<Vector3f> v1(data); // uses v1 as a Vector3f object
Map<ArrayXf> v2(data,3); // uses v2 as a ArrayXf object
Map<Array22f> m1(data); // uses m1 as a Array22f object
Map<MatrixXf> m2(data,2,2); // uses m2 as a MatrixXf object
\endcode</td>
</tr>
<tr>
<td>Typical usage of strides</td>
<td>\code
float data[] = {1,2,3,4,5,6,7,8,9};
Map<VectorXf,0,InnerStride<2> > v1(data,3); // == [1,3,5]
Map<VectorXf,0,InnerStride<> > v2(data,3,InnerStride<>(3)); // == [1,4,7]
Map<MatrixXf,0,OuterStride<> > m1(data,2,3,OuterStride<>(3)); // == |1,4,7|
Map<MatrixXf,0,OuterStride<3> > m2(data,2,3); // |2,5,8|
\endcode</td>
</tr>
</table>
<a href="#" class="top">top</a>
\section QuickRef_ArithmeticOperators Arithmetic Operators
<table class="tutorial_code">
<tr><td>
add/subtract</td><td>\code
mat3 = mat1 + mat2; mat3 += mat1;
mat3 = mat1 - mat2; mat3 -= mat1;\endcode
</td></tr>
<tr><td>
scalar product</td><td>\code
mat3 = mat1 * s1; mat3 = s1 * mat1; mat3 *= s1;
mat3 = mat1 / s1; mat3 /= s1;\endcode
</td></tr>
<tr><td>
matrix/vector product \matrixworld</td><td>\code
col2 = mat1 * col1;
row2 = row1 * mat1; row1 *= mat1;
mat3 = mat1 * mat2; mat3 *= mat1; \endcode
</td></tr>
<tr><td>
\link MatrixBase::dot() dot \endlink \& inner products \matrixworld</td><td>\code
scalar = col1.adjoint() * col2;
scalar = (col1.adjoint() * col2).value();
scalar = vec1.dot(vec2);\endcode
</td></tr>
<tr><td>
outer product \matrixworld</td><td>\code
mat = col1 * col2.transpose();\endcode
</td></tr>
<tr><td>
\link MatrixBase::cross() cross product \endlink \matrixworld</td><td>\code
#include <Eigen/Geometry>
vec3 = vec1.cross(vec2);\endcode</td></tr>
</table>
<a href="#" class="top">top</a>
\section QuickRef_Coeffwise Coefficient-wise \& Array operators
Coefficient-wise operators for matrices and vectors:
<table class="tutorial_code">
<tr><td>Matrix API \matrixworld</td><td>Via Array conversions</td></tr>
<tr><td>\code
mat1.cwiseMin(mat2)
mat1.cwiseMax(mat2)
mat1.cwiseAbs2()
mat1.cwiseAbs()
mat1.cwiseSqrt()
mat1.cwiseProduct(mat2)
mat1.cwiseQuotient(mat2)\endcode
</td><td>\code
mat1.array().min(mat2.array())
mat1.array().max(mat2.array())
mat1.array().abs2()
mat1.array().abs()
mat1.array().sqrt()
mat1.array() * mat2.array()
mat1.array() / mat2.array()
\endcode</td></tr>
</table>
Array operators:\arrayworld
<table class="tutorial_code">
<tr><td>Arithmetic operators</td><td>\code
array1 * array2 array1 / array2 array1 *= array2 array1 /= array2
array1 + scalar array1 - scalar array1 += scalar array1 -= scalar
\endcode</td></tr>
<tr><td>Comparisons</td><td>\code
array1 < array2 array1 > array2 array1 < scalar array1 > scalar
array1 <= array2 array1 >= array2 array1 <= scalar array1 >= scalar
array1 == array2 array1 != array2 array1 == scalar array1 != scalar
\endcode</td></tr>
<tr><td>Special functions \n and STL variants</td><td>\code
array1.min(array2) std::min(array1,array2)
array1.max(array2) std::max(array1,array2)
array1.abs2()
array1.abs() std::abs(array1)
array1.sqrt() std::sqrt(array1)
array1.log() std::log(array1)
array1.exp() std::exp(array1)
array1.pow(exponent) std::pow(array1,exponent)
array1.square()
array1.cube()
array1.inverse()
array1.sin() std::sin(array1)
array1.cos() std::cos(array1)
array1.tan() std::tan(array1)
\endcode
</td></tr>
</table>
*/
// FIXME I stopped here
/**
<a href="#" class="top">top</a>
\section TutorialCoreReductions Reductions
Eigen provides several reduction methods such as:
\link DenseBase::minCoeff() minCoeff() \endlink, \link DenseBase::maxCoeff() maxCoeff() \endlink,
\link DenseBase::sum() sum() \endlink, \link MatrixBase::trace() trace() \endlink \matrixworld,
\link MatrixBase::norm() norm() \endlink \matrixworld, \link MatrixBase::squaredNorm() squaredNorm() \endlink \matrixworld,
\link DenseBase::all() all() \endlink \redstar,and \link DenseBase::any() any() \endlink \redstar.
All reduction operations can be done matrix-wise,
\link DenseBase::colwise() column-wise \endlink \redstar or
\link DenseBase::rowwise() row-wise \endlink \redstar. Usage example:
<table class="tutorial_code">
<tr><td rowspan="3" style="border-right-style:dashed">\code
5 3 1
mat = 2 7 8
9 4 6 \endcode
</td> <td>\code mat.minCoeff(); \endcode</td><td>\code 1 \endcode</td></tr>
<tr><td>\code mat.colwise().minCoeff(); \endcode</td><td>\code 2 3 1 \endcode</td></tr>
<tr><td>\code mat.rowwise().minCoeff(); \endcode</td><td>\code
1
2
4
\endcode</td></tr>
</table>
Also note that maxCoeff and minCoeff can takes optional arguments returning the coordinates of the respective min/max coeff: \link DenseBase::maxCoeff(int*,int*) const maxCoeff(int* i, int* j) \endlink, \link DenseBase::minCoeff(int*,int*) const minCoeff(int* i, int* j) \endlink.
<span class="note">\b Side \b note: The all() and any() functions are especially useful in combination with coeff-wise comparison operators.</span>
<a href="#" class="top">top</a>\section TutorialCoreMatrixBlocks Matrix blocks
Read-write access to a \link DenseBase::col(int) column \endlink
or a \link DenseBase::row(int) row \endlink of a matrix (or array):
\code
mat1.row(i) = mat2.col(j);
mat1.col(j1).swap(mat1.col(j2));
\endcode
Read-write access to sub-vectors:
<table class="tutorial_code">
<tr>
<td>Default versions</td>
<td>Optimized versions when the size \n is known at compile time</td></tr>
<td></td>
<tr><td>\code vec1.head(n)\endcode</td><td>\code vec1.head<n>()\endcode</td><td>the first \c n coeffs </td></tr>
<tr><td>\code vec1.tail(n)\endcode</td><td>\code vec1.tail<n>()\endcode</td><td>the last \c n coeffs </td></tr>
<tr><td>\code vec1.segment(pos,n)\endcode</td><td>\code vec1.segment<n>(pos)\endcode</td>
<td>the \c size coeffs in \n the range [\c pos : \c pos + \c n [</td></tr>
<tr style="border-style: dashed none dashed none;"><td>
Read-write access to sub-matrices:</td><td></td><td></td></tr>
<tr>
<td>\code mat1.block(i,j,rows,cols)\endcode
\link DenseBase::block(int,int,int,int) (more) \endlink</td>
<td>\code mat1.block<rows,cols>(i,j)\endcode
\link DenseBase::block(int,int) (more) \endlink</td>
<td>the \c rows x \c cols sub-matrix \n starting from position (\c i,\c j)</td></tr><tr>
<td>\code
mat1.topLeftCorner(rows,cols)
mat1.topRightCorner(rows,cols)
mat1.bottomLeftCorner(rows,cols)
mat1.bottomRightCorner(rows,cols)\endcode
<td>\code
mat1.topLeftCorner<rows,cols>()
mat1.topRightCorner<rows,cols>()
mat1.bottomLeftCorner<rows,cols>()
mat1.bottomRightCorner<rows,cols>()\endcode
<td>the \c rows x \c cols sub-matrix \n taken in one of the four corners</td></tr>
</table>
<a href="#" class="top">top</a>\section TutorialCoreDiagonalMatrices Diagonal matrices
\matrixworld
<table class="tutorial_code">
<tr><td>
\link MatrixBase::asDiagonal() make a diagonal matrix \endlink from a vector \n
<em class="note">this product is automatically optimized !</em></td><td>\code
mat3 = mat1 * vec2.asDiagonal();\endcode
</td></tr>
<tr><td>Access \link MatrixBase::diagonal() the diagonal of a matrix \endlink as a vector (read/write)</td>
<td>\code
vec1 = mat1.diagonal();
mat1.diagonal() = vec1;
\endcode
</td>
</tr>
</table>
<a href="#" class="top">top</a>
\section TutorialCoreTransposeAdjoint Transpose and Adjoint operations
<table class="tutorial_code">
<tr><td>
\link DenseBase::transpose() transposition \endlink (read-write)</td><td>\code
mat3 = mat1.transpose() * mat2;
mat3.transpose() = mat1 * mat2.transpose();
\endcode
</td></tr>
<tr><td>
\link MatrixBase::adjoint() adjoint \endlink (read only) \matrixworld\n</td><td>\code
mat3 = mat1.adjoint() * mat2;
\endcode
</td></tr>
</table>
<a href="#" class="top">top</a>
\section TutorialCoreDotNorm Dot-product, vector norm, normalization \matrixworld
<table class="tutorial_code">
<tr><td>
\link MatrixBase::dot() Dot-product \endlink of two vectors
</td><td>\code vec1.dot(vec2);\endcode
</td></tr>
<tr><td>
\link MatrixBase::norm() norm \endlink of a vector \n
\link MatrixBase::squaredNorm() squared norm \endlink of a vector
</td><td>\code vec.norm(); \endcode \n \code vec.squaredNorm() \endcode
</td></tr>
<tr><td>
returns a \link MatrixBase::normalized() normalized \endlink vector \n
\link MatrixBase::normalize() normalize \endlink a vector
</td><td>\code
vec3 = vec1.normalized();
vec1.normalize();\endcode
</td></tr>
</table>
<a href="#" class="top">top</a>
\section TutorialCoreTriangularMatrix Dealing with triangular matrices \matrixworld
Currently, Eigen does not provide any explicit triangular matrix, with storage class. Instead, we
can reference a triangular part of a square matrix or expression to perform special treatment on it.
This is achieved by the class TriangularView and the MatrixBase::triangularView template function.
Note that the opposite triangular part of the matrix is never referenced, and so it can, e.g., store
a second triangular matrix.
<table class="tutorial_code">
<tr><td>
Reference a read/write triangular part of a given \n
matrix (or expression) m with optional unit diagonal:
</td><td>\code
m.triangularView<Eigen::UpperTriangular>()
m.triangularView<Eigen::UnitUpperTriangular>()
m.triangularView<Eigen::LowerTriangular>()
m.triangularView<Eigen::UnitLowerTriangular>()\endcode
</td></tr>
<tr><td>
Writing to a specific triangular part:\n (only the referenced triangular part is evaluated)
</td><td>\code
m1.triangularView<Eigen::LowerTriangular>() = m2 + m3 \endcode
</td></tr>
<tr><td>
Conversion to a dense matrix setting the opposite triangular part to zero:
</td><td>\code
m2 = m1.triangularView<Eigen::UnitUpperTriangular>()\endcode
</td></tr>
<tr><td>
Products:
</td><td>\code
m3 += s1 * m1.adjoint().triangularView<Eigen::UnitUpperTriangular>() * m2
m3 -= s1 * m2.conjugate() * m1.adjoint().triangularView<Eigen::LowerTriangular>() \endcode
</td></tr>
<tr><td>
Solving linear equations:\n(\f$ m_2 := m_1^{-1} m_2 \f$)
</td><td>\code
m1.triangularView<Eigen::UnitLowerTriangular>().solveInPlace(m2)
m1.adjoint().triangularView<Eigen::UpperTriangular>().solveInPlace(m2)\endcode
</td></tr>
</table>
<a href="#" class="top">top</a>
\section TutorialCoreSelfadjointMatrix Dealing with symmetric/selfadjoint matrices \matrixworld
Just as for triangular matrix, you can reference any triangular part of a square matrix to see it a selfadjoint
matrix to perform special and optimized operations. Again the opposite triangular is never referenced and can be
used to store other information.
<table class="tutorial_code">
<tr><td>
Conversion to a dense matrix:
</td><td>\code
m2 = m.selfadjointView<Eigen::LowerTriangular>();\endcode
</td></tr>
<tr><td>
Product with another general matrix or vector:
</td><td>\code
m3 = s1 * m1.conjugate().selfadjointView<Eigen::UpperTriangular>() * m3;
m3 -= s1 * m3.adjoint() * m1.selfadjointView<Eigen::UpperTriangular>();\endcode
</td></tr>
<tr><td>
Rank 1 and rank K update:
</td><td>\code
// fast version of m1 += s1 * m2 * m2.adjoint():
m1.selfadjointView<Eigen::UpperTriangular>().rankUpdate(m2,s1);
// fast version of m1 -= m2.adjoint() * m2:
m1.selfadjointView<Eigen::LowerTriangular>().rankUpdate(m2.adjoint(),-1); \endcode
</td></tr>
<tr><td>
Rank 2 update: (\f$ m += s u v^* + s v u^* \f$)
</td><td>\code
m.selfadjointView<Eigen::UpperTriangular>().rankUpdate(u,v,s);
\endcode
</td></tr>
<tr><td>
Solving linear equations:\n(\f$ m_2 := m_1^{-1} m_2 \f$)
</td><td>\code
// via a standard Cholesky factorization
m1.selfadjointView<Eigen::UpperTriangular>().llt().solveInPlace(m2);
// via a Cholesky factorization with pivoting
m1.selfadjointView<Eigen::UpperTriangular>().ldlt().solveInPlace(m2);
\endcode
</td></tr>
</table>
<a href="#" class="top">top</a>
\section TutorialCoreSpecialTopics Special Topics
\ref TopicLazyEvaluation "Lazy Evaluation and Aliasing": Thanks to expression templates, Eigen is able to apply lazy evaluation wherever that is beneficial.
*/
}