\section QuickRef_Headers Modules and Header files
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.
<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>
<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>
<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>
<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>
<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>
<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>
<tr ><td></td><td>\code#include <Eigen/Eigen>\endcode</td><td>Includes %Dense and %Sparse header files (the whole Eigen library)</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:
<tr><td>Trigo, power, and \n misc functions \n and the STL-like variants</td><td>\code
array1.abs2()
array1.abs() abs(array1)
array1.sqrt() sqrt(array1)
array1.log() log(array1)
array1.log10() log10(array1)
array1.exp() exp(array1)
array1.pow(array2) pow(array1,array2)
array1.pow(scalar) pow(array1,scalar)
pow(scalar,array2)
array1.square()
array1.cube()
array1.inverse()
array1.sin() sin(array1)
array1.cos() cos(array1)
array1.tan() tan(array1)
array1.asin() asin(array1)
array1.acos() acos(array1)
array1.atan() atan(array1)
array1.sinh() sinh(array1)
array1.cosh() cosh(array1)
array1.tanh() tanh(array1)
array1.arg() arg(array1)
array1.floor() floor(array1)
array1.ceil() ceil(array1)
array1.round() round(aray1)
array1.isFinite() isfinite(array1)
array1.isInf() isinf(array1)
array1.isNaN() isnan(array1)
\endcode
</td></tr>
</table>
The following coefficient-wise operators are available for all kind of expressions (matrices, vectors, and arrays), and for both real or complex scalar types:
The main difference between the two API is that the one based on cwise* methods returns an expression in the matrix world,
while the second one (based on .array()) returns an array expression.
Recall that .array() has no cost, it only changes the available API and interpretation of the data.
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, 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):
\code
mat1.unaryExpr(std::ptr_fun(foo));
mat1.unaryExpr(std::ref(foo));
mat1.unaryExpr([](double x) { return foo(x); });
\endcode
Please note that it's not possible to pass a raw function pointer to \c unaryExpr, so please warp it as shown above.
Special versions of \link DenseBase::minCoeff(IndexType*,IndexType*) const minCoeff \endlink and \link DenseBase::maxCoeff(IndexType*,IndexType*) const maxCoeff \endlink:
\code
int i, j;
s = vector.minCoeff(&i); // s == vector[i]
s = matrix.maxCoeff(&i, &j); // s == matrix(i,j)
\endcode
Typical use cases of all() and any():
\code
if((array1 > 0).all()) ... // if all coefficients of array1 are greater than 0 ...
if((array1 < array2).any()) ... // if there exist a pair i,j such that array1(i,j) < array2(i,j) ...
<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>
<td>\code
vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal
vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal
vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal
vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal
vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal
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.
\note The .triangularView() template member function requires the \c template keyword if it is used on an
object of a type that depends on a template parameter; see \ref TopicTemplateKeyword for details.