template<typename Derived>
Eigen::DenseBase class

Base class for all dense matrices, vectors, and arrays.

Template parameters
Derived is the derived type, e.g., a matrix type or an expression.

This class is the base that is inherited by all dense objects (matrix, vector, arrays, and related expression types). The common Eigen API for dense objects is contained in this class.

This class can be extended with the help of the plugin mechanism described on the page Extending MatrixBase (and other classes) by defining the preprocessor symbol EIGEN_DENSEBASE_PLUGIN.

Base classes

template<typename Derived>
class DenseCoeffsBase<Derived, DirectWriteAccessors>
Base class providing direct read/write coefficient access to matrices and arrays.

Derived classes

template<typename Derived>
class ArrayBase
Base class for all 1D and 2D array, and related expressions.
template<typename Derived>
class ArrayBase
Base class for all 1D and 2D array, and related expressions.
template<typename Derived>
class MatrixBase
Base class for all dense matrices, vectors, and expressions.
template<typename Derived>
class MatrixBase
Base class for all dense matrices, vectors, and expressions.
template<typename Derived>
class MatrixBase
Base class for all dense matrices, vectors, and expressions.
template<typename Derived>
class MatrixBase
Base class for all dense matrices, vectors, and expressions.

Public types

enum (anonymous) { RowsAtCompileTime = internal::traits<Derived>::RowsAtCompileTime, ColsAtCompileTime = internal::traits<Derived>::ColsAtCompileTime, SizeAtCompileTime = (internal::size_at_compile_time<internal::traits<Derived>::RowsAtCompileTime, internal::traits<Derived>::ColsAtCompileTime>::ret), MaxRowsAtCompileTime = internal::traits<Derived>::MaxRowsAtCompileTime, MaxColsAtCompileTime = internal::traits<Derived>::MaxColsAtCompileTime, MaxSizeAtCompileTime = (internal::size_at_compile_time<internal::traits<Derived>::MaxRowsAtCompileTime, internal::traits<Derived>::MaxColsAtCompileTime>::ret), IsVectorAtCompileTime = internal::traits<Derived>::RowsAtCompileTime == 1 || internal::traits<Derived>::ColsAtCompileTime == 1, NumDimensions = int(MaxSizeAtCompileTime) == 1 ? 0 : bool(IsVectorAtCompileTime) ? 1 : 2, Flags = internal::traits<Derived>::Flags, IsRowMajor = int(Flags) & RowMajorBit, InnerSizeAtCompileTime = int(IsVectorAtCompileTime) ? int(SizeAtCompileTime) : int(IsRowMajor) ? int(ColsAtCompileTime) : int(RowsAtCompileTime), InnerStrideAtCompileTime = internal::inner_stride_at_compile_time<Derived>::ret, OuterStrideAtCompileTime = internal::outer_stride_at_compile_time<Derived>::ret }
using const_iterator = random_access_iterator_type
using iterator = random_access_iterator_type
using PlainArray = Array<typename internal::traits<Derived>::Scalar, internal::traits<Derived>::RowsAtCompileTime, internal::traits<Derived>::ColsAtCompileTime, AutoAlign|(internal::traits<Derived>::Flags&RowMajorBit ? RowMajor :ColMajor), internal::traits<Derived>::MaxRowsAtCompileTime, internal::traits<Derived>::MaxColsAtCompileTime>
using PlainMatrix = Matrix<typename internal::traits<Derived>::Scalar, internal::traits<Derived>::RowsAtCompileTime, internal::traits<Derived>::ColsAtCompileTime, AutoAlign|(internal::traits<Derived>::Flags&RowMajorBit ? RowMajor :ColMajor), internal::traits<Derived>::MaxRowsAtCompileTime, internal::traits<Derived>::MaxColsAtCompileTime>
using PlainObject = internal::conditional<internal::is_same<typename internal::traits<Derived>::XprKind, MatrixXpr>::value, PlainMatrix, PlainArray>::type
The plain matrix or array type corresponding to this expression.
using Scalar = internal::traits<Derived>::Scalar
using StorageIndex = internal::traits<Derived>::StorageIndex
The type used to store indices.
using value_type = Scalar

Public static functions

static auto Constant(Index rows, Index cols, const Scalar& value) -> const ConstantReturnType
static auto Constant(Index size, const Scalar& value) -> const ConstantReturnType
static auto Constant(const Scalar& value) -> const ConstantReturnType
static auto LinSpaced(Sequential_t, Index size, const Scalar& low, const Scalar& high) -> const SequentialLinSpacedReturnType deprecated
static auto LinSpaced(Index size, const Scalar& low, const Scalar& high) -> const RandomAccessLinSpacedReturnType
Sets a linearly spaced vector.
static auto LinSpaced(Sequential_t, const Scalar& low, const Scalar& high) -> const SequentialLinSpacedReturnType deprecated
static auto LinSpaced(const Scalar& low, const Scalar& high) -> const RandomAccessLinSpacedReturnType
Sets a linearly spaced vector.
template<typename CustomNullaryOp>
static auto NullaryExpr(Index rows, Index cols, const CustomNullaryOp& func) -> const CwiseNullaryOp<CustomNullaryOp, PlainObject>
template<typename CustomNullaryOp>
static auto NullaryExpr(Index size, const CustomNullaryOp& func) -> const CwiseNullaryOp<CustomNullaryOp, PlainObject>
template<typename CustomNullaryOp>
static auto NullaryExpr(const CustomNullaryOp& func) -> const CwiseNullaryOp<CustomNullaryOp, PlainObject>
static auto Ones(Index rows, Index cols) -> const ConstantReturnType
static auto Ones(Index size) -> const ConstantReturnType
static auto Ones() -> const ConstantReturnType
static auto Random(Index rows, Index cols) -> const RandomReturnType
static auto Random(Index size) -> const RandomReturnType
static auto Random() -> const RandomReturnType
static auto Zero(Index rows, Index cols) -> const ConstantReturnType
static auto Zero(Index size) -> const ConstantReturnType
static auto Zero() -> const ConstantReturnType

Constructors, destructors, conversion operators

DenseBase() protected

Public functions

