| 1D objects | 2D objects | Notes |
Constructors |
\code
Vector4d v4;
Vector2f v1(x, y);
Array3i v2(x, y, z);
Vector4d v3(x, y, z, w);
VectorXf v5; // empty object
ArrayXf v6(size);
\endcode | \code
Matrix4f m1;
MatrixXf m5; // empty object
MatrixXf m6(nb_rows, nb_columns);
\endcode |
By default, the coefficients \n are left uninitialized |
Comma initializer |
\code
Vector3f v1; v1 << x, y, z;
ArrayXf v2(4); v2 << 1, 2, 3, 4;
\endcode | \code
Matrix3f m1; m1 << 1, 2, 3,
4, 5, 6,
7, 8, 9;
\endcode | |
Comma initializer (bis) |
\include Tutorial_commainit_02.cpp
|
output:
\verbinclude Tutorial_commainit_02.out
|
Runtime info |
\code
vector.size();
vector.innerStride();
vector.data();
\endcode | \code
matrix.rows(); matrix.cols();
matrix.innerSize(); matrix.outerSize();
matrix.innerStride(); matrix.outerStride();
matrix.data();
\endcode | Inner/Outer* are storage order dependent |
Compile-time info |
\code
ObjectType::Scalar ObjectType::RowsAtCompileTime
ObjectType::RealScalar ObjectType::ColsAtCompileTime
ObjectType::Index ObjectType::SizeAtCompileTime
\endcode | |
Resizing |
\code
vector.resize(size);
vector.resizeLike(other_vector);
vector.conservativeResize(size);
\endcode | \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 | no-op if the new sizes match, otherwise data are lost
resizing with data preservation |
Coeff access with \n range checking |
\code
vector(i) vector.x()
vector[i] vector.y()
vector.z()
vector.w()
\endcode | \code
matrix(i,j)
\endcode | Range checking is disabled if \n NDEBUG or #EIGEN_NO_DEBUG is defined |
Coeff access without \n range checking |
\code
vector.coeff(i)
vector.coeffRef(i)
\endcode | \code
matrix.coeff(i,j)
matrix.coeffRef(i,j)
\endcode | |
Assignment/copy |
\code
object = expression;
object_of_float = expression_of_double.cast();
\endcode | the destination is automatically resized (if possible) |
\subsection QuickRef_PredefMat Predefined Matrices
Fixed-size matrix or vector |
Dynamic-size matrix |
Dynamic-size vector |
\code
typedef {Matrix3f|Array33f} FixedXD;
FixedXD x;
x = FixedXD::Zero();
x = FixedXD::Ones();
x = FixedXD::Constant(value);
x = FixedXD::Random();
x = FixedXD::LinSpaced(low, high, size);
x.setZero();
x.setOnes();
x.setConstant(value);
x.setRandom();
x.setLinSpaced(low, high, size);
\endcode
|
\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);
N/A
x.setZero(rows, cols);
x.setOnes(rows, cols);
x.setConstant(rows, cols, value);
x.setRandom(rows, cols);
N/A
\endcode
|
\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 = Dynamic1D::LinSpaced(low, high, size);
x.setZero(size);
x.setOnes(size);
x.setConstant(size, value);
x.setRandom(size);
x.setLinSpaced(low, high, size);
\endcode
|
Identity and \link MatrixBase::Unit basis vectors \endlink \matrixworld |
\code
x = FixedXD::Identity();
x.setIdentity();
Vector3f::UnitX() // 1 0 0
Vector3f::UnitY() // 0 1 0
Vector3f::UnitZ() // 0 0 1
\endcode
|
\code
x = Dynamic2D::Identity(rows, cols);
x.setIdentity(rows, cols);
N/A
\endcode
|
\code
N/A
VectorXf::Unit(size,i)
VectorXf::Unit(4,1) == Vector4f(0,1,0,0)
== Vector4f::UnitY()
\endcode
|
\subsection QuickRef_Map Mapping external arrays
Contiguous \n memory |
\code
float data[] = {1,2,3,4};
Map v1(data); // uses v1 as a Vector3f object
Map v2(data,3); // uses v2 as a ArrayXf object
Map m1(data); // uses m1 as a Array22f object
Map m2(data,2,2); // uses m2 as a MatrixXf object
\endcode |
Typical usage \n of strides |
\code
float data[] = {1,2,3,4,5,6,7,8,9};
Map > v1(data,3); // = [1,3,5]
Map > v2(data,3,InnerStride<>(3)); // = [1,4,7]
Map > m2(data,2,3); // both lines |1,4,7|
Map > m1(data,2,3,OuterStride<>(3)); // are equal to: |2,5,8|
\endcode |
add \n subtract | \code
mat3 = mat1 + mat2; mat3 += mat1;
mat3 = mat1 - mat2; mat3 -= mat1;\endcode
|
scalar product | \code
mat3 = mat1 * s1; mat3 *= s1; mat3 = s1 * mat1;
mat3 = mat1 / s1; mat3 /= s1;\endcode
|
matrix/vector \n products \matrixworld | \code
col2 = mat1 * col1;
row2 = row1 * mat1; row1 *= mat1;
mat3 = mat1 * mat2; mat3 *= mat1; \endcode
|
transposition \n adjoint \matrixworld | \code
mat1 = mat2.transpose(); mat1.transposeInPlace();
mat1 = mat2.adjoint(); mat1.adjointInPlace();
\endcode
|
\link MatrixBase::dot() dot \endlink product \n inner product \matrixworld | \code
scalar = vec1.dot(vec2);
scalar = col1.adjoint() * col2;
scalar = (col1.adjoint() * col2).value();\endcode
|
outer product \matrixworld | \code
mat = col1 * col2.transpose();\endcode
|
\link MatrixBase::norm() norm \endlink \n \link MatrixBase::normalized() normalization \endlink \matrixworld | \code
scalar = vec1.norm(); scalar = vec1.squaredNorm()
vec2 = vec1.normalized(); vec1.normalize(); // inplace \endcode
|
\link MatrixBase::cross() cross product \endlink \matrixworld | \code
#include
vec3 = vec1.cross(vec2);\endcode |
Arithmetic operators | \code
array1 * array2 array1 / array2 array1 *= array2 array1 /= array2
array1 + scalar array1 - scalar array1 += scalar array1 -= scalar
\endcode |
Comparisons | \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 |
Trigo, power, and \n misc functions \n and the STL variants | \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
|
Default versions |
Optimized versions when the size \n is known at compile time |
|
\code vec1.head(n)\endcode | \code vec1.head()\endcode | the first \c n coeffs |
\code vec1.tail(n)\endcode | \code vec1.tail()\endcode | the last \c n coeffs |
\code vec1.segment(pos,n)\endcode | \code vec1.segment(pos)\endcode |
the \c n coeffs in \n the range [\c pos : \c pos + \c n [ |
Read-write access to sub-matrices: |
\code mat1.block(i,j,rows,cols)\endcode
\link DenseBase::block(Index,Index,Index,Index) (more) \endlink |
\code mat1.block(i,j)\endcode
\link DenseBase::block(Index,Index) (more) \endlink |
the \c rows x \c cols sub-matrix \n starting from position (\c i,\c j) |
\code
mat1.topLeftCorner(rows,cols)
mat1.topRightCorner(rows,cols)
mat1.bottomLeftCorner(rows,cols)
mat1.bottomRightCorner(rows,cols)\endcode
| \code
mat1.topLeftCorner()
mat1.topRightCorner()
mat1.bottomLeftCorner()
mat1.bottomRightCorner()\endcode
the \c rows x \c cols sub-matrix \n taken in one of the four corners | |
\code
mat1.topRows(rows)
mat1.bottomRows(rows)
mat1.leftCols(cols)
mat1.rightCols(cols)\endcode
| \code
mat1.topRows()
mat1.bottomRows()
mat1.leftCols()
mat1.rightCols()\endcode
specialized versions of block() \n when the block fit two corners | |
Operation | Code |
view a vector \link MatrixBase::asDiagonal() as a diagonal matrix \endlink \n | \code
mat1 = vec1.asDiagonal();\endcode
|
Declare a diagonal matrix | \code
DiagonalMatrix diag1(size);
diag1.diagonal() = vector;\endcode
|
Access the \link MatrixBase::diagonal() diagonal \endlink and \link MatrixBase::diagonal(Index) super/sub diagonals \endlink of a matrix as a vector (read/write) |
\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
\endcode |
Optimized products and inverse |
\code
mat3 = scalar * diag1 * mat1;
mat3 += scalar * mat1 * vec1.asDiagonal();
mat3 = vec1.asDiagonal().inverse() * mat1
mat3 = mat1 * diag1.inverse()
\endcode |
\subsection QuickRef_TriangularView Triangular views
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.
Operation | Code |
Reference to a triangular with optional \n
unit or null diagonal (read/write):
| \code
m.triangularView()
\endcode \n
\c Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower
|
Writing to a specific triangular part:\n (only the referenced triangular part is evaluated)
| \code
m1.triangularView() = m2 + m3 \endcode
|
Conversion to a dense matrix setting the opposite triangular part to zero:
| \code
m2 = m1.triangularView()\endcode
|
Products:
| \code
m3 += s1 * m1.adjoint().triangularView() * m2
m3 -= s1 * m2.conjugate() * m1.adjoint().triangularView() \endcode
|
Solving linear equations:\n
\f$ M_2 := L_1^{-1} M_2 \f$ \n
\f$ M_3 := {L_1^*}^{-1} M_3 \f$ \n
\f$ M_4 := M_4 U_1^{-1} \f$
| \n \code
L1.triangularView().solveInPlace(M2)
L1.triangularView().adjoint().solveInPlace(M3)
U1.triangularView().solveInPlace(M4)\endcode
|
\subsection QuickRef_SelfadjointMatrix Symmetric/selfadjoint views
Just as for triangular matrix, you can reference any triangular part of a square matrix to see it as a selfadjoint
matrix and perform special and optimized operations. Again the opposite triangular part is never referenced and can be
used to store other information.
Operation | Code |
Conversion to a dense matrix:
| \code
m2 = m.selfadjointView();\endcode
|
Product with another general matrix or vector:
| \code
m3 = s1 * m1.conjugate().selfadjointView() * m3;
m3 -= s1 * m3.adjoint() * m1.selfadjointView();\endcode
|
Rank 1 and rank K update: \n
\f$ upper(M_1) += s1 M_2^* M_2 \f$ \n
\f$ lower(M_1) -= M_2 M_2^* \f$
| \n \code
M1.selfadjointView().rankUpdate(M2,s1);
m1.selfadjointView().rankUpdate(m2.adjoint(),-1); \endcode
|
Rank 2 update: (\f$ M += s u v^* + s v u^* \f$)
| \code
M.selfadjointView().rankUpdate(u,v,s);
\endcode
|
Solving linear equations:\n(\f$ M_2 := M_1^{-1} M_2 \f$)
| \code
// via a standard Cholesky factorization
m2 = m1.selfadjointView().llt().solve(m2);
// via a Cholesky factorization with pivoting
m2 = m1.selfadjointView().ldlt().solve(m2);
\endcode
|
*/
/*
\link MatrixBase::asDiagonal() make a diagonal matrix \endlink \n from a vector | \code
mat1 = vec1.asDiagonal();\endcode
|
Declare a diagonal matrix | \code
DiagonalMatrix diag1(size);
diag1.diagonal() = vector;\endcode
|
Access \link MatrixBase::diagonal() the diagonal and super/sub diagonals of a matrix \endlink as a vector (read/write) |
\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
\endcode |
View on a triangular part of a matrix (read/write) |
\code
mat2 = mat1.triangularView();
// Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower
mat1.triangularView() = mat2 + mat3; // only the upper part is evaluated and referenced
\endcode |
View a triangular part as a symmetric/self-adjoint matrix (read/write) |
\code
mat2 = mat1.selfadjointView(); // Xxx = Upper or Lower
mat1.selfadjointView() = mat2 + mat2.adjoint(); // evaluated and write to the upper triangular part only
\endcode |
Optimized products:
\code
mat3 += scalar * vec1.asDiagonal() * mat1
mat3 += scalar * mat1 * vec1.asDiagonal()
mat3.noalias() += scalar * mat1.triangularView