auto all() const -> bool
auto allFinite() const -> bool
auto any() const -> bool
auto begin() -> iterator
auto begin() const -> const_iterator
template<typename NRowsType, typename NColsType>
auto block(Index startRow, Index startCol, NRowsType blockRows, NColsType blockCols) -> FixedBlockXpr<...,...>::Type
template<typename NRowsType, typename NColsType>
auto block(Index startRow, Index startCol, NRowsType blockRows, NColsType blockCols) const -> const ConstFixedBlockXpr<...,...>::Type
This is the const version of block(Index,Index,NRowsType,NColsType)
template<int NRows, int NCols>
auto block(Index startRow, Index startCol) -> FixedBlockXpr<NRows, NCols>::Type
template<int NRows, int NCols>
auto block(Index startRow, Index startCol) const -> const ConstFixedBlockXpr<NRows, NCols>::Type
This is the const version of block<>(Index, Index). */.
template<int NRows, int NCols>
auto block(Index startRow, Index startCol, Index blockRows, Index blockCols) -> FixedBlockXpr<NRows, NCols>::Type
template<int NRows, int NCols>
auto block(Index startRow, Index startCol, Index blockRows, Index blockCols) const -> const ConstFixedBlockXpr<NRows, NCols>::Type
This is the const version of block<>(Index, Index, Index, Index).
template<typename NRowsType, typename NColsType>
auto bottomLeftCorner(NRowsType cRows, NColsType cCols) -> FixedBlockXpr<...,...>::Type
template<typename NRowsType, typename NColsType>
auto bottomLeftCorner(NRowsType cRows, NColsType cCols) const -> ConstFixedBlockXpr<...,...>::Type
This is the const version of bottomLeftCorner(NRowsType, NColsType).
template<int CRows, int CCols>
auto bottomLeftCorner() -> FixedBlockXpr<CRows, CCols>::Type
template<int CRows, int CCols>
auto bottomLeftCorner() const -> const ConstFixedBlockXpr<CRows, CCols>::Type
This is the const version of bottomLeftCorner<int, int>().
template<int CRows, int CCols>
auto bottomLeftCorner(Index cRows, Index cCols) -> FixedBlockXpr<CRows, CCols>::Type
template<int CRows, int CCols>
auto bottomLeftCorner(Index cRows, Index cCols) const -> const ConstFixedBlockXpr<CRows, CCols>::Type
This is the const version of bottomLeftCorner<int, int>(Index, Index).
template<typename NRowsType, typename NColsType>
auto bottomRightCorner(NRowsType cRows, NColsType cCols) -> FixedBlockXpr<...,...>::Type
template<typename NRowsType, typename NColsType>
auto bottomRightCorner(NRowsType cRows, NColsType cCols) const -> const ConstFixedBlockXpr<...,...>::Type
This is the const version of bottomRightCorner(NRowsType, NColsType).
template<int CRows, int CCols>
auto bottomRightCorner() -> FixedBlockXpr<CRows, CCols>::Type
template<int CRows, int CCols>
auto bottomRightCorner() const -> const ConstFixedBlockXpr<CRows, CCols>::Type
This is the const version of bottomRightCorner<int, int>().
template<int CRows, int CCols>
auto bottomRightCorner(Index cRows, Index cCols) -> FixedBlockXpr<CRows, CCols>::Type
template<int CRows, int CCols>
auto bottomRightCorner(Index cRows, Index cCols) const -> const ConstFixedBlockXpr<CRows, CCols>::Type
This is the const version of bottomRightCorner<int, int>(Index, Index).
template<typename NRowsType>
auto bottomRows(NRowsType n) -> NRowsBlockXpr<...>::Type
template<typename NRowsType>
auto bottomRows(NRowsType n) const -> const ConstNRowsBlockXpr<...>::Type
This is the const version of bottomRows(NRowsType).
template<int N>
auto bottomRows(Index n = N) -> NRowsBlockXpr<N>::Type
template<int N>
auto bottomRows(Index n = N) const -> ConstNRowsBlockXpr<N>::Type
This is the const version of bottomRows<int>().
template<typename NewType>
auto cast() const -> CastXpr<NewType>::Type
auto cbegin() const -> const_iterator
auto cend() const -> const_iterator
auto col(Index i) -> ColXpr
auto col(Index i) const -> ConstColXpr
This is the const version of col().
auto colwise() const -> ConstColwiseReturnType
auto colwise() -> ColwiseReturnType
auto conjugate() const -> ConjugateReturnType
template<bool Cond>
auto conjugateIf() const -> internal::conditional<Cond, ConjugateReturnType, const Derived&>::type
auto count() const -> Index
auto end() -> iterator
auto end() const -> const_iterator
auto eval() const -> EvalReturnType
void fill(const Scalar& value)
template<unsigned int Added, unsigned int Removed>
auto flagged() const -> EIGEN_DEPRECATED const Derived& deprecated
auto format(const IOFormat& fmt) const -> const WithFormat<Derived>
auto hasNaN() const -> bool
template<typename NType>
auto head(NType n) -> FixedSegmentReturnType<...>::Type
template<typename NType>
auto head(NType n) const -> const ConstFixedSegmentReturnType<...>::Type
This is the const version of head(NType).
template<int N>
auto head(Index n = N) -> FixedSegmentReturnType<N>::Type
template<int N>
auto head(Index n = N) const -> ConstFixedSegmentReturnType<N>::Type
This is the const version of head<int>().
auto imag() const -> const ImagReturnType
auto imag() -> NonConstImagReturnType
auto innerSize() const -> Index
auto innerVector(Index outer) -> InnerVectorReturnType
auto innerVector(Index outer) const -> const ConstInnerVectorReturnType
auto innerVectors(Index outerStart, Index outerSize) -> InnerVectorsReturnType
auto innerVectors(Index outerStart, Index outerSize) const -> const ConstInnerVectorsReturnType
template<typename OtherDerived>
auto isApprox(const DenseBase<OtherDerived>& other, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
auto isApproxToConstant(const Scalar& value, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
auto isConstant(const Scalar& value, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
template<typename Derived>
auto isMuchSmallerThan(const typename NumTraits<Scalar>::Real& other, const RealScalar& prec) const -> bool
template<typename OtherDerived>
auto isMuchSmallerThan(const DenseBase<OtherDerived>& other, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
auto isOnes(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
auto isZero(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
template<typename NColsType>
auto leftCols(NColsType n) -> NColsBlockXpr<...>::Type
template<typename NColsType>
auto leftCols(NColsType n) const -> const ConstNColsBlockXpr<...>::Type
This is the const version of leftCols(NColsType).
template<int N>
auto leftCols(Index n = N) -> NColsBlockXpr<N>::Type
template<int N>
auto leftCols(Index n = N) const -> ConstNColsBlockXpr<N>::Type
This is the const version of leftCols<int>().
auto maxCoeff() const -> internal::traits<Derived>::Scalar
template<typename IndexType>
auto maxCoeff(IndexType* row, IndexType* col) const -> internal::traits<Derived>::Scalar
template<typename IndexType>
auto maxCoeff(IndexType* index) const -> internal::traits<Derived>::Scalar
auto mean() const -> Scalar
template<typename NColsType>
auto middleCols(Index startCol, NColsType numCols) -> NColsBlockXpr<...>::Type
template<typename NColsType>
auto middleCols(Index startCol, NColsType numCols) const -> const ConstNColsBlockXpr<...>::Type
This is the const version of middleCols(Index,NColsType).
template<int N>
auto middleCols(Index startCol, Index n = N) -> NColsBlockXpr<N>::Type
template<int N>
auto middleCols(Index startCol, Index n = N) const -> ConstNColsBlockXpr<N>::Type
This is the const version of middleCols<int>().
template<typename NRowsType>
auto middleRows(Index startRow, NRowsType n) -> NRowsBlockXpr<...>::Type
template<typename NRowsType>
auto middleRows(Index startRow, NRowsType n) const -> const ConstNRowsBlockXpr<...>::Type
This is the const version of middleRows(Index,NRowsType).
template<int N>
auto middleRows(Index startRow, Index n = N) -> NRowsBlockXpr<N>::Type
template<int N>
auto middleRows(Index startRow, Index n = N) const -> ConstNRowsBlockXpr<N>::Type
This is the const version of middleRows<int>().
auto minCoeff() const -> internal::traits<Derived>::Scalar
template<typename IndexType>
auto minCoeff(IndexType* row, IndexType* col) const -> internal::traits<Derived>::Scalar
template<typename IndexType>
auto minCoeff(IndexType* index) const -> internal::traits<Derived>::Scalar
auto nestByValue() const -> const NestByValue<Derived>
auto nonZeros() const -> Index
template<typename RowIndices, typename ColIndices>
auto operator()(const RowIndices& rowIndices, const ColIndices& colIndices) -> IndexedView_or_Block
template<typename Indices>
auto operator()(const Indices& indices) -> IndexedView_or_VectorBlock
auto operator-() const -> const NegativeReturnType
auto operator<<(const Scalar& s) -> CommaInitializer<Derived>
template<typename OtherDerived>
auto operator<<(const DenseBase<OtherDerived>& other) -> CommaInitializer<Derived>
template<typename OtherDerived>
auto operator=(const DenseBase<OtherDerived>& other) -> Derived&
auto operator=(const DenseBase& other) -> Derived&
template<typename OtherDerived>
auto operator=(const EigenBase<OtherDerived>& other) -> Derived&
Copies the generic expression other into *this.
auto outerSize() const -> Index
auto prod() const -> Scalar
auto real() const -> RealReturnType
auto real() -> NonConstRealReturnType
template<typename Func>
auto redux(const Func& func) const -> internal::traits<Derived>::Scalar
template<int RowFactor, int ColFactor>
auto replicate() const -> const Replicate<Derived, RowFactor, ColFactor>
auto replicate(Index rowFactor, Index colFactor) const -> const Replicate<Derived, Dynamic, Dynamic>
template<int Order = ColMajor, typename NRowsType, typename NColsType>
auto reshaped(NRowsType nRows, NColsType nCols) -> Reshaped<Derived,...>
template<int Order = ColMajor, typename NRowsType, typename NColsType>
auto reshaped(NRowsType nRows, NColsType nCols) const -> const Reshaped<const Derived,...>
This is the const version of reshaped(NRowsType,NColsType).
template<int Order = ColMajor>
auto reshaped() -> Reshaped<Derived,...>
template<int Order = ColMajor>
auto reshaped() const -> const Reshaped<const Derived,...>
This is the const version of reshaped().
void resize(Index newSize)
void resize(Index rows, Index cols)
auto reverse() -> ReverseReturnType
auto reverse() const -> ConstReverseReturnType
void reverseInPlace()
template<typename NColsType>
auto rightCols(NColsType n) -> NColsBlockXpr<...>::Type
template<typename NColsType>
auto rightCols(NColsType n) const -> const ConstNColsBlockXpr<...>::Type
This is the const version of rightCols(NColsType).
template<int N>
auto rightCols(Index n = N) -> NColsBlockXpr<N>::Type
template<int N>
auto rightCols(Index n = N) const -> ConstNColsBlockXpr<N>::Type
This is the const version of rightCols<int>().
auto row(Index i) -> RowXpr
auto row(Index i) const -> ConstRowXpr
This is the const version of row(). */.
auto rowwise() const -> ConstRowwiseReturnType
auto rowwise() -> RowwiseReturnType
template<typename NType>
auto segment(Index start, NType n) -> FixedSegmentReturnType<...>::Type
template<typename NType>
auto segment(Index start, NType n) const -> const ConstFixedSegmentReturnType<...>::Type
This is the const version of segment(Index,NType).
template<int N>
auto segment(Index start, Index n = N) -> FixedSegmentReturnType<N>::Type
template<int N>
auto segment(Index start, Index n = N) const -> ConstFixedSegmentReturnType<N>::Type
This is the const version of segment<int>(Index).
template<typename ThenDerived, typename ElseDerived>
auto select(const DenseBase<ThenDerived>& thenMatrix, const DenseBase<ElseDerived>& elseMatrix) const -> const Select<Derived, ThenDerived, ElseDerived>
template<typename ThenDerived>
auto select(const DenseBase<ThenDerived>& thenMatrix, const typename ThenDerived::Scalar& elseScalar) const -> const Select<Derived, ThenDerived, typename ThenDerived::ConstantReturnType>
template<typename ElseDerived>
auto select(const typename ElseDerived::Scalar& thenScalar, const DenseBase<ElseDerived>& elseMatrix) const -> const Select<Derived, typename ElseDerived::ConstantReturnType, ElseDerived>
auto setConstant(const Scalar& value) -> Derived&
auto setLinSpaced(Index size, const Scalar& low, const Scalar& high) -> Derived&
Sets a linearly spaced vector.
auto setLinSpaced(const Scalar& low, const Scalar& high) -> Derived&
Sets a linearly spaced vector.
auto setOnes() -> Derived&
auto setRandom() -> Derived&
auto setZero() -> Derived&
template<DirectionType Direction>
auto subVector(Index i) -> internal::conditional<Direction==Vertical, ColXpr, RowXpr>::type
template<DirectionType Direction>
auto subVector(Index i) const -> internal::conditional<Direction==Vertical, ConstColXpr, ConstRowXpr>::type
template<DirectionType Direction>
auto subVectors() const -> Index
auto sum() const -> Scalar
template<typename OtherDerived>
void swap(const DenseBase<OtherDerived>& other)
template<typename OtherDerived>
void swap(PlainObjectBase<OtherDerived>& other)
template<typename NType>
auto tail(NType n) -> FixedSegmentReturnType<...>::Type
template<typename NType>
auto tail(NType n) const -> const ConstFixedSegmentReturnType<...>::Type
This is the const version of tail(Index).
template<int N>
auto tail(Index n = N) -> FixedSegmentReturnType<N>::Type
template<int N>
auto tail(Index n = N) const -> ConstFixedSegmentReturnType<N>::Type
This is the const version of tail<int>.
template<typename NRowsType, typename NColsType>
auto topLeftCorner(NRowsType cRows, NColsType cCols) -> FixedBlockXpr<...,...>::Type
template<typename NRowsType, typename NColsType>
auto topLeftCorner(NRowsType cRows, NColsType cCols) const -> const ConstFixedBlockXpr<...,...>::Type
This is the const version of topLeftCorner(Index, Index).
template<int CRows, int CCols>
auto topLeftCorner() -> FixedBlockXpr<CRows, CCols>::Type
template<int CRows, int CCols>
auto topLeftCorner() const -> const ConstFixedBlockXpr<CRows, CCols>::Type
This is the const version of topLeftCorner<int, int>().
template<int CRows, int CCols>
auto topLeftCorner(Index cRows, Index cCols) -> FixedBlockXpr<CRows, CCols>::Type
template<int CRows, int CCols>
auto topLeftCorner(Index cRows, Index cCols) const -> const ConstFixedBlockXpr<CRows, CCols>::Type
This is the const version of topLeftCorner<int, int>(Index, Index).
template<typename NRowsType, typename NColsType>
auto topRightCorner(NRowsType cRows, NColsType cCols) -> FixedBlockXpr<...,...>::Type
template<typename NRowsType, typename NColsType>
auto topRightCorner(NRowsType cRows, NColsType cCols) const -> const ConstFixedBlockXpr<...,...>::Type
This is the const version of topRightCorner(NRowsType, NColsType).
template<int CRows, int CCols>
auto topRightCorner() -> FixedBlockXpr<CRows, CCols>::Type
template<int CRows, int CCols>
auto topRightCorner() const -> const ConstFixedBlockXpr<CRows, CCols>::Type
This is the const version of topRightCorner<int, int>().
template<int CRows, int CCols>
auto topRightCorner(Index cRows, Index cCols) -> FixedBlockXpr<CRows, CCols>::Type
template<int CRows, int CCols>
auto topRightCorner(Index cRows, Index cCols) const -> const ConstFixedBlockXpr<CRows, CCols>::Type
This is the const version of topRightCorner<int, int>(Index, Index).
template<typename NRowsType>
auto topRows(NRowsType n) -> NRowsBlockXpr<...>::Type
template<typename NRowsType>
auto topRows(NRowsType n) const -> const ConstNRowsBlockXpr<...>::Type
This is the const version of topRows(NRowsType).
template<int N>
auto topRows(Index n = N) -> NRowsBlockXpr<N>::Type
template<int N>
auto topRows(Index n = N) const -> ConstNRowsBlockXpr<N>::Type
This is the const version of topRows<int>().
auto transpose() -> TransposeReturnType
auto transpose() const -> ConstTransposeReturnType
void transposeInPlace()
template<typename CustomUnaryOp>
auto unaryExpr(const CustomUnaryOp& func = CustomUnaryOp()) const -> const CwiseUnaryOp<CustomUnaryOp, const Derived>
Apply a unary operator coefficient-wise.
template<typename CustomViewOp>
auto unaryViewExpr(const CustomViewOp& func = CustomViewOp()) const -> const CwiseUnaryView<CustomViewOp, const Derived>
auto value() const -> CoeffReturnType
template<typename Visitor>
void visit(Visitor& func) const

Enum documentation

template<typename Derived>
enum Eigen::DenseBase<Derived>::(anonymous)

Enumerators
RowsAtCompileTime

The number of rows at compile-time. This is just a copy of the value provided by the Derived type. If a value is not known at compile-time, it is set to the Dynamic constant.

ColsAtCompileTime

The number of columns at compile-time. This is just a copy of the value provided by the Derived type. If a value is not known at compile-time, it is set to the Dynamic constant.

SizeAtCompileTime

This is equal to the number of coefficients, i.e. the number of rows times the number of columns, or to Dynamic if this is not known at compile-time.

MaxRowsAtCompileTime

This value is equal to the maximum possible number of rows that this expression might have. If this expression might have an arbitrarily high number of rows, this value is set to Dynamic.

This value is useful to know when evaluating an expression, in order to determine whether it is possible to avoid doing a dynamic memory allocation.

MaxColsAtCompileTime

This value is equal to the maximum possible number of columns that this expression might have. If this expression might have an arbitrarily high number of columns, this value is set to Dynamic.

This value is useful to know when evaluating an expression, in order to determine whether it is possible to avoid doing a dynamic memory allocation.

MaxSizeAtCompileTime

This value is equal to the maximum possible number of coefficients that this expression might have. If this expression might have an arbitrarily high number of coefficients, this value is set to Dynamic.

This value is useful to know when evaluating an expression, in order to determine whether it is possible to avoid doing a dynamic memory allocation.

IsVectorAtCompileTime

This is set to true if either the number of rows or the number of columns is known at compile-time to be equal to 1. Indeed, in that case, we are dealing with a column-vector (if there is only one column) or with a row-vector (if there is only one row).

NumDimensions

This value is equal to Tensor::NumDimensions, i.e. 0 for scalars, 1 for vectors, and 2 for matrices.

Flags

This stores expression Flags flags which may or may not be inherited by new expressions constructed from this one. See the list of flags.

IsRowMajor

True if this expression has row-major storage order.

InnerSizeAtCompileTime
InnerStrideAtCompileTime
OuterStrideAtCompileTime

Typedef documentation

template<typename Derived>
typedef random_access_iterator_type Eigen::DenseBase<Derived>::const_iterator

This is the const version of iterator (aka read-only)

template<typename Derived>
typedef random_access_iterator_type Eigen::DenseBase<Derived>::iterator

STL-like RandomAccessIterator iterator type as returned by the begin() and end() methods.

template<typename Derived>
typedef Array<typename internal::traits<Derived>::Scalar, internal::traits<Derived>::RowsAtCompileTime, internal::traits<Derived>::ColsAtCompileTime, AutoAlign|(internal::traits<Derived>::Flags&RowMajorBit ? RowMajor :ColMajor), internal::traits<Derived>::MaxRowsAtCompileTime, internal::traits<Derived>::MaxColsAtCompileTime> Eigen::DenseBase<Derived>::PlainArray

The plain array type corresponding to this expression.

template<typename Derived>
typedef Matrix<typename internal::traits<Derived>::Scalar, internal::traits<Derived>::RowsAtCompileTime, internal::traits<Derived>::ColsAtCompileTime, AutoAlign|(internal::traits<Derived>::Flags&RowMajorBit ? RowMajor :ColMajor), internal::traits<Derived>::MaxRowsAtCompileTime, internal::traits<Derived>::MaxColsAtCompileTime> Eigen::DenseBase<Derived>::PlainMatrix

The plain matrix type corresponding to this expression.

template<typename Derived>
typedef internal::conditional<internal::is_same<typename internal::traits<Derived>::XprKind, MatrixXpr>::value, PlainMatrix, PlainArray>::type Eigen::DenseBase<Derived>::PlainObject

The plain matrix or array type corresponding to this expression.

This is not necessarily exactly the return type of eval(). In the case of plain matrices, the return type of eval() is a const reference to a matrix, not a matrix! It is however guaranteed that the return type of eval() is either PlainObject or const PlainObject&.

template<typename Derived>
typedef internal::traits<Derived>::Scalar Eigen::DenseBase<Derived>::Scalar

The numeric type of the expression' coefficients, e.g. float, double, int or std::complex<float>, etc.

template<typename Derived>
typedef internal::traits<Derived>::StorageIndex Eigen::DenseBase<Derived>::StorageIndex

The type used to store indices.

This typedef is relevant for types that store multiple indices such as PermutationMatrix or Transpositions, otherwise it defaults to Eigen::Index

template<typename Derived>
typedef Scalar Eigen::DenseBase<Derived>::value_type

The numeric type of the expression' coefficients, e.g. float, double, int or std::complex<float>, etc.

It is an alias for the Scalar type

Function documentation

template<typename Derived>
static const ConstantReturnType Eigen::DenseBase<Derived>::Constant(Index rows, Index cols, const Scalar& value)

Returns an expression of a constant matrix of value value

The parameters rows and cols are the number of rows and of columns of the returned matrix. Must be compatible with this DenseBase type.

This variant is meant to be used for dynamic-size matrix types. For fixed-size types, it is redundant to pass rows and cols as arguments, so Zero() should be used instead.

The template parameter CustomNullaryOp is the type of the functor.

template<typename Derived>
static const ConstantReturnType Eigen::DenseBase<Derived>::Constant(Index size, const Scalar& value)

Returns an expression of a constant matrix of value value

The parameter size is the size of the returned vector. Must be compatible with this DenseBase type.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

This variant is meant to be used for dynamic-size vector types. For fixed-size types, it is redundant to pass size as argument, so Zero() should be used instead.

The template parameter CustomNullaryOp is the type of the functor.

template<typename Derived>
static const ConstantReturnType Eigen::DenseBase<Derived>::Constant(const Scalar& value)

Returns an expression of a constant matrix of value value

This variant is only for fixed-size DenseBase types. For dynamic-size types, you need to use the variants taking size arguments.

The template parameter CustomNullaryOp is the type of the functor.

template<typename Derived>
static const SequentialLinSpacedReturnType Eigen::DenseBase<Derived>::LinSpaced(Sequential_t, Index size, const Scalar& low, const Scalar& high)

template<typename Derived>
static const RandomAccessLinSpacedReturnType Eigen::DenseBase<Derived>::LinSpaced(Index size, const Scalar& low, const Scalar& high)

Sets a linearly spaced vector.

The function generates 'size' equally spaced values in the closed interval [low,high]. When size is set to 1, a vector of length 1 containing 'high' is returned.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

Example:

cout << VectorXi::LinSpaced(4,7,10).transpose() << endl;
cout << VectorXd::LinSpaced(5,0.0,1.0).transpose() << endl;

Output:

 7  8  9 10
   0 0.25  0.5 0.75    1

For integer scalar types, an even spacing is possible if and only if the length of the range, i.e., high-low is a scalar multiple of size-1, or if size is a scalar multiple of the number of values high-low+1 (meaning each value can be repeated the same number of time). If one of these two considions is not satisfied, then high is lowered to the largest value satisfying one of this constraint. Here are some examples:

Example:

cout << "Even spacing inputs:" << endl;
cout << VectorXi::LinSpaced(8,1,4).transpose() << endl;
cout << VectorXi::LinSpaced(8,1,8).transpose() << endl;
cout << VectorXi::LinSpaced(8,1,15).transpose() << endl;
cout << "Uneven spacing inputs:" << endl;
cout << VectorXi::LinSpaced(8,1,7).transpose() << endl;
cout << VectorXi::LinSpaced(8,1,9).transpose() << endl;
cout << VectorXi::LinSpaced(8,1,16).transpose() << endl;

Output:

Even spacing inputs:
1 1 2 2 3 3 4 4
1 2 3 4 5 6 7 8
 1  3  5  7  9 11 13 15
Uneven spacing inputs:
1 1 2 2 3 3 4 4
1 2 3 4 5 6 7 8
 1  3  5  7  9 11 13 15

template<typename Derived>
static const SequentialLinSpacedReturnType Eigen::DenseBase<Derived>::LinSpaced(Sequential_t, const Scalar& low, const Scalar& high)

template<typename Derived>
static const RandomAccessLinSpacedReturnType Eigen::DenseBase<Derived>::LinSpaced(const Scalar& low, const Scalar& high)

Sets a linearly spaced vector.

The function generates 'size' equally spaced values in the closed interval [low,high]. When size is set to 1, a vector of length 1 containing 'high' is returned.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

Example:

cout << VectorXi::LinSpaced(4,7,10).transpose() << endl;
cout << VectorXd::LinSpaced(5,0.0,1.0).transpose() << endl;

Output:

 7  8  9 10
   0 0.25  0.5 0.75    1

For integer scalar types, an even spacing is possible if and only if the length of the range, i.e., high-low is a scalar multiple of size-1, or if size is a scalar multiple of the number of values high-low+1 (meaning each value can be repeated the same number of time). If one of these two considions is not satisfied, then high is lowered to the largest value satisfying one of this constraint. Here are some examples:

Example:

cout << "Even spacing inputs:" << endl;
cout << VectorXi::LinSpaced(8,1,4).transpose() << endl;
cout << VectorXi::LinSpaced(8,1,8).transpose() << endl;
cout << VectorXi::LinSpaced(8,1,15).transpose() << endl;
cout << "Uneven spacing inputs:" << endl;
cout << VectorXi::LinSpaced(8,1,7).transpose() << endl;
cout << VectorXi::LinSpaced(8,1,9).transpose() << endl;
cout << VectorXi::LinSpaced(8,1,16).transpose() << endl;

Output:

Even spacing inputs:
1 1 2 2 3 3 4 4
1 2 3 4 5 6 7 8
 1  3  5  7  9 11 13 15
Uneven spacing inputs:
1 1 2 2 3 3 4 4
1 2 3 4 5 6 7 8
 1  3  5  7  9 11 13 15

template<typename Derived> template<typename CustomNullaryOp>
static const CwiseNullaryOp<CustomNullaryOp, PlainObject> Eigen::DenseBase<Derived>::NullaryExpr(Index rows, Index cols, const CustomNullaryOp& func)

Returns an expression of a matrix defined by a custom functor func

The parameters rows and cols are the number of rows and of columns of the returned matrix. Must be compatible with this MatrixBase type.

This variant is meant to be used for dynamic-size matrix types. For fixed-size types, it is redundant to pass rows and cols as arguments, so Zero() should be used instead.

The template parameter CustomNullaryOp is the type of the functor.

template<typename Derived> template<typename CustomNullaryOp>
static const CwiseNullaryOp<CustomNullaryOp, PlainObject> Eigen::DenseBase<Derived>::NullaryExpr(Index size, const CustomNullaryOp& func)

Returns an expression of a matrix defined by a custom functor func

The parameter size is the size of the returned vector. Must be compatible with this MatrixBase type.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

This variant is meant to be used for dynamic-size vector types. For fixed-size types, it is redundant to pass size as argument, so Zero() should be used instead.

The template parameter CustomNullaryOp is the type of the functor.

Here is an example with C++11 random generators:

#include <Eigen/Core>
#include <iostream>
#include <random>

using namespace Eigen;

int main() {
  std::default_random_engine generator;
  std::poisson_distribution<int> distribution(4.1);
  auto poisson = [&] () {return distribution(generator);};

  RowVectorXi v = RowVectorXi::NullaryExpr(10, poisson );
  std::cout << v << "\n";
}

Output:

2 3 1 4 3 4 4 3 2 3

template<typename Derived> template<typename CustomNullaryOp>
static const CwiseNullaryOp<CustomNullaryOp, PlainObject> Eigen::DenseBase<Derived>::NullaryExpr(const CustomNullaryOp& func)

Returns an expression of a matrix defined by a custom functor func

This variant is only for fixed-size DenseBase types. For dynamic-size types, you need to use the variants taking size arguments.

The template parameter CustomNullaryOp is the type of the functor.

template<typename Derived>
static const ConstantReturnType Eigen::DenseBase<Derived>::Ones(Index rows, Index cols)

Returns an expression of a matrix where all coefficients equal one.

The parameters rows and cols are the number of rows and of columns of the returned matrix. Must be compatible with this MatrixBase type.

This variant is meant to be used for dynamic-size matrix types. For fixed-size types, it is redundant to pass rows and cols as arguments, so Ones() should be used instead.

Example:

cout << MatrixXi::Ones(2,3) << endl;

Output:

1 1 1
1 1 1

template<typename Derived>
static const ConstantReturnType Eigen::DenseBase<Derived>::Ones(Index size)

Returns an expression of a vector where all coefficients equal one.

The parameter newSize is the size of the returned vector. Must be compatible with this MatrixBase type.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

This variant is meant to be used for dynamic-size vector types. For fixed-size types, it is redundant to pass size as argument, so Ones() should be used instead.

Example:

cout << 6 * RowVectorXi::Ones(4) << endl;
cout << VectorXf::Ones(2) << endl;

Output:

6 6 6 6
1
1

template<typename Derived>
static const ConstantReturnType Eigen::DenseBase<Derived>::Ones()

Returns an expression of a fixed-size matrix or vector where all coefficients equal one.

This variant is only for fixed-size MatrixBase types. For dynamic-size types, you need to use the variants taking size arguments.

Example:

cout << Matrix2d::Ones() << endl;
cout << 6 * RowVector4i::Ones() << endl;

Output:

1 1
1 1
6 6 6 6

template<typename Derived>
static const RandomReturnType Eigen::DenseBase<Derived>::Random(Index rows, Index cols)

Returns a random matrix expression

Numbers are uniformly spread through their whole definition range for integer types, and in the [-1:1] range for floating point scalar types.

The parameters rows and cols are the number of rows and of columns of the returned matrix. Must be compatible with this MatrixBase type.

This variant is meant to be used for dynamic-size matrix types. For fixed-size types, it is redundant to pass rows and cols as arguments, so Random() should be used instead.

Example:

cout << MatrixXi::Random(2,3) << endl;

Output:

 7  6  9
-2  6 -6

This expression has the "evaluate before nesting" flag so that it will be evaluated into a temporary matrix whenever it is nested in a larger expression. This prevents unexpected behavior with expressions involving random matrices.

See DenseBase::NullaryExpr(Index, const CustomNullaryOp&) for an example using C++11 random generators.

template<typename Derived>
static const RandomReturnType Eigen::DenseBase<Derived>::Random(Index size)

Returns a random vector expression

Numbers are uniformly spread through their whole definition range for integer types, and in the [-1:1] range for floating point scalar types.

The parameter size is the size of the returned vector. Must be compatible with this MatrixBase type.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

This variant is meant to be used for dynamic-size vector types. For fixed-size types, it is redundant to pass size as argument, so Random() should be used instead.

Example:

cout << VectorXi::Random(2) << endl;

Output:

 7
-2

This expression has the "evaluate before nesting" flag so that it will be evaluated into a temporary vector whenever it is nested in a larger expression. This prevents unexpected behavior with expressions involving random matrices.

template<typename Derived>
static const RandomReturnType Eigen::DenseBase<Derived>::Random()

Returns a fixed-size random matrix or vector expression

Numbers are uniformly spread through their whole definition range for integer types, and in the [-1:1] range for floating point scalar types.

This variant is only for fixed-size MatrixBase types. For dynamic-size types, you need to use the variants taking size arguments.

Example:

cout << 100 * Matrix2i::Random() << endl;

Output:

 700  600
-200  600

This expression has the "evaluate before nesting" flag so that it will be evaluated into a temporary matrix whenever it is nested in a larger expression. This prevents unexpected behavior with expressions involving random matrices.

template<typename Derived>
static const ConstantReturnType Eigen::DenseBase<Derived>::Zero(Index rows, Index cols)

Returns an expression of a zero matrix.

The parameters rows and cols are the number of rows and of columns of the returned matrix. Must be compatible with this MatrixBase type.

This variant is meant to be used for dynamic-size matrix types. For fixed-size types, it is redundant to pass rows and cols as arguments, so Zero() should be used instead.

Example:

cout << MatrixXi::Zero(2,3) << endl;

Output:

0 0 0
0 0 0

template<typename Derived>
static const ConstantReturnType Eigen::DenseBase<Derived>::Zero(Index size)

Returns an expression of a zero vector.

The parameter size is the size of the returned vector. Must be compatible with this MatrixBase type.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

This variant is meant to be used for dynamic-size vector types. For fixed-size types, it is redundant to pass size as argument, so Zero() should be used instead.

Example:

cout << RowVectorXi::Zero(4) << endl;
cout << VectorXf::Zero(2) << endl;

Output:

0 0 0 0
0
0

template<typename Derived>
static const ConstantReturnType Eigen::DenseBase<Derived>::Zero()

Returns an expression of a fixed-size zero matrix or vector.

This variant is only for fixed-size MatrixBase types. For dynamic-size types, you need to use the variants taking size arguments.

Example:

cout << Matrix2d::Zero() << endl;
cout << RowVector4i::Zero() << endl;

Output:

0 0
0 0
0 0 0 0

template<typename Derived>
Eigen::DenseBase<Derived>::DenseBase() protected

Default constructor. Do nothing.

template<typename Derived>
bool Eigen::DenseBase<Derived>::all() const

Returns true if all coefficients are true

Example:

Vector3f boxMin(Vector3f::Zero()), boxMax(Vector3f::Ones());
Vector3f p0 = Vector3f::Random(), p1 = Vector3f::Random().cwiseAbs();
// let's check if p0 and p1 are inside the axis aligned box defined by the corners boxMin,boxMax:
cout << "Is (" << p0.transpose() << ") inside the box: "
     << ((boxMin.array()<p0.array()).all() && (boxMax.array()>p0.array()).all()) << endl;
cout << "Is (" << p1.transpose() << ") inside the box: "
     << ((boxMin.array()<p1.array()).all() && (boxMax.array()>p1.array()).all()) << endl;

Output:

Is (  0.68 -0.211  0.566) inside the box: 0
Is (0.597 0.823 0.605) inside the box: 1

template<typename Derived>
bool Eigen::DenseBase<Derived>::allFinite() const

Returns true if *this contains only finite numbers, i.e., no NaN and no +/-INF values.

template<typename Derived>
bool Eigen::DenseBase<Derived>::any() const

Returns true if at least one coefficient is true

template<typename Derived>
iterator Eigen::DenseBase<Derived>::begin()

returns an iterator to the first element of the 1D vector or array This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

template<typename Derived>
const_iterator Eigen::DenseBase<Derived>::begin() const

const version of begin()

template<typename Derived> template<typename NRowsType, typename NColsType>
FixedBlockXpr<...,...>::Type Eigen::DenseBase<Derived>::block(Index startRow, Index startCol, NRowsType blockRows, NColsType blockCols)

Template parameters
NRowsType the type of the value handling the number of rows in the block, typically Index.
NColsType the type of the value handling the number of columns in the block, typically Index.
Parameters
startRow the first row in the block
startCol the first column in the block
blockRows number of rows in the block, specified at either run-time or compile-time
blockCols number of columns in the block, specified at either run-time or compile-time
Returns an expression of a block in *this with either dynamic or fixed sizes.

Example using runtime (aka dynamic) sizes:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.block(1, 1, 2, 2):" << endl << m.block(1, 1, 2, 2) << endl;
m.block(1, 1, 2, 2).setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is m.block(1, 1, 2, 2):
-6  1
-3  0
Now the matrix m is:
 7  9 -5 -3
-2  0  0  0
 6  0  0  9
 6  6  3  9

New in Eigen 3.4.:

The number of rows blockRows and columns blockCols can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. In the later case, n plays the role of a runtime fallback value in case N equals Eigen::Dynamic. Here is an example with a fixed number of rows NRows and dynamic number of columns cols: mat.block(i,j,fix<NRows>,cols)

This function thus fully covers the features offered by the following overloads block<NRows,NCols>(Index, Index), and block<NRows,NCols>(Index, Index, Index, Index) that are thus obsolete. Indeed, this generic version avoids redundancy, it preserves the argument order, and prevents the need to rely on the template keyword in templated code.

but with less redundancy and more consistency as it does not modify the argument order and seamlessly enable hybrid fixed/dynamic sizes.

template<typename Derived> template<typename NRowsType, typename NColsType>
const ConstFixedBlockXpr<...,...>::Type Eigen::DenseBase<Derived>::block(Index startRow, Index startCol, NRowsType blockRows, NColsType blockCols) const

This is the const version of block(Index,Index,NRowsType,NColsType)

template<typename Derived> template<int NRows, int NCols>
FixedBlockXpr<NRows, NCols>::Type Eigen::DenseBase<Derived>::block(Index startRow, Index startCol)

Parameters
startRow the first row in the block
startCol the first column in the block
Returns a fixed-size expression of a block of *this.

The template parameters NRows and NCols are the number of rows and columns in the block.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.block<2,2>(1,1):" << endl << m.block<2,2>(1,1) << endl;
m.block<2,2>(1,1).setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is m.block<2,2>(1,1):
-6  1
-3  0
Now the matrix m is:
 7  9 -5 -3
-2  0  0  0
 6  0  0  9
 6  6  3  9

template<typename Derived> template<int NRows, int NCols>
const ConstFixedBlockXpr<NRows, NCols>::Type Eigen::DenseBase<Derived>::block(Index startRow, Index startCol) const

This is the const version of block<>(Index, Index). */.

template<typename Derived> template<int NRows, int NCols>
FixedBlockXpr<NRows, NCols>::Type Eigen::DenseBase<Derived>::block(Index startRow, Index startCol, Index blockRows, Index blockCols)

Template parameters
NRows number of rows in block as specified at compile-time
NCols number of columns in block as specified at compile-time
Parameters
startRow the first row in the block
startCol the first column in the block
blockRows number of rows in block as specified at run-time
blockCols number of columns in block as specified at run-time
Returns an expression of a block of *this.

This function is mainly useful for blocks where the number of rows is specified at compile-time and the number of columns is specified at run-time, or vice versa. The compile-time and run-time information should not contradict. In other words, blockRows should equal NRows unless NRows is Dynamic, and the same for the number of columns.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the block:" << endl << m.block<2, Dynamic>(1, 1, 2, 3) << endl;
m.block<2, Dynamic>(1, 1, 2, 3).setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is the block:
-6  1  0
-3  0  9
Now the matrix m is:
 7  9 -5 -3
-2  0  0  0
 6  0  0  0
 6  6  3  9

template<typename Derived> template<int NRows, int NCols>
const ConstFixedBlockXpr<NRows, NCols>::Type Eigen::DenseBase<Derived>::block(Index startRow, Index startCol, Index blockRows, Index blockCols) const

This is the const version of block<>(Index, Index, Index, Index).

template<typename Derived> template<typename NRowsType, typename NColsType>
FixedBlockXpr<...,...>::Type Eigen::DenseBase<Derived>::bottomLeftCorner(NRowsType cRows, NColsType cCols)

Template parameters
NRowsType the type of the value handling the number of rows in the block, typically Index.
NColsType the type of the value handling the number of columns in the block, typically Index.
Parameters
cRows the number of rows in the corner
cCols the number of columns in the corner
Returns an expression of a bottom-left corner of *this with either dynamic or fixed sizes.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.bottomLeftCorner(2, 2):" << endl;
cout << m.bottomLeftCorner(2, 2) << endl;
m.bottomLeftCorner(2, 2).setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is m.bottomLeftCorner(2, 2):
 6 -3
 6  6
Now the matrix m is:
 7  9 -5 -3
-2 -6  1  0
 0  0  0  9
 0  0  3  9

The number of rows blockRows and columns blockCols can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NRowsType, typename NColsType>
ConstFixedBlockXpr<...,...>::Type Eigen::DenseBase<Derived>::bottomLeftCorner(NRowsType cRows, NColsType cCols) const

This is the const version of bottomLeftCorner(NRowsType, NColsType).

template<typename Derived> template<int CRows, int CCols>
FixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::bottomLeftCorner()

Returns an expression of a fixed-size bottom-left corner of *this.

The template parameters CRows and CCols are the number of rows and columns in the corner.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.bottomLeftCorner<2,2>():" << endl;
cout << m.bottomLeftCorner<2,2>() << endl;
m.bottomLeftCorner<2,2>().setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

 

block(Index,Index,NRowsType,NColsType), class Block

template<typename Derived> template<int CRows, int CCols>
const ConstFixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::bottomLeftCorner() const

This is the const version of bottomLeftCorner<int, int>().

template<typename Derived> template<int CRows, int CCols>
FixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::bottomLeftCorner(Index cRows, Index cCols)

Template parameters
CRows number of rows in corner as specified at compile-time
CCols number of columns in corner as specified at compile-time
Parameters
cRows number of rows in corner as specified at run-time
cCols number of columns in corner as specified at run-time
Returns an expression of a bottom-left corner of *this.

This function is mainly useful for corners where the number of rows is specified at compile-time and the number of columns is specified at run-time, or vice versa. The compile-time and run-time information should not contradict. In other words, cRows should equal CRows unless CRows is Dynamic, and the same for the number of columns.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.bottomLeftCorner<2,Dynamic>(2,2):" << endl;
cout << m.bottomLeftCorner<2,Dynamic>(2,2) << endl;
m.bottomLeftCorner<2,Dynamic>(2,2).setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

class Block

template<typename Derived> template<int CRows, int CCols>
const ConstFixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::bottomLeftCorner(Index cRows, Index cCols) const

This is the const version of bottomLeftCorner<int, int>(Index, Index).

template<typename Derived> template<typename NRowsType, typename NColsType>
FixedBlockXpr<...,...>::Type Eigen::DenseBase<Derived>::bottomRightCorner(NRowsType cRows, NColsType cCols)

Template parameters
NRowsType the type of the value handling the number of rows in the block, typically Index.
NColsType the type of the value handling the number of columns in the block, typically Index.
Parameters
cRows the number of rows in the corner
cCols the number of columns in the corner
Returns an expression of a bottom-right corner of *this with either dynamic or fixed sizes.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.bottomRightCorner(2, 2):" << endl;
cout << m.bottomRightCorner(2, 2) << endl;
m.bottomRightCorner(2, 2).setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is m.bottomRightCorner(2, 2):
0 9
3 9
Now the matrix m is:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  0
 6  6  0  0

The number of rows blockRows and columns blockCols can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NRowsType, typename NColsType>
const ConstFixedBlockXpr<...,...>::Type Eigen::DenseBase<Derived>::bottomRightCorner(NRowsType cRows, NColsType cCols) const

This is the const version of bottomRightCorner(NRowsType, NColsType).

template<typename Derived> template<int CRows, int CCols>
FixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::bottomRightCorner()

Returns an expression of a fixed-size bottom-right corner of *this.

The template parameters CRows and CCols are the number of rows and columns in the corner.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.bottomRightCorner<2,2>():" << endl;
cout << m.bottomRightCorner<2,2>() << endl;
m.bottomRightCorner<2,2>().setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

 

block(Index,Index,NRowsType,NColsType), class Block

template<typename Derived> template<int CRows, int CCols>
const ConstFixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::bottomRightCorner() const

This is the const version of bottomRightCorner<int, int>().

template<typename Derived> template<int CRows, int CCols>
FixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::bottomRightCorner(Index cRows, Index cCols)

Template parameters
CRows number of rows in corner as specified at compile-time
CCols number of columns in corner as specified at compile-time
Parameters
cRows number of rows in corner as specified at run-time
cCols number of columns in corner as specified at run-time
Returns an expression of a bottom-right corner of *this.

This function is mainly useful for corners where the number of rows is specified at compile-time and the number of columns is specified at run-time, or vice versa. The compile-time and run-time information should not contradict. In other words, cRows should equal CRows unless CRows is Dynamic, and the same for the number of columns.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.bottomRightCorner<2,Dynamic>(2,2):" << endl;
cout << m.bottomRightCorner<2,Dynamic>(2,2) << endl;
m.bottomRightCorner<2,Dynamic>(2,2).setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

class Block

template<typename Derived> template<int CRows, int CCols>
const ConstFixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::bottomRightCorner(Index cRows, Index cCols) const

This is the const version of bottomRightCorner<int, int>(Index, Index).

template<typename Derived> template<typename NRowsType>
NRowsBlockXpr<...>::Type Eigen::DenseBase<Derived>::bottomRows(NRowsType n)

Template parameters
NRowsType the type of the value handling the number of rows in the block, typically Index.
Parameters
n the number of rows in the block
Returns a block consisting of the bottom rows of *this.

Example:

Array44i a = Array44i::Random();
cout << "Here is the array a:" << endl << a << endl;
cout << "Here is a.bottomRows(2):" << endl;
cout << a.bottomRows(2) << endl;
a.bottomRows(2).setZero();
cout << "Now the array a is:" << endl << a << endl;

Output:

Here is the array a:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is a.bottomRows(2):
 6 -3  0  9
 6  6  3  9
Now the array a is:
 7  9 -5 -3
-2 -6  1  0
 0  0  0  0
 0  0  0  0

The number of rows n can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NRowsType>
const ConstNRowsBlockXpr<...>::Type Eigen::DenseBase<Derived>::bottomRows(NRowsType n) const

This is the const version of bottomRows(NRowsType).

template<typename Derived> template<int N>
NRowsBlockXpr<N>::Type Eigen::DenseBase<Derived>::bottomRows(Index n = N)

Template parameters
N the number of rows in the block as specified at compile-time
Parameters
n the number of rows in the block as specified at run-time
Returns a block consisting of the bottom rows of *this.

The compile-time and run-time information should not contradict. In other words, n should equal N unless N is Dynamic.

Example:

Array44i a = Array44i::Random();
cout << "Here is the array a:" << endl << a << endl;
cout << "Here is a.bottomRows<2>():" << endl;
cout << a.bottomRows<2>() << endl;
a.bottomRows<2>().setZero();
cout << "Now the array a is:" << endl << a << endl;

Output:

 

block(Index,Index,NRowsType,NColsType), class Block

template<typename Derived> template<int N>
ConstNRowsBlockXpr<N>::Type Eigen::DenseBase<Derived>::bottomRows(Index n = N) const

This is the const version of bottomRows<int>().

template<typename Derived> template<typename NewType>
CastXpr<NewType>::Type Eigen::DenseBase<Derived>::cast() const

Returns an expression of *this with the Scalar type casted to NewScalar.

The template parameter NewScalar is the type we are casting the scalars to.

template<typename Derived>
const_iterator Eigen::DenseBase<Derived>::cbegin() const

returns a read-only const_iterator to the first element of the 1D vector or array This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

template<typename Derived>
const_iterator Eigen::DenseBase<Derived>::cend() const

returns a read-only const_iterator to the element following the last element of the 1D vector or array This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

template<typename Derived>
ColXpr Eigen::DenseBase<Derived>::col(Index i)

Returns an expression of the i-th column of *this. Note that the numbering starts at 0.

Example:

Matrix3d m = Matrix3d::Identity();
m.col(1) = Vector3d(4,5,6);
cout << m << endl;

Output:

 

row(), class Block

template<typename Derived>
ConstColXpr Eigen::DenseBase<Derived>::col(Index i) const

This is the const version of col().

template<typename Derived>
ConstColwiseReturnType Eigen::DenseBase<Derived>::colwise() const

Returns a VectorwiseOp wrapper of *this broadcasting and partial reductions

Example:

Matrix3d m = Matrix3d::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the sum of each column:" << endl << m.colwise().sum() << endl;
cout << "Here is the maximum absolute value of each column:"
     << endl << m.cwiseAbs().colwise().maxCoeff() << endl;

Output:

Here is the matrix m:
  0.68  0.597  -0.33
-0.211  0.823  0.536
 0.566 -0.605 -0.444
Here is the sum of each column:
  1.04  0.815 -0.238
Here is the maximum absolute value of each column:
 0.68 0.823 0.536

template<typename Derived>
ColwiseReturnType Eigen::DenseBase<Derived>::colwise()

Returns a writable VectorwiseOp wrapper of *this providing additional partial reduction operations

template<typename Derived>
ConjugateReturnType Eigen::DenseBase<Derived>::conjugate() const

Returns an expression of the complex conjugate of *this.

template<typename Derived> template<bool Cond>
internal::conditional<Cond, ConjugateReturnType, const Derived&>::type Eigen::DenseBase<Derived>::conjugateIf() const

Returns an expression of the complex conjugate of *this if Cond==true, returns derived() otherwise.

template<typename Derived>
Index Eigen::DenseBase<Derived>::count() const

Returns the number of coefficients which evaluate to true

template<typename Derived>
iterator Eigen::DenseBase<Derived>::end()

returns an iterator to the element following the last element of the 1D vector or array This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

template<typename Derived>
const_iterator Eigen::DenseBase<Derived>::end() const

const version of end()

template<typename Derived>
EvalReturnType Eigen::DenseBase<Derived>::eval() const

Returns the matrix or vector obtained by evaluating this expression.

Notice that in the case of a plain matrix or vector (not an expression) this function just returns a const reference, in order to avoid a useless copy.

template<typename Derived>
void Eigen::DenseBase<Derived>::fill(const Scalar& value)

Alias for setConstant(): sets all coefficients in this expression to val.

template<typename Derived> template<unsigned int Added, unsigned int Removed>
EIGEN_DEPRECATED const Derived& Eigen::DenseBase<Derived>::flagged() const

template<typename Derived>
const WithFormat<Derived> Eigen::DenseBase<Derived>::format(const IOFormat& fmt) const

Returns a WithFormat proxy object allowing to print a matrix the with given format fmt.

See class IOFormat for some examples.

template<typename Derived>
bool Eigen::DenseBase<Derived>::hasNaN() const

Returns true is *this contains at least one Not A Number (NaN).

template<typename Derived> template<typename NType>
FixedSegmentReturnType<...>::Type Eigen::DenseBase<Derived>::head(NType n)

Template parameters
NType the type of the value handling the number of coefficients in the segment, typically Index.
Parameters
n the number of coefficients in the segment
Returns an expression of the first coefficients of *this with either dynamic or fixed sizes.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

Example:

RowVector4i v = RowVector4i::Random();
cout << "Here is the vector v:" << endl << v << endl;
cout << "Here is v.head(2):" << endl << v.head(2) << endl;
v.head(2).setZero();
cout << "Now the vector v is:" << endl << v << endl;

Output:

Here is the vector v:
 7 -2  6  6
Here is v.head(2):
 7 -2
Now the vector v is:
0 0 6 6

The number of coefficients n can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NType>
const ConstFixedSegmentReturnType<...>::Type Eigen::DenseBase<Derived>::head(NType n) const

This is the const version of head(NType).

template<typename Derived> template<int N>
FixedSegmentReturnType<N>::Type Eigen::DenseBase<Derived>::head(Index n = N)

Template parameters
N the number of coefficients in the segment as specified at compile-time
Parameters
n the number of coefficients in the segment as specified at run-time
Returns a fixed-size expression of the first coefficients of *this.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

The compile-time and run-time information should not contradict. In other words, n should equal N unless N is Dynamic.

Example:

RowVector4i v = RowVector4i::Random();
cout << "Here is the vector v:" << endl << v << endl;
cout << "Here is v.head(2):" << endl << v.head<2>() << endl;
v.head<2>().setZero();
cout << "Now the vector v is:" << endl << v << endl;

Output:

Here is the vector v:
 7 -2  6  6
Here is v.head(2):
 7 -2
Now the vector v is:
0 0 6 6

template<typename Derived> template<int N>
ConstFixedSegmentReturnType<N>::Type Eigen::DenseBase<Derived>::head(Index n = N) const

This is the const version of head<int>().

template<typename Derived>
const ImagReturnType Eigen::DenseBase<Derived>::imag() const

Returns an read-only expression of the imaginary part of *this.

template<typename Derived>
NonConstImagReturnType Eigen::DenseBase<Derived>::imag()

Returns a non const expression of the imaginary part of *this.

template<typename Derived>
Index Eigen::DenseBase<Derived>::innerSize() const

Returns the inner size.

template<typename Derived>
InnerVectorReturnType Eigen::DenseBase<Derived>::innerVector(Index outer)

Returns the outer -th column (resp. row) of the matrix *this if *this is col-major (resp. row-major).

template<typename Derived>
const ConstInnerVectorReturnType Eigen::DenseBase<Derived>::innerVector(Index outer) const

Returns the outer -th column (resp. row) of the matrix *this if *this is col-major (resp. row-major). Read-only.

template<typename Derived>
InnerVectorsReturnType Eigen::DenseBase<Derived>::innerVectors(Index outerStart, Index outerSize)

Returns the outer -th column (resp. row) of the matrix *this if *this is col-major (resp. row-major).

template<typename Derived>
const ConstInnerVectorsReturnType Eigen::DenseBase<Derived>::innerVectors(Index outerStart, Index outerSize) const

Returns the outer -th column (resp. row) of the matrix *this if *this is col-major (resp. row-major). Read-only.

template<typename Derived> template<typename OtherDerived>
bool Eigen::DenseBase<Derived>::isApprox(const DenseBase<OtherDerived>& other, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const

Returns true if *this is approximately equal to other, within the precision determined by prec.

template<typename Derived>
bool Eigen::DenseBase<Derived>::isApproxToConstant(const Scalar& value, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const

Returns true if all coefficients in this matrix are approximately equal to val, to within precision prec

template<typename Derived>
bool Eigen::DenseBase<Derived>::isConstant(const Scalar& value, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const

Returns true if all coefficients in this matrix are approximately equal to value, to within precision prec

This is just an alias for isApproxToConstant().

template<typename Derived> template<typename Derived>
bool Eigen::DenseBase<Derived>::isMuchSmallerThan(const typename NumTraits<Scalar>::Real& other, const RealScalar& prec) const

Returns true if the norm of *this is much smaller than other, within the precision determined by prec.

For matrices, the comparison is done using the Hilbert-Schmidt norm. For this reason, the value of the reference scalar other should come from the Hilbert-Schmidt norm of a reference matrix of same dimensions.

template<typename Derived> template<typename OtherDerived>
bool Eigen::DenseBase<Derived>::isMuchSmallerThan(const DenseBase<OtherDerived>& other, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const

Returns true if the norm of *this is much smaller than the norm of other, within the precision determined by prec.

template<typename Derived>
bool Eigen::DenseBase<Derived>::isOnes(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const

Returns true if *this is approximately equal to the matrix where all coefficients are equal to 1, within the precision given by prec.

Example:

Matrix3d m = Matrix3d::Ones();
m(0,2) += 1e-4;
cout << "Here's the matrix m:" << endl << m << endl;
cout << "m.isOnes() returns: " << m.isOnes() << endl;
cout << "m.isOnes(1e-3) returns: " << m.isOnes(1e-3) << endl;

Output:

Here's the matrix m:
1 1 1
1 1 1
1 1 1
m.isOnes() returns: 0
m.isOnes(1e-3) returns: 1

template<typename Derived>
bool Eigen::DenseBase<Derived>::isZero(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const

Returns true if *this is approximately equal to the zero matrix, within the precision given by prec.

Example:

Matrix3d m = Matrix3d::Zero();
m(0,2) = 1e-4;
cout << "Here's the matrix m:" << endl << m << endl;
cout << "m.isZero() returns: " << m.isZero() << endl;
cout << "m.isZero(1e-3) returns: " << m.isZero(1e-3) << endl;

Output:

Here's the matrix m:
     0      0 0.0001
     0      0      0
     0      0      0
m.isZero() returns: 0
m.isZero(1e-3) returns: 1

template<typename Derived> template<typename NColsType>
NColsBlockXpr<...>::Type Eigen::DenseBase<Derived>::leftCols(NColsType n)

Template parameters
NColsType the type of the value handling the number of columns in the block, typically Index.
Parameters
n the number of columns in the block
Returns a block consisting of the left columns of *this.

Example:

Array44i a = Array44i::Random();
cout << "Here is the array a:" << endl << a << endl;
cout << "Here is a.leftCols(2):" << endl;
cout << a.leftCols(2) << endl;
a.leftCols(2).setZero();
cout << "Now the array a is:" << endl << a << endl;

Output:

Here is the array a:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is a.leftCols(2):
 7  9
-2 -6
 6 -3
 6  6
Now the array a is:
 0  0 -5 -3
 0  0  1  0
 0  0  0  9
 0  0  3  9

The number of columns n can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NColsType>
const ConstNColsBlockXpr<...>::Type Eigen::DenseBase<Derived>::leftCols(NColsType n) const

This is the const version of leftCols(NColsType).

template<typename Derived> template<int N>
NColsBlockXpr<N>::Type Eigen::DenseBase<Derived>::leftCols(Index n = N)

Template parameters
N the number of columns in the block as specified at compile-time
Parameters
n the number of columns in the block as specified at run-time
Returns a block consisting of the left columns of *this.

The compile-time and run-time information should not contradict. In other words, n should equal N unless N is Dynamic.

Example:

Array44i a = Array44i::Random();
cout << "Here is the array a:" << endl << a << endl;
cout << "Here is a.leftCols<2>():" << endl;
cout << a.leftCols<2>() << endl;
a.leftCols<2>().setZero();
cout << "Now the array a is:" << endl << a << endl;

Output:

 

block(Index,Index,NRowsType,NColsType), class Block

template<typename Derived> template<int N>
ConstNColsBlockXpr<N>::Type Eigen::DenseBase<Derived>::leftCols(Index n = N) const

This is the const version of leftCols<int>().

template<typename Derived>
internal::traits<Derived>::Scalar Eigen::DenseBase<Derived>::maxCoeff() const

Returns the maximum of all coefficients of *this.

template<typename Derived> template<typename IndexType>
internal::traits<Derived>::Scalar Eigen::DenseBase<Derived>::maxCoeff(IndexType* row, IndexType* col) const

Returns the maximum of all coefficients of *this and puts in *row and *col its location.

template<typename Derived> template<typename IndexType>
internal::traits<Derived>::Scalar Eigen::DenseBase<Derived>::maxCoeff(IndexType* index) const

Returns the maximum of all coefficients of *this and puts in *index its location.

template<typename Derived>
Scalar Eigen::DenseBase<Derived>::mean() const

Returns the mean of all coefficients of *this

template<typename Derived> template<typename NColsType>
NColsBlockXpr<...>::Type Eigen::DenseBase<Derived>::middleCols(Index startCol, NColsType numCols)

Template parameters
NColsType the type of the value handling the number of columns in the block, typically Index.
Parameters
startCol the index of the first column in the block
numCols the number of columns in the block
Returns a block consisting of a range of columns of *this.

Example:

#include <Eigen/Core>
#include <iostream>

using namespace Eigen;
using namespace std;

int main(void)
{
    int const N = 5;
    MatrixXi A(N,N);
    A.setRandom();
    cout << "A =\n" << A << '\n' << endl;
    cout << "A(1..3,:) =\n" << A.middleCols(1,3) << endl;
    return 0;
}

Output:

A =
  7  -6   0   9 -10
 -2  -3   3   3  -5
  6   6  -3   5  -8
  6  -5   0  -8   6
  9   1   9   2  -7

A(1..3,:) =
-6  0  9
-3  3  3
 6 -3  5
-5  0 -8
 1  9  2

The number of columns n can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NColsType>
const ConstNColsBlockXpr<...>::Type Eigen::DenseBase<Derived>::middleCols(Index startCol, NColsType numCols) const

This is the const version of middleCols(Index,NColsType).

template<typename Derived> template<int N>
NColsBlockXpr<N>::Type Eigen::DenseBase<Derived>::middleCols(Index startCol, Index n = N)

Template parameters
N the number of columns in the block as specified at compile-time
Parameters
startCol the index of the first column in the block
n the number of columns in the block as specified at run-time
Returns a block consisting of a range of columns of *this.

The compile-time and run-time information should not contradict. In other words, n should equal N unless N is Dynamic.

Example:

#include <Eigen/Core>
#include <iostream>

using namespace Eigen;
using namespace std;

int main(void)
{
    int const N = 5;
    MatrixXi A(N,N);
    A.setRandom();
    cout << "A =\n" << A << '\n' << endl;
    cout << "A(:,1..3) =\n" << A.middleCols<3>(1) << endl;
    return 0;
}

Output:

 

block(Index,Index,NRowsType,NColsType), class Block

template<typename Derived> template<int N>
ConstNColsBlockXpr<N>::Type Eigen::DenseBase<Derived>::middleCols(Index startCol, Index n = N) const

This is the const version of middleCols<int>().

template<typename Derived> template<typename NRowsType>
NRowsBlockXpr<...>::Type Eigen::DenseBase<Derived>::middleRows(Index startRow, NRowsType n)

Template parameters
NRowsType the type of the value handling the number of rows in the block, typically Index.
Parameters
startRow the index of the first row in the block
n the number of rows in the block
Returns a block consisting of a range of rows of *this.

Example:

#include <Eigen/Core>
#include <iostream>

using namespace Eigen;
using namespace std;

int main(void)
{
    int const N = 5;
    MatrixXi A(N,N);
    A.setRandom();
    cout << "A =\n" << A << '\n' << endl;
    cout << "A(2..3,:) =\n" << A.middleRows(2,2) << endl;
    return 0;
}

Output:

A =
  7  -6   0   9 -10
 -2  -3   3   3  -5
  6   6  -3   5  -8
  6  -5   0  -8   6
  9   1   9   2  -7

A(2..3,:) =
 6  6 -3  5 -8
 6 -5  0 -8  6

The number of rows n can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NRowsType>
const ConstNRowsBlockXpr<...>::Type Eigen::DenseBase<Derived>::middleRows(Index startRow, NRowsType n) const

This is the const version of middleRows(Index,NRowsType).

template<typename Derived> template<int N>
NRowsBlockXpr<N>::Type Eigen::DenseBase<Derived>::middleRows(Index startRow, Index n = N)

Template parameters
N the number of rows in the block as specified at compile-time
Parameters
startRow the index of the first row in the block
n the number of rows in the block as specified at run-time
Returns a block consisting of a range of rows of *this.

The compile-time and run-time information should not contradict. In other words, n should equal N unless N is Dynamic.

Example:

#include <Eigen/Core>
#include <iostream>

using namespace Eigen;
using namespace std;

int main(void)
{
    int const N = 5;
    MatrixXi A(N,N);
    A.setRandom();
    cout << "A =\n" << A << '\n' << endl;
    cout << "A(1..3,:) =\n" << A.middleRows<3>(1) << endl;
    return 0;
}

Output:

 

block(Index,Index,NRowsType,NColsType), class Block

template<typename Derived> template<int N>
ConstNRowsBlockXpr<N>::Type Eigen::DenseBase<Derived>::middleRows(Index startRow, Index n = N) const

This is the const version of middleRows<int>().

template<typename Derived>
internal::traits<Derived>::Scalar Eigen::DenseBase<Derived>::minCoeff() const

Returns the minimum of all coefficients of *this.

template<typename Derived> template<typename IndexType>
internal::traits<Derived>::Scalar Eigen::DenseBase<Derived>::minCoeff(IndexType* row, IndexType* col) const

Returns the minimum of all coefficients of *this and puts in *row and *col its location.

template<typename Derived> template<typename IndexType>
internal::traits<Derived>::Scalar Eigen::DenseBase<Derived>::minCoeff(IndexType* index) const

Returns the minimum of all coefficients of *this and puts in *index its location.

template<typename Derived>
const NestByValue<Derived> Eigen::DenseBase<Derived>::nestByValue() const

Returns an expression of the temporary version of *this.

template<typename Derived>
Index Eigen::DenseBase<Derived>::nonZeros() const

Returns the number of nonzero coefficients which is in practice the number of stored coefficients.

template<typename Derived> template<typename RowIndices, typename ColIndices>
IndexedView_or_Block Eigen::DenseBase<Derived>::operator()(const RowIndices& rowIndices, const ColIndices& colIndices)

Returns a generic submatrix view defined by the rows and columns indexed rowIndices and colIndices respectively.

Each parameter must either be:

  • An integer indexing a single row or column
  • Eigen::all indexing the full set of respective rows or columns in increasing order
  • An ArithmeticSequence as returned by the Eigen::seq and Eigen::seqN functions
  • Any Eigen's vector/array of integers or expressions
  • Plain C arrays: int[N]
  • And more generally any type exposing the following two member functions:

    <integral type> operator[](<integral type>) const;
    <integral type> size() const;

    where <integral type> stands for any integer type compatible with Eigen::Index (i.e. std::ptrdiff_t).

The last statement implies compatibility with std::vector, std::valarray, std::array, many of the Range-v3's ranges, etc.

If the submatrix can be represented using a starting position (i,j) and positive sizes (rows,columns), then this method will returns a Block object after extraction of the relevant information from the passed arguments. This is the case when all arguments are either:

  • An integer
  • Eigen::all
  • An ArithmeticSequence with compile-time increment strictly equal to 1, as returned by Eigen::seq(a,b), and Eigen::seqN(a,N).

Otherwise a more general IndexedView<Derived,RowIndices',ColIndices'> object will be returned, after conversion of the inputs to more suitable types RowIndices' and ColIndices'.

For 1D vectors and arrays, you better use the operator()(const Indices&) overload, which behave the same way but taking a single parameter.

See also this question and its answer for an example of how to duplicate coefficients.

template<typename Derived> template<typename Indices>
IndexedView_or_VectorBlock Eigen::DenseBase<Derived>::operator()(const Indices& indices)

This is an overload of operator()(const RowIndices&, const ColIndices&) for 1D vectors or arrays

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

template<typename Derived>
const NegativeReturnType Eigen::DenseBase<Derived>::operator-() const

Returns an expression of the opposite of *this

template<typename Derived>
CommaInitializer<Derived> Eigen::DenseBase<Derived>::operator<<(const Scalar& s)

Convenient operator to set the coefficients of a matrix.

The coefficients must be provided in a row major order and exactly match the size of the matrix. Otherwise an assertion is raised.

Example:

Matrix3i m1;
m1 << 1, 2, 3,
      4, 5, 6,
      7, 8, 9;
cout << m1 << endl << endl;
Matrix3i m2 = Matrix3i::Identity();
m2.block(0,0, 2,2) << 10, 11, 12, 13;
cout << m2 << endl << endl;
Vector2i v1;
v1 << 14, 15;
m2 << v1.transpose(), 16,
      v1, m1.block(1,1,2,2);
cout << m2 << endl;

Output:

1 2 3
4 5 6
7 8 9

10 11  0
12 13  0
 0  0  1

14 15 16
14  5  6
15  8  9

template<typename Derived> template<typename OtherDerived>
CommaInitializer<Derived> Eigen::DenseBase<Derived>::operator<<(const DenseBase<OtherDerived>& other)

template<typename Derived> template<typename OtherDerived>
Derived& Eigen::DenseBase<Derived>::operator=(const DenseBase<OtherDerived>& other)

Returns a reference to *this.

Copies other into *this.

template<typename Derived>
Derived& Eigen::DenseBase<Derived>::operator=(const DenseBase& other)

Special case of the template operator=, in order to prevent the compiler from generating a default operator= (issue hit with g++ 4.1)

template<typename Derived> template<typename OtherDerived>
Derived& Eigen::DenseBase<Derived>::operator=(const EigenBase<OtherDerived>& other)

Copies the generic expression other into *this.

Returns a reference to *this.

The expression must provide a (templated) evalTo(Derived& dst) const function which does the actual job. In practice, this allows any user to write its own special matrix without having to modify MatrixBase

template<typename Derived>
Index Eigen::DenseBase<Derived>::outerSize() const

Returns the outer size.

template<typename Derived>
Scalar Eigen::DenseBase<Derived>::prod() const

Returns the product of all coefficients of *this

Example:

Matrix3d m = Matrix3d::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the product of all the coefficients:" << endl << m.prod() << endl;

Output:

Here is the matrix m:
  0.68  0.597  -0.33
-0.211  0.823  0.536
 0.566 -0.605 -0.444
Here is the product of all the coefficients:
0.0019

template<typename Derived>
RealReturnType Eigen::DenseBase<Derived>::real() const

Returns a read-only expression of the real part of *this.

template<typename Derived>
NonConstRealReturnType Eigen::DenseBase<Derived>::real()

Returns a non const expression of the real part of *this.

template<typename Derived> template<typename Func>
internal::traits<Derived>::Scalar Eigen::DenseBase<Derived>::redux(const Func& func) const

Returns the result of a full redux operation on the whole matrix or vector using func

The template parameter BinaryOp is the type of the functor func which must be an associative operator. Both current C++98 and C++11 functor styles are handled.

template<typename Derived> template<int RowFactor, int ColFactor>
const Replicate<Derived, RowFactor, ColFactor> Eigen::DenseBase<Derived>::replicate() const

Returns an expression of the replication of *this

Example:

MatrixXi m = MatrixXi::Random(2,3);
cout << "Here is the matrix m:" << endl << m << endl;
cout << "m.replicate<3,2>() = ..." << endl;
cout << m.replicate<3,2>() << endl;

Output:

Here is the matrix m:
 7  6  9
-2  6 -6
m.replicate<3,2>() = ...
 7  6  9  7  6  9
-2  6 -6 -2  6 -6
 7  6  9  7  6  9
-2  6 -6 -2  6 -6
 7  6  9  7  6  9
-2  6 -6 -2  6 -6

template<typename Derived>
const Replicate<Derived, Dynamic, Dynamic> Eigen::DenseBase<Derived>::replicate(Index rowFactor, Index colFactor) const

Returns an expression of the replication of *this

Example:

Vector3i v = Vector3i::Random();
cout << "Here is the vector v:" << endl << v << endl;
cout << "v.replicate(2,5) = ..." << endl;
cout << v.replicate(2,5) << endl;

Output:

Here is the vector v:
 7
-2
 6
v.replicate(2,5) = ...
 7  7  7  7  7
-2 -2 -2 -2 -2
 6  6  6  6  6
 7  7  7  7  7
-2 -2 -2 -2 -2
 6  6  6  6  6

template<typename Derived> template<int Order = ColMajor, typename NRowsType, typename NColsType>
Reshaped<Derived,...> Eigen::DenseBase<Derived>::reshaped(NRowsType nRows, NColsType nCols)

Template parameters
Order specifies whether the coefficients should be processed in column-major-order (ColMajor), in row-major-order (RowMajor), or follows the natural order of the nested expression (AutoOrder). The default is ColMajor.
NRowsType the type of the value handling the number of rows, typically Index.
NColsType the type of the value handling the number of columns, typically Index.
Parameters
nRows the number of rows in the reshaped expression, specified at either run-time or compile-time, or AutoSize
nCols the number of columns in the reshaped expression, specified at either run-time or compile-time, or AutoSize
Returns an expression of *this with reshaped sizes.

Dynamic size example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.reshaped(2, 8):" << endl << m.reshaped(2, 8) << endl;

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is m.reshaped(2, 8):
 7  6  9 -3 -5  0 -3  9
-2  6 -6  6  1  3  0  9

The number of rows nRows and columns nCols can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. In the later case, n plays the role of a runtime fallback value in case N equals Eigen::Dynamic. Here is an example with a fixed number of rows and columns:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.reshaped(fix<2>,fix<8>):" << endl << m.reshaped(fix<2>,fix<8>) << endl;

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is m.reshaped(fix<2>,fix<8>):
 7  6  9 -3 -5  0 -3  9
-2  6 -6  6  1  3  0  9

Finally, one of the sizes parameter can be automatically deduced from the other one by passing AutoSize as in the following example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.reshaped(2, AutoSize):" << endl << m.reshaped(2, AutoSize) << endl;
cout << "Here is m.reshaped<RowMajor>(AutoSize, fix<8>):" << endl << m.reshaped<RowMajor>(AutoSize, fix<8>) << endl;

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is m.reshaped(2, AutoSize):
 7  6  9 -3 -5  0 -3  9
-2  6 -6  6  1  3  0  9
Here is m.reshaped<RowMajor>(AutoSize, fix<8>):
 7  9 -5 -3 -2 -6  1  0
 6 -3  0  9  6  6  3  9

AutoSize does preserve compile-time sizes when possible, i.e., when the sizes of the input are known at compile time and that the other size is passed at compile-time using Eigen::fix<N> as above.

template<typename Derived> template<int Order = ColMajor, typename NRowsType, typename NColsType>
const Reshaped<const Derived,...> Eigen::DenseBase<Derived>::reshaped(NRowsType nRows, NColsType nCols) const

This is the const version of reshaped(NRowsType,NColsType).

template<typename Derived> template<int Order = ColMajor>
Reshaped<Derived,...> Eigen::DenseBase<Derived>::reshaped()

Template parameters
Order specifies whether the coefficients should be processed in column-major-order (ColMajor), in row-major-order (RowMajor), or follows the natural order of the nested expression (AutoOrder). The default is ColMajor.
Returns an expression of *this with columns (or rows) stacked to a linear column vector

This overloads is essentially a shortcut for A.reshaped<Order>(AutoSize,fix<1>).

  • If Order==ColMajor (the default), then it returns a column-vector from the stacked columns of *this.
  • If Order==RowMajor, then it returns a column-vector from the stacked rows of *this.
  • If Order==AutoOrder, then it returns a column-vector with elements stacked following the storage order of *this. This mode is the recommended one when the particular ordering of the element is not relevant.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.reshaped().transpose():" << endl << m.reshaped().transpose() << endl;
cout << "Here is m.reshaped<RowMajor>().transpose():  " << endl << m.reshaped<RowMajor>().transpose() << endl;

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is m.reshaped().transpose():
 7 -2  6  6  9 -6 -3  6 -5  1  0  3 -3  0  9  9
Here is m.reshaped<RowMajor>().transpose():  
 7  9 -5 -3 -2 -6  1  0  6 -3  0  9  6  6  3  9

If you want more control, you can still fall back to reshaped(NRowsType,NColsType).

template<typename Derived> template<int Order = ColMajor>
const Reshaped<const Derived,...> Eigen::DenseBase<Derived>::reshaped() const

This is the const version of reshaped().

template<typename Derived>
void Eigen::DenseBase<Derived>::resize(Index newSize)

Only plain matrices/arrays, not expressions, may be resized; therefore the only useful resize methods are Matrix::resize() and Array::resize(). The present method only asserts that the new size equals the old size, and does nothing else.

template<typename Derived>
void Eigen::DenseBase<Derived>::resize(Index rows, Index cols)

Only plain matrices/arrays, not expressions, may be resized; therefore the only useful resize methods are Matrix::resize() and Array::resize(). The present method only asserts that the new size equals the old size, and does nothing else.

template<typename Derived>
ReverseReturnType Eigen::DenseBase<Derived>::reverse()

Returns an expression of the reverse of *this.

Example:

MatrixXi m = MatrixXi::Random(3,4);
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the reverse of m:" << endl << m.reverse() << endl;
cout << "Here is the coefficient (1,0) in the reverse of m:" << endl
     << m.reverse()(1,0) << endl;
cout << "Let us overwrite this coefficient with the value 4." << endl;
m.reverse()(1,0) = 4;
cout << "Now the matrix m is:" << endl << m << endl;

Output:

Here is the matrix m:
 7  6 -3  1
-2  9  6  0
 6 -6 -5  3
Here is the reverse of m:
 3 -5 -6  6
 0  6  9 -2
 1 -3  6  7
Here is the coefficient (1,0) in the reverse of m:
0
Let us overwrite this coefficient with the value 4.
Now the matrix m is:
 7  6 -3  1
-2  9  6  4
 6 -6 -5  3

template<typename Derived>
ConstReverseReturnType Eigen::DenseBase<Derived>::reverse() const

This is the const version of reverse().

template<typename Derived>
void Eigen::DenseBase<Derived>::reverseInPlace()

This is the "in place" version of reverse: it reverses *this.

In most cases it is probably better to simply use the reversed expression of a matrix. However, when reversing the matrix data itself is really needed, then this "in-place" version is probably the right choice because it provides the following additional benefits:

  • less error prone: doing the same operation with .reverse() requires special care: m = m.reverse().eval();
  • this API enables reverse operations without the need for a temporary
  • it allows future optimizations (cache friendliness, etc.)

template<typename Derived> template<typename NColsType>
NColsBlockXpr<...>::Type Eigen::DenseBase<Derived>::rightCols(NColsType n)

Template parameters
NColsType the type of the value handling the number of columns in the block, typically Index.
Parameters
n the number of columns in the block
Returns a block consisting of the right columns of *this.

Example:

Array44i a = Array44i::Random();
cout << "Here is the array a:" << endl << a << endl;
cout << "Here is a.rightCols(2):" << endl;
cout << a.rightCols(2) << endl;
a.rightCols(2).setZero();
cout << "Now the array a is:" << endl << a << endl;

Output:

Here is the array a:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is a.rightCols(2):
-5 -3
 1  0
 0  9
 3  9
Now the array a is:
 7  9  0  0
-2 -6  0  0
 6 -3  0  0
 6  6  0  0

The number of columns n can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NColsType>
const ConstNColsBlockXpr<...>::Type Eigen::DenseBase<Derived>::rightCols(NColsType n) const

This is the const version of rightCols(NColsType).

template<typename Derived> template<int N>
NColsBlockXpr<N>::Type Eigen::DenseBase<Derived>::rightCols(Index n = N)

Template parameters
N the number of columns in the block as specified at compile-time
Parameters
n the number of columns in the block as specified at run-time
Returns a block consisting of the right columns of *this.

The compile-time and run-time information should not contradict. In other words, n should equal N unless N is Dynamic.

Example:

Array44i a = Array44i::Random();
cout << "Here is the array a:" << endl << a << endl;
cout << "Here is a.rightCols<2>():" << endl;
cout << a.rightCols<2>() << endl;
a.rightCols<2>().setZero();
cout << "Now the array a is:" << endl << a << endl;

Output:

 

block(Index,Index,NRowsType,NColsType), class Block

template<typename Derived> template<int N>
ConstNColsBlockXpr<N>::Type Eigen::DenseBase<Derived>::rightCols(Index n = N) const

This is the const version of rightCols<int>().

template<typename Derived>
RowXpr Eigen::DenseBase<Derived>::row(Index i)

Returns an expression of the i-th row of *this. Note that the numbering starts at 0.

Example:

Matrix3d m = Matrix3d::Identity();
m.row(1) = Vector3d(4,5,6);
cout << m << endl;

Output:

 

col(), class Block

template<typename Derived>
ConstRowXpr Eigen::DenseBase<Derived>::row(Index i) const

This is the const version of row(). */.

template<typename Derived>
ConstRowwiseReturnType Eigen::DenseBase<Derived>::rowwise() const

Returns a VectorwiseOp wrapper of *this for broadcasting and partial reductions

Example:

Matrix3d m = Matrix3d::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the sum of each row:" << endl << m.rowwise().sum() << endl;
cout << "Here is the maximum absolute value of each row:"
     << endl << m.cwiseAbs().rowwise().maxCoeff() << endl;

Output:

Here is the matrix m:
  0.68  0.597  -0.33
-0.211  0.823  0.536
 0.566 -0.605 -0.444
Here is the sum of each row:
 0.948
  1.15
-0.483
Here is the maximum absolute value of each row:
 0.68
0.823
0.605

template<typename Derived>
RowwiseReturnType Eigen::DenseBase<Derived>::rowwise()

Returns a writable VectorwiseOp wrapper of *this providing additional partial reduction operations

template<typename Derived> template<typename NType>
FixedSegmentReturnType<...>::Type Eigen::DenseBase<Derived>::segment(Index start, NType n)

Template parameters
NType the type of the value handling the number of coefficients in the segment, typically Index.
Parameters
start the first coefficient in the segment
n the number of coefficients in the segment
Returns an expression of a segment (i.e. a vector block) in *this with either dynamic or fixed sizes.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

Example:

RowVector4i v = RowVector4i::Random();
cout << "Here is the vector v:" << endl << v << endl;
cout << "Here is v.segment(1, 2):" << endl << v.segment(1, 2) << endl;
v.segment(1, 2).setZero();
cout << "Now the vector v is:" << endl << v << endl;

Output:

Here is the vector v:
 7 -2  6  6
Here is v.segment(1, 2):
-2  6
Now the vector v is:
7 0 0 6

The number of coefficients n can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NType>
const ConstFixedSegmentReturnType<...>::Type Eigen::DenseBase<Derived>::segment(Index start, NType n) const

This is the const version of segment(Index,NType).

template<typename Derived> template<int N>
FixedSegmentReturnType<N>::Type Eigen::DenseBase<Derived>::segment(Index start, Index n = N)

Template parameters
N the number of coefficients in the segment as specified at compile-time
Parameters
start the index of the first element in the segment
n the number of coefficients in the segment as specified at compile-time
Returns a fixed-size expression of a segment (i.e. a vector block) in *this

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

The compile-time and run-time information should not contradict. In other words, n should equal N unless N is Dynamic.

Example:

RowVector4i v = RowVector4i::Random();
cout << "Here is the vector v:" << endl << v << endl;
cout << "Here is v.segment<2>(1):" << endl << v.segment<2>(1) << endl;
v.segment<2>(2).setZero();
cout << "Now the vector v is:" << endl << v << endl;

Output:

Here is the vector v:
 7 -2  6  6
Here is v.segment<2>(1):
-2  6
Now the vector v is:
 7 -2  0  0

template<typename Derived> template<int N>
ConstFixedSegmentReturnType<N>::Type Eigen::DenseBase<Derived>::segment(Index start, Index n = N) const

This is the const version of segment<int>(Index).

template<typename Derived> template<typename ThenDerived, typename ElseDerived>
const Select<Derived, ThenDerived, ElseDerived> Eigen::DenseBase<Derived>::select(const DenseBase<ThenDerived>& thenMatrix, const DenseBase<ElseDerived>& elseMatrix) const

Returns a matrix where each coefficient (i,j) is equal to thenMatrix(i,j) if *this(i,j), and elseMatrix(i,j) otherwise.

Example:

MatrixXi m(3, 3);
m << 1, 2, 3,
     4, 5, 6,
     7, 8, 9;
m = (m.array() >= 5).select(-m, m);
cout << m << endl;

Output:

 1  2  3
 4 -5 -6
-7 -8 -9

template<typename Derived> template<typename ThenDerived>
const Select<Derived, ThenDerived, typename ThenDerived::ConstantReturnType> Eigen::DenseBase<Derived>::select(const DenseBase<ThenDerived>& thenMatrix, const typename ThenDerived::Scalar& elseScalar) const

Version of DenseBase::select(const DenseBase&, const DenseBase&) with the else expression being a scalar value.

template<typename Derived> template<typename ElseDerived>
const Select<Derived, typename ElseDerived::ConstantReturnType, ElseDerived> Eigen::DenseBase<Derived>::select(const typename ElseDerived::Scalar& thenScalar, const DenseBase<ElseDerived>& elseMatrix) const

Version of DenseBase::select(const DenseBase&, const DenseBase&) with the then expression being a scalar value.

template<typename Derived>
Derived& Eigen::DenseBase<Derived>::setConstant(const Scalar& value)

Sets all coefficients in this expression to value val.

template<typename Derived>
Derived& Eigen::DenseBase<Derived>::setLinSpaced(Index size, const Scalar& low, const Scalar& high)

Sets a linearly spaced vector.

The function generates 'size' equally spaced values in the closed interval [low,high]. When size is set to 1, a vector of length 1 containing 'high' is returned.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

Example:

VectorXf v;
v.setLinSpaced(5,0.5f,1.5f);
cout << v << endl;

Output:

 0.5
0.75
   1
1.25
 1.5

For integer scalar types, do not miss the explanations on the definition of even spacing.

template<typename Derived>
Derived& Eigen::DenseBase<Derived>::setLinSpaced(const Scalar& low, const Scalar& high)

Sets a linearly spaced vector.

The function fills *this with equally spaced values in the closed interval [low,high]. When size is set to 1, a vector of length 1 containing 'high' is returned.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

For integer scalar types, do not miss the explanations on the definition of even spacing.

template<typename Derived>
Derived& Eigen::DenseBase<Derived>::setOnes()

Sets all coefficients in this expression to one.

Example:

Matrix4i m = Matrix4i::Random();
m.row(1).setOnes();
cout << m << endl;

Output:

 7  9 -5 -3
 1  1  1  1
 6 -3  0  9
 6  6  3  9

template<typename Derived>
Derived& Eigen::DenseBase<Derived>::setRandom()

Sets all coefficients in this expression to random values.

Numbers are uniformly spread through their whole definition range for integer types, and in the [-1:1] range for floating point scalar types.

Example:

Matrix4i m = Matrix4i::Zero();
m.col(1).setRandom();
cout << m << endl;

Output:

 0  7  0  0
 0 -2  0  0
 0  6  0  0
 0  6  0  0

template<typename Derived>
Derived& Eigen::DenseBase<Derived>::setZero()

Sets all coefficients in this expression to zero.

Example:

Matrix4i m = Matrix4i::Random();
m.row(1).setZero();
cout << m << endl;

Output:

 7  9 -5 -3
 0  0  0  0
 6 -3  0  9
 6  6  3  9

template<typename Derived> template<DirectionType Direction>
internal::conditional<Direction==Vertical, ColXpr, RowXpr>::type Eigen::DenseBase<Derived>::subVector(Index i)

Returns the i-th subvector (column or vector) according to the Direction

template<typename Derived> template<DirectionType Direction>
internal::conditional<Direction==Vertical, ConstColXpr, ConstRowXpr>::type Eigen::DenseBase<Derived>::subVector(Index i) const

This is the const version of subVector(Index)

template<typename Derived> template<DirectionType Direction>
Index Eigen::DenseBase<Derived>::subVectors() const

Returns the number of subvectors (rows or columns) in the direction Direction

template<typename Derived>
Scalar Eigen::DenseBase<Derived>::sum() const

Returns the sum of all coefficients of *this

If *this is empty, then the value 0 is returned.

template<typename Derived> template<typename OtherDerived>
void Eigen::DenseBase<Derived>::swap(const DenseBase<OtherDerived>& other)

swaps *this with the expression other.

template<typename Derived> template<typename OtherDerived>
void Eigen::DenseBase<Derived>::swap(PlainObjectBase<OtherDerived>& other)

swaps *this with the matrix or array other.

template<typename Derived> template<typename NType>
FixedSegmentReturnType<...>::Type Eigen::DenseBase<Derived>::tail(NType n)

Template parameters
NType the type of the value handling the number of coefficients in the segment, typically Index.
Parameters
n the number of coefficients in the segment
Returns an expression of a last coefficients of *this with either dynamic or fixed sizes.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

Example:

RowVector4i v = RowVector4i::Random();
cout << "Here is the vector v:" << endl << v << endl;
cout << "Here is v.tail(2):" << endl << v.tail(2) << endl;
v.tail(2).setZero();
cout << "Now the vector v is:" << endl << v << endl;

Output:

Here is the vector v:
 7 -2  6  6
Here is v.tail(2):
6 6
Now the vector v is:
 7 -2  0  0

The number of coefficients n can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NType>
const ConstFixedSegmentReturnType<...>::Type Eigen::DenseBase<Derived>::tail(NType n) const

This is the const version of tail(Index).

template<typename Derived> template<int N>
FixedSegmentReturnType<N>::Type Eigen::DenseBase<Derived>::tail(Index n = N)

Template parameters
N the number of coefficients in the segment as specified at compile-time
Parameters
n the number of coefficients in the segment as specified at run-time
Returns a fixed-size expression of the last coefficients of *this.

This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.

The compile-time and run-time information should not contradict. In other words, n should equal N unless N is Dynamic.

Example:

RowVector4i v = RowVector4i::Random();
cout << "Here is the vector v:" << endl << v << endl;
cout << "Here is v.tail(2):" << endl << v.tail<2>() << endl;
v.tail<2>().setZero();
cout << "Now the vector v is:" << endl << v << endl;

Output:

Here is the vector v:
 7 -2  6  6
Here is v.tail(2):
6 6
Now the vector v is:
 7 -2  0  0

template<typename Derived> template<int N>
ConstFixedSegmentReturnType<N>::Type Eigen::DenseBase<Derived>::tail(Index n = N) const

This is the const version of tail<int>.

template<typename Derived> template<typename NRowsType, typename NColsType>
FixedBlockXpr<...,...>::Type Eigen::DenseBase<Derived>::topLeftCorner(NRowsType cRows, NColsType cCols)

Template parameters
NRowsType the type of the value handling the number of rows in the block, typically Index.
NColsType the type of the value handling the number of columns in the block, typically Index.
Parameters
cRows the number of rows in the corner
cCols the number of columns in the corner
Returns an expression of a top-left corner of *this with either dynamic or fixed sizes.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.topLeftCorner(2, 2):" << endl;
cout << m.topLeftCorner(2, 2) << endl;
m.topLeftCorner(2, 2).setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is m.topLeftCorner(2, 2):
 7  9
-2 -6
Now the matrix m is:
 0  0 -5 -3
 0  0  1  0
 6 -3  0  9
 6  6  3  9

The number of rows blockRows and columns blockCols can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NRowsType, typename NColsType>
const ConstFixedBlockXpr<...,...>::Type Eigen::DenseBase<Derived>::topLeftCorner(NRowsType cRows, NColsType cCols) const

This is the const version of topLeftCorner(Index, Index).

template<typename Derived> template<int CRows, int CCols>
FixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::topLeftCorner()

Returns an expression of a fixed-size top-left corner of *this.

The template parameters CRows and CCols are the number of rows and columns in the corner.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.topLeftCorner<2,2>():" << endl;
cout << m.topLeftCorner<2,2>() << endl;
m.topLeftCorner<2,2>().setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

 

block(Index,Index,NRowsType,NColsType), class Block

template<typename Derived> template<int CRows, int CCols>
const ConstFixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::topLeftCorner() const

This is the const version of topLeftCorner<int, int>().

template<typename Derived> template<int CRows, int CCols>
FixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::topLeftCorner(Index cRows, Index cCols)

Template parameters
CRows number of rows in corner as specified at compile-time
CCols number of columns in corner as specified at compile-time
Parameters
cRows number of rows in corner as specified at run-time
cCols number of columns in corner as specified at run-time
Returns an expression of a top-left corner of *this.

This function is mainly useful for corners where the number of rows is specified at compile-time and the number of columns is specified at run-time, or vice versa. The compile-time and run-time information should not contradict. In other words, cRows should equal CRows unless CRows is Dynamic, and the same for the number of columns.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.topLeftCorner<2,Dynamic>(2,2):" << endl;
cout << m.topLeftCorner<2,Dynamic>(2,2) << endl;
m.topLeftCorner<2,Dynamic>(2,2).setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

class Block

template<typename Derived> template<int CRows, int CCols>
const ConstFixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::topLeftCorner(Index cRows, Index cCols) const

This is the const version of topLeftCorner<int, int>(Index, Index).

template<typename Derived> template<typename NRowsType, typename NColsType>
FixedBlockXpr<...,...>::Type Eigen::DenseBase<Derived>::topRightCorner(NRowsType cRows, NColsType cCols)

Template parameters
NRowsType the type of the value handling the number of rows in the block, typically Index.
NColsType the type of the value handling the number of columns in the block, typically Index.
Parameters
cRows the number of rows in the corner
cCols the number of columns in the corner
Returns a expression of a top-right corner of *this with either dynamic or fixed sizes.

Example with dynamic sizes:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.topRightCorner(2, 2):" << endl;
cout << m.topRightCorner(2, 2) << endl;
m.topRightCorner(2, 2).setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is m.topRightCorner(2, 2):
-5 -3
 1  0
Now the matrix m is:
 7  9  0  0
-2 -6  0  0
 6 -3  0  9
 6  6  3  9

The number of rows blockRows and columns blockCols can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NRowsType, typename NColsType>
const ConstFixedBlockXpr<...,...>::Type Eigen::DenseBase<Derived>::topRightCorner(NRowsType cRows, NColsType cCols) const

This is the const version of topRightCorner(NRowsType, NColsType).

template<typename Derived> template<int CRows, int CCols>
FixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::topRightCorner()

Template parameters
CRows the number of rows in the corner
CCols the number of columns in the corner
Returns an expression of a fixed-size top-right corner of *this.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.topRightCorner<2,2>():" << endl;
cout << m.topRightCorner<2,2>() << endl;
m.topRightCorner<2,2>().setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

class Block, block<int,int>(Index,Index)

template<typename Derived> template<int CRows, int CCols>
const ConstFixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::topRightCorner() const

This is the const version of topRightCorner<int, int>().

template<typename Derived> template<int CRows, int CCols>
FixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::topRightCorner(Index cRows, Index cCols)

Template parameters
CRows number of rows in corner as specified at compile-time
CCols number of columns in corner as specified at compile-time
Parameters
cRows number of rows in corner as specified at run-time
cCols number of columns in corner as specified at run-time
Returns an expression of a top-right corner of *this.

This function is mainly useful for corners where the number of rows is specified at compile-time and the number of columns is specified at run-time, or vice versa. The compile-time and run-time information should not contradict. In other words, cRows should equal CRows unless CRows is Dynamic, and the same for the number of columns.

Example:

Matrix4i m = Matrix4i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is m.topRightCorner<2,Dynamic>(2,2):" << endl;
cout << m.topRightCorner<2,Dynamic>(2,2) << endl;
m.topRightCorner<2,Dynamic>(2,2).setZero();
cout << "Now the matrix m is:" << endl << m << endl;

Output:

class Block

template<typename Derived> template<int CRows, int CCols>
const ConstFixedBlockXpr<CRows, CCols>::Type Eigen::DenseBase<Derived>::topRightCorner(Index cRows, Index cCols) const

This is the const version of topRightCorner<int, int>(Index, Index).

template<typename Derived> template<typename NRowsType>
NRowsBlockXpr<...>::Type Eigen::DenseBase<Derived>::topRows(NRowsType n)

Template parameters
NRowsType the type of the value handling the number of rows in the block, typically Index.
Parameters
n the number of rows in the block
Returns a block consisting of the top rows of *this.

Example:

Array44i a = Array44i::Random();
cout << "Here is the array a:" << endl << a << endl;
cout << "Here is a.topRows(2):" << endl;
cout << a.topRows(2) << endl;
a.topRows(2).setZero();
cout << "Now the array a is:" << endl << a << endl;

Output:

Here is the array a:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here is a.topRows(2):
 7  9 -5 -3
-2 -6  1  0
Now the array a is:
 0  0  0  0
 0  0  0  0
 6 -3  0  9
 6  6  3  9

The number of rows n can also be specified at compile-time by passing Eigen::fix<N>, or Eigen::fix<N>(n) as arguments. See block() for the details.

template<typename Derived> template<typename NRowsType>
const ConstNRowsBlockXpr<...>::Type Eigen::DenseBase<Derived>::topRows(NRowsType n) const

This is the const version of topRows(NRowsType).

template<typename Derived> template<int N>
NRowsBlockXpr<N>::Type Eigen::DenseBase<Derived>::topRows(Index n = N)

Template parameters
N the number of rows in the block as specified at compile-time
Parameters
n the number of rows in the block as specified at run-time
Returns a block consisting of the top rows of *this.

The compile-time and run-time information should not contradict. In other words, n should equal N unless N is Dynamic.

Example:

Array44i a = Array44i::Random();
cout << "Here is the array a:" << endl << a << endl;
cout << "Here is a.topRows<2>():" << endl;
cout << a.topRows<2>() << endl;
a.topRows<2>().setZero();
cout << "Now the array a is:" << endl << a << endl;

Output:

 

block(Index,Index,NRowsType,NColsType), class Block

template<typename Derived> template<int N>
ConstNRowsBlockXpr<N>::Type Eigen::DenseBase<Derived>::topRows(Index n = N) const

This is the const version of topRows<int>().

template<typename Derived>
TransposeReturnType Eigen::DenseBase<Derived>::transpose()

Returns an expression of the transpose of *this.

Example:

Matrix2i m = Matrix2i::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the transpose of m:" << endl << m.transpose() << endl;
cout << "Here is the coefficient (1,0) in the transpose of m:" << endl
     << m.transpose()(1,0) << endl;
cout << "Let us overwrite this coefficient with the value 0." << endl;
m.transpose()(1,0) = 0;
cout << "Now the matrix m is:" << endl << m << endl;

Output:

Here is the matrix m:
 7  6
-2  6
Here is the transpose of m:
 7 -2
 6  6
Here is the coefficient (1,0) in the transpose of m:
6
Let us overwrite this coefficient with the value 0.
Now the matrix m is:
 7  0
-2  6

template<typename Derived>
ConstTransposeReturnType Eigen::DenseBase<Derived>::transpose() const

This is the const version of transpose().

Make sure you read the warning for transpose() !

template<typename Derived>
void Eigen::DenseBase<Derived>::transposeInPlace()

This is the "in place" version of transpose(): it replaces *this by its own transpose. Thus, doing m.transposeInPlace(); has the same effect on m as doing m = m.transpose().eval(); and is faster and also safer because in the latter line of code, forgetting the eval() results in a bug caused by aliasing.

Notice however that this method is only useful if you want to replace a matrix by its own transpose. If you just need the transpose of a matrix, use transpose().

template<typename Derived> template<typename CustomUnaryOp>
const CwiseUnaryOp<CustomUnaryOp, const Derived> Eigen::DenseBase<Derived>::unaryExpr(const CustomUnaryOp& func = CustomUnaryOp()) const

Apply a unary operator coefficient-wise.

Template parameters
CustomUnaryOp Type of func
Parameters
func in Functor implementing the unary operator
Returns An expression of a custom coefficient-wise unary operator func of *this

The function ptr_fun() from the C++ standard library can be used to make functors out of normal functions.

Example:

#include <Eigen/Core>
#include <iostream>
using namespace Eigen;
using namespace std;

// define function to be applied coefficient-wise
double ramp(double x)
{
  if (x > 0)
    return x;
  else 
    return 0;
}

int main(int, char**)
{
  Matrix4d m1 = Matrix4d::Random();
  cout << m1 << endl << "becomes: " << endl << m1.unaryExpr(ptr_fun(ramp)) << endl;
  return 0;
}

Output:

   0.68   0.823  -0.444   -0.27
 -0.211  -0.605   0.108  0.0268
  0.566   -0.33 -0.0452   0.904
  0.597   0.536   0.258   0.832
becomes: 
  0.68  0.823      0      0
     0      0  0.108 0.0268
 0.566      0      0  0.904
 0.597  0.536  0.258  0.832

Genuine functors allow for more possibilities, for instance it may contain a state.

Example:

#include <Eigen/Core>
#include <iostream>
using namespace Eigen;
using namespace std;

// define a custom template unary functor
template<typename Scalar>
struct CwiseClampOp {
  CwiseClampOp(const Scalar& inf, const Scalar& sup) : m_inf(inf), m_sup(sup) {}
  const Scalar operator()(const Scalar& x) const { return x<m_inf ? m_inf : (x>m_sup ? m_sup : x); }
  Scalar m_inf, m_sup;
};

int main(int, char**)
{
  Matrix4d m1 = Matrix4d::Random();
  cout << m1 << endl << "becomes: " << endl << m1.unaryExpr(CwiseClampOp<double>(-0.5,0.5)) << endl;
  return 0;
}

Output:

unaryViewExpr, binaryExpr, class CwiseUnaryOp

template<typename Derived> template<typename CustomViewOp>
const CwiseUnaryView<CustomViewOp, const Derived> Eigen::DenseBase<Derived>::unaryViewExpr(const CustomViewOp& func = CustomViewOp()) const

Returns an expression of a custom coefficient-wise unary operator func of *this

The template parameter CustomUnaryOp is the type of the functor of the custom unary operator.

Example:

#include <Eigen/Core>
#include <iostream>
using namespace Eigen;
using namespace std;

// define a custom template unary functor
template<typename Scalar>
struct CwiseClampOp {
  CwiseClampOp(const Scalar& inf, const Scalar& sup) : m_inf(inf), m_sup(sup) {}
  const Scalar operator()(const Scalar& x) const { return x<m_inf ? m_inf : (x>m_sup ? m_sup : x); }
  Scalar m_inf, m_sup;
};

int main(int, char**)
{
  Matrix4d m1 = Matrix4d::Random();
  cout << m1 << endl << "becomes: " << endl << m1.unaryExpr(CwiseClampOp<double>(-0.5,0.5)) << endl;
  return 0;
}

Output:

unaryExpr, binaryExpr class CwiseUnaryOp

template<typename Derived>
CoeffReturnType Eigen::DenseBase<Derived>::value() const

Returns the unique coefficient of a 1x1 expression

template<typename Derived> template<typename Visitor>
void Eigen::DenseBase<Derived>::visit(Visitor& func) const

Applies the visitor visitor to the whole coefficients of the matrix or vector.

The template parameter Visitor is the type of the visitor and provides the following interface:

struct MyVisitor {
  // called for the first coefficient
  void init(const Scalar& value, Index i, Index j);
  // called for all other coefficients
  void operator() (const Scalar& value, Index i, Index j);
};

template<typename Derived> template<typename Derived>
std::ostream& operator<<(std::ostream& s, const DenseBase<Derived>& m)

Outputs the matrix, to the given stream.

If you wish to print the matrix with a format different than the default, use DenseBase::format().

It is also possible to change the default format by defining EIGEN_DEFAULT_IO_FORMAT before including Eigen headers. If not defined, this will automatically be defined to Eigen::IOFormat(), that is the Eigen::IOFormat with default parameters.