Eigen  3.4.90 (git rev 9589cc4e7fd8e4538bedef80dd36c7738977a8be)
 
Loading...
Searching...
No Matches
Eigen::DenseBase< Derived > Class Template Reference

#include <Eigen/src/Core/DenseBase.h>

Detailed Description

template<typename Derived>
class Eigen::DenseBase< Derived >

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

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.

Template Parameters
Derivedis the derived type, e.g., a matrix type or an expression.

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.

See also
The class hierarchy
+ Inheritance diagram for Eigen::DenseBase< Derived >:

Public Types

enum  {
  RowsAtCompileTime ,
  ColsAtCompileTime ,
  SizeAtCompileTime ,
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime ,
  MaxSizeAtCompileTime ,
  IsVectorAtCompileTime ,
  NumDimensions ,
  Flags ,
  IsRowMajor ,
  InnerSizeAtCompileTime ,
  InnerStrideAtCompileTime ,
  OuterStrideAtCompileTime
}
 
typedef random_access_iterator_type const_iterator
 
typedef Eigen::InnerIterator< Derived > InnerIterator
 
typedef random_access_iterator_type iterator
 
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 >::MaxColsAtCompileTimePlainArray
 
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 >::MaxColsAtCompileTimePlainMatrix
 
typedef std::conditional_t< internal::is_same< typename internal::traits< Derived >::XprKind, MatrixXpr >::value, PlainMatrix, PlainArrayPlainObject
 The plain matrix or array type corresponding to this expression.
 
typedef internal::traits< Derived >::Scalar Scalar
 
typedef internal::traits< Derived >::StorageIndex StorageIndex
 The type used to store indices.
 
typedef Scalar value_type
 

Public Member Functions

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

Static Public Member Functions

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

Protected Member Functions

constexpr DenseBase ()=default
 

Related Symbols

(Note that these are not member symbols.)

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

Member Typedef Documentation

◆ const_iterator

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

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

◆ InnerIterator

template<typename Derived>
typedef Eigen::InnerIterator<Derived> Eigen::DenseBase< Derived >::InnerIterator

Inner iterator type to iterate over the coefficients of a row or column.

See also
class InnerIterator

◆ iterator

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.

◆ PlainArray

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.

See also
PlainObject

◆ PlainMatrix

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.

See also
PlainObject

◆ PlainObject

template<typename Derived>
typedef std::conditional_t<internal::is_same<typename internal::traits<Derived>::XprKind, MatrixXpr>::value, PlainMatrix, PlainArray> 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&.

◆ Scalar

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.

◆ StorageIndex

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

See also
Preprocessor directives, Eigen::Index, SparseMatrixBase.

◆ value_type

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

Member Enumeration Documentation

◆ anonymous enum

template<typename Derived>
anonymous enum
Enumerator
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.

See also
MatrixBase::rows(), MatrixBase::cols(), ColsAtCompileTime, SizeAtCompileTime
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.

See also
MatrixBase::rows(), MatrixBase::cols(), RowsAtCompileTime, SizeAtCompileTime
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.

See also
RowsAtCompileTime, ColsAtCompileTime
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.

See also
RowsAtCompileTime, MaxColsAtCompileTime, MaxSizeAtCompileTime
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.

See also
ColsAtCompileTime, MaxRowsAtCompileTime, MaxSizeAtCompileTime
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.

See also
SizeAtCompileTime, MaxRowsAtCompileTime, MaxColsAtCompileTime
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.

Constructor & Destructor Documentation

◆ DenseBase()

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

Default constructor. Do nothing.

Member Function Documentation

◆ all()

template<typename Derived>
bool Eigen::DenseBase< Derived >::all ( ) const
inline
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;
Matrix< float, 3, 1 > Vector3f
3×1 vector of type float.
Definition Matrix.h:478

Output:

Is ( -0.824   0.924 -0.0532) inside the box: 0
Is (0.199 0.237 0.146) inside the box: 1
See also
any(), Cwise::operator<()

◆ allFinite()

template<typename Derived>
bool Eigen::DenseBase< Derived >::allFinite ( ) const
inline
Returns
true if *this contains only finite numbers, i.e., no NaN and no +/-INF values.
See also
hasNaN()

◆ any()

template<typename Derived>
bool Eigen::DenseBase< Derived >::any ( ) const
inline
Returns
true if at least one coefficient is true
See also
all()

◆ begin() [1/2]

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

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.

See also
end(), cbegin()

◆ begin() [2/2]

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

const version of begin()

◆ block() [1/3]

template<typename Derived>
template<int NRows, int NCols>
FixedBlockXpr< NRows, NCols >::Type Eigen::DenseBase< Derived >::block ( Index startRow,
Index startCol )
inline
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.

Parameters
startRowthe first row in the block
startColthe first column 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;
Derived & setZero()
Definition CwiseNullaryOp.h:554
Matrix< int, 4, 4 > Matrix4i
4×4 matrix of type int.
Definition Matrix.h:477

Output:

Here is the matrix m:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.block<2,2>(1,1):
-1122281286   304089172
-1364114958    35005211
Now the matrix m is:
 1804289383 -1550966999  1365180540   336465782
 -465790871           0           0 -1868760786
 -189735855           0           0    -2309581
  719885386  2044897763 -1852781081  1101513929
Note
The usage of of this overload is discouraged from Eigen 3.4, better used the generic block(Index,Index,NRowsType,NColsType), here is the one-to-one equivalence:
mat.template block<NRows,NCols>(i,j) <--> mat.block(i,j,fix<NRows>,fix<NCols>)
FixedBlockXpr<...,... >::Type block(Index startRow, Index startCol, NRowsType blockRows, NColsType blockCols)
Definition DenseBase.h:122
static const auto fix()
since block is a templated member, the keyword template has to be used if the matrix type is also a template parameter:
m.template block<3,3>(1,1);
See also
block(Index,Index,NRowsType,NColsType), class Block

◆ block() [2/3]

template<typename Derived>
template<int NRows, int NCols>
FixedBlockXpr< NRows, NCols >::Type Eigen::DenseBase< Derived >::block ( Index startRow,
Index startCol,
Index blockRows,
Index blockCols )
inline
Returns
an expression of a block of *this.
Template Parameters
NRowsnumber of rows in block as specified at compile-time
NColsnumber of columns in block as specified at compile-time
Parameters
startRowthe first row in the block
startColthe first column in the block
blockRowsnumber of rows in block as specified at run-time
blockColsnumber of columns in block as specified at run-time

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;
const int Dynamic
Definition Constants.h:25

Output:

Here is the matrix m:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is the block:
-1122281286   304089172 -1868760786
-1364114958    35005211    -2309581
Now the matrix m is:
 1804289383 -1550966999  1365180540   336465782
 -465790871           0           0           0
 -189735855           0           0           0
  719885386  2044897763 -1852781081  1101513929
Note
The usage of of this overload is discouraged from Eigen 3.4, better used the generic block(Index,Index,NRowsType,NColsType), here is the one-to-one complete equivalence:
mat.template block<NRows,NCols>(i,j,rows,cols) <--> mat.block(i,j,fix<NRows>(rows),fix<NCols>(cols))
If we known that, e.g., NRows==Dynamic and NCols!=Dynamic, then the equivalence becomes:
mat.template block<Dynamic,NCols>(i,j,rows,NCols) <--> mat.block(i,j,rows,fix<NCols>)
See also
block(Index,Index,NRowsType,NColsType), class Block

◆ block() [3/3]

template<typename Derived>
template<typename NRowsType, typename NColsType>
FixedBlockXpr<...,... >::Type Eigen::DenseBase< Derived >::block ( Index startRow,
Index startCol,
NRowsType blockRows,
NColsType blockCols )
inline
Returns
an expression of a block in *this with either dynamic or fixed sizes.
Parameters
startRowthe first row in the block
startColthe first column in the block
blockRowsnumber of rows in the block, specified at either run-time or compile-time
blockColsnumber of columns in the block, specified at either run-time or compile-time
Template Parameters
NRowsTypethe type of the value handling the number of rows in the block, typically Index.
NColsTypethe type of the value handling the number of columns in the block, typically Index.

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;
Derived & setZero(Index size)
Definition CwiseNullaryOp.h:569

Output:

Here is the matrix m:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.block(1, 1, 2, 2):
-1122281286   304089172
-1364114958    35005211
Now the matrix m is:
 1804289383 -1550966999  1365180540   336465782
 -465790871           0           0 -1868760786
 -189735855           0           0    -2309581
  719885386  2044897763 -1852781081  1101513929

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.

Note
Even in the case that the returned expression has dynamic size, in the case when it is applied to a fixed-size matrix, it inherits a fixed maximal size, which means that evaluating it does not cause a dynamic memory allocation.
See also
class Block, fix, fix<N>(int)

◆ bottomLeftCorner() [1/3]

template<typename Derived>
template<int CRows, int CCols>
FixedBlockXpr< CRows, CCols >::Type Eigen::DenseBase< Derived >::bottomLeftCorner ( )
inline
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:

Here is the matrix m:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.bottomLeftCorner<2,2>():
 -189735855 -1364114958
  719885386  2044897763
Now the matrix m is:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
          0           0    35005211    -2309581
          0           0 -1852781081  1101513929
See also
block(Index,Index,NRowsType,NColsType), class Block

◆ bottomLeftCorner() [2/3]

template<typename Derived>
template<int CRows, int CCols>
FixedBlockXpr< CRows, CCols >::Type Eigen::DenseBase< Derived >::bottomLeftCorner ( Index cRows,
Index cCols )
inline
Returns
an expression of a bottom-left corner of *this.
Template Parameters
CRowsnumber of rows in corner as specified at compile-time
CColsnumber of columns in corner as specified at compile-time
Parameters
cRowsnumber of rows in corner as specified at run-time
cColsnumber of columns in corner as specified at run-time

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:

Here is the matrix m:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.bottomLeftCorner<2,Dynamic>(2,2):
 -189735855 -1364114958
  719885386  2044897763
Now the matrix m is:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
          0           0    35005211    -2309581
          0           0 -1852781081  1101513929
See also
class Block

◆ bottomLeftCorner() [3/3]

template<typename Derived>
template<typename NRowsType, typename NColsType>
FixedBlockXpr<...,... >::Type Eigen::DenseBase< Derived >::bottomLeftCorner ( NRowsType cRows,
NColsType cCols )
inline
Returns
an expression of a bottom-left corner of *this with either dynamic or fixed sizes.
Parameters
cRowsthe number of rows in the corner
cColsthe number of columns in the corner
Template Parameters
NRowsTypethe type of the value handling the number of rows in the block, typically Index.
NColsTypethe type of the value handling the number of columns in the block, typically Index.

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:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.bottomLeftCorner(2, 2):
 -189735855 -1364114958
  719885386  2044897763
Now the matrix m is:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
          0           0    35005211    -2309581
          0           0 -1852781081  1101513929

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.

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

◆ bottomRightCorner() [1/3]

template<typename Derived>
template<int CRows, int CCols>
FixedBlockXpr< CRows, CCols >::Type Eigen::DenseBase< Derived >::bottomRightCorner ( )
inline
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:

Here is the matrix m:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.bottomRightCorner<2,2>():
   35005211    -2309581
-1852781081  1101513929
Now the matrix m is:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958           0           0
  719885386  2044897763           0           0
See also
block(Index,Index,NRowsType,NColsType), class Block

◆ bottomRightCorner() [2/3]

template<typename Derived>
template<int CRows, int CCols>
FixedBlockXpr< CRows, CCols >::Type Eigen::DenseBase< Derived >::bottomRightCorner ( Index cRows,
Index cCols )
inline
Returns
an expression of a bottom-right corner of *this.
Template Parameters
CRowsnumber of rows in corner as specified at compile-time
CColsnumber of columns in corner as specified at compile-time
Parameters
cRowsnumber of rows in corner as specified at run-time
cColsnumber of columns in corner as specified at run-time

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:

Here is the matrix m:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.bottomRightCorner<2,Dynamic>(2,2):
   35005211    -2309581
-1852781081  1101513929
Now the matrix m is:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958           0           0
  719885386  2044897763           0           0
See also
class Block

◆ bottomRightCorner() [3/3]

template<typename Derived>
template<typename NRowsType, typename NColsType>
FixedBlockXpr<...,... >::Type Eigen::DenseBase< Derived >::bottomRightCorner ( NRowsType cRows,
NColsType cCols )
inline
Returns
an expression of a bottom-right corner of *this with either dynamic or fixed sizes.
Parameters
cRowsthe number of rows in the corner
cColsthe number of columns in the corner
Template Parameters
NRowsTypethe type of the value handling the number of rows in the block, typically Index.
NColsTypethe type of the value handling the number of columns in the block, typically Index.

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:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.bottomRightCorner(2, 2):
   35005211    -2309581
-1852781081  1101513929
Now the matrix m is:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958           0           0
  719885386  2044897763           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.

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

◆ bottomRows() [1/2]

template<typename Derived>
template<int N>
NRowsBlockXpr< N >::Type Eigen::DenseBase< Derived >::bottomRows ( Index n = N)
inline
Returns
a block consisting of the bottom rows of *this.
Template Parameters
Nthe number of rows in the block as specified at compile-time
Parameters
nthe number of rows in the block as specified at run-time

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:

Here is the array a:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is a.bottomRows<2>():
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Now the array a is:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
          0           0           0           0
          0           0           0           0
See also
block(Index,Index,NRowsType,NColsType), class Block

◆ bottomRows() [2/2]

template<typename Derived>
template<typename NRowsType>
NRowsBlockXpr<... >::Type Eigen::DenseBase< Derived >::bottomRows ( NRowsType n)
inline
Returns
a block consisting of the bottom rows of *this.
Parameters
nthe number of rows in the block
Template Parameters
NRowsTypethe type of the value handling the number of rows in the block, typically Index.

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:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is a.bottomRows(2):
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Now the array a is:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
          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.

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

◆ cast()

template<typename Derived>
template<typename NewType>
CastXpr< NewType >::Type Eigen::DenseBase< Derived >::cast ( ) const
inline
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.

See also
class CwiseUnaryOp

◆ cbegin()

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

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.

See also
cend(), begin()

◆ cend()

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

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.

See also
begin(), cend()

◆ col()

template<typename Derived>
ColXpr Eigen::DenseBase< Derived >::col ( Index i)
inline
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;
Matrix< double, 3, 1 > Vector3d
3×1 vector of type double.
Definition Matrix.h:479
Matrix< double, 3, 3 > Matrix3d
3×3 matrix of type double.
Definition Matrix.h:479

Output:

1 4 0
0 5 0
0 6 1
See also
row(), class Block

◆ colwise() [1/2]

template<typename Derived>
DenseBase< Derived >::ColwiseReturnType Eigen::DenseBase< Derived >::colwise ( )
inline
Returns
a writable VectorwiseOp wrapper of *this providing additional partial reduction operations
See also
rowwise(), class VectorwiseOp, Reductions, visitors and broadcasting

◆ colwise() [2/2]

template<typename Derived>
ConstColwiseReturnType Eigen::DenseBase< Derived >::colwise ( ) const
inline
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.696  0.334  0.445
 0.205  -0.47 -0.633
-0.415  0.928 0.0241
Here is the sum of each column:
 0.487  0.792 -0.163
Here is the maximum absolute value of each column:
0.696 0.928 0.633
See also
rowwise(), class VectorwiseOp, Reductions, visitors and broadcasting

◆ conjugate()

template<typename Derived>
ConjugateReturnType Eigen::DenseBase< Derived >::conjugate ( ) const
inline
Returns
an expression of the complex conjugate of *this.
See also
Math functions, MatrixBase\adjoint()

◆ conjugateIf()

template<typename Derived>
template<bool Cond>
std::conditional_t< Cond, ConjugateReturnType, const Derived & > Eigen::DenseBase< Derived >::conjugateIf ( ) const
inline
Returns
an expression of the complex conjugate of *this if Cond==true, returns derived() otherwise.
See also
conjugate()

◆ Constant() [1/3]

template<typename Derived>
const DenseBase< Derived >::ConstantReturnType Eigen::DenseBase< Derived >::Constant ( const Scalar & value)
inlinestatic
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.

See also
class CwiseNullaryOp

◆ Constant() [2/3]

template<typename Derived>
const DenseBase< Derived >::ConstantReturnType Eigen::DenseBase< Derived >::Constant ( Index rows,
Index cols,
const Scalar & value )
inlinestatic
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.

See also
class CwiseNullaryOp

◆ Constant() [3/3]

template<typename Derived>
const DenseBase< Derived >::ConstantReturnType Eigen::DenseBase< Derived >::Constant ( Index size,
const Scalar & value )
inlinestatic
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.

See also
class CwiseNullaryOp

◆ count()

template<typename Derived>
Index Eigen::DenseBase< Derived >::count ( ) const
Returns
the number of coefficients which evaluate to true
See also
all(), any()

◆ end() [1/2]

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

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.

See also
begin(), cend()

◆ end() [2/2]

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

const version of end()

◆ eval()

template<typename Derived>
EvalReturnType Eigen::DenseBase< Derived >::eval ( ) const
inline
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.

Warning
Be careful with eval() and the auto C++ keyword, as detailed in this page .

◆ fill()

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

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

See also
setConstant(), Constant(), class CwiseNullaryOp

◆ flagged()

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

◆ format()

template<typename Derived>
const WithFormat< Derived > Eigen::DenseBase< Derived >::format ( const IOFormat & fmt) const
inline
Returns
a WithFormat proxy object allowing to print a matrix the with given format fmt.

See class IOFormat for some examples.

See also
class IOFormat, class WithFormat

◆ head() [1/2]

template<typename Derived>
template<int N>
FixedSegmentReturnType< N >::Type Eigen::DenseBase< Derived >::head ( Index n = N)
inline
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.

Template Parameters
Nthe number of coefficients in the segment as specified at compile-time
Parameters
nthe number of coefficients in the segment as specified at run-time

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;
Matrix< int, 1, 4 > RowVector4i
1×4 vector of type int.
Definition Matrix.h:477

Output:

Here is the vector v:
1804289383 -465790871 -189735855  719885386
Here is v.head(2):
1804289383 -465790871
Now the vector v is:
         0          0 -189735855  719885386
See also
head(NType), class Block

◆ head() [2/2]

template<typename Derived>
template<typename NType>
FixedSegmentReturnType<... >::Type Eigen::DenseBase< Derived >::head ( NType n)
inline
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.

Parameters
nthe number of coefficients in the segment
Template Parameters
NTypethe type of the value handling the number of coefficients in the segment, typically Index.

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:
1804289383 -465790871 -189735855  719885386
Here is v.head(2):
1804289383 -465790871
Now the vector v is:
         0          0 -189735855  719885386

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.

Note
Even in the case that the returned expression has dynamic size, in the case when it is applied to a fixed-size vector, it inherits a fixed maximal size, which means that evaluating it does not cause a dynamic memory allocation.
See also
class Block, block(Index,Index)

◆ imag() [1/2]

template<typename Derived>
NonConstImagReturnType Eigen::DenseBase< Derived >::imag ( )
inline
Returns
a non const expression of the imaginary part of *this.
See also
real()

◆ imag() [2/2]

template<typename Derived>
const ImagReturnType Eigen::DenseBase< Derived >::imag ( ) const
inline
Returns
an read-only expression of the imaginary part of *this.
See also
real()

◆ innerSize()

template<typename Derived>
EIGEN_CONSTEXPR Index Eigen::DenseBase< Derived >::innerSize ( ) const
inline
Returns
the inner size.
Note
For a vector, this is just the size. For a matrix (non-vector), this is the minor dimension with respect to the storage order, i.e., the number of rows for a column-major matrix, and the number of columns for a row-major matrix.

◆ innerVector() [1/2]

template<typename Derived>
InnerVectorReturnType Eigen::DenseBase< Derived >::innerVector ( Index outer)
inline
Returns
the outer -th column (resp. row) of the matrix *this if *this is col-major (resp. row-major).

◆ innerVector() [2/2]

template<typename Derived>
const ConstInnerVectorReturnType Eigen::DenseBase< Derived >::innerVector ( Index outer) const
inline
Returns
the outer -th column (resp. row) of the matrix *this if *this is col-major (resp. row-major). Read-only.

◆ innerVectors() [1/2]

template<typename Derived>
InnerVectorsReturnType Eigen::DenseBase< Derived >::innerVectors ( Index outerStart,
Index outerSize )
inline
Returns
the outer -th column (resp. row) of the matrix *this if *this is col-major (resp. row-major).

◆ innerVectors() [2/2]

template<typename Derived>
const ConstInnerVectorsReturnType Eigen::DenseBase< Derived >::innerVectors ( Index outerStart,
Index outerSize ) const
inline
Returns
the outer -th column (resp. row) of the matrix *this if *this is col-major (resp. row-major). Read-only.

◆ isApprox()

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.
Note
The fuzzy compares are done multiplicatively. Two vectors $ v $ and $ w $ are considered to be approximately equal within precision $ p $ if

\[ \Vert v - w \Vert \leqslant p\,\min(\Vert v\Vert, \Vert w\Vert). \]

For matrices, the comparison is done using the Hilbert-Schmidt norm (aka Frobenius norm L2 norm).
Because of the multiplicativeness of this comparison, one can't use this function to check whether *this is approximately equal to the zero matrix or vector. Indeed, isApprox(zero) returns false unless *this itself is exactly the zero matrix or vector. If you want to test whether *this is zero, use internal::isMuchSmallerThan(const RealScalar&, RealScalar) instead.
See also
internal::isMuchSmallerThan(const RealScalar&, RealScalar) const

◆ isApproxToConstant()

template<typename Derived>
bool Eigen::DenseBase< Derived >::isApproxToConstant ( const Scalar & val,
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

◆ isConstant()

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

This is just an alias for isApproxToConstant().

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

◆ isMuchSmallerThan() [1/2]

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.
Note
The fuzzy compares are done multiplicatively. A vector $ v $ is considered to be much smaller than a vector $ w $ within precision $ p $ if

\[ \Vert v \Vert \leqslant p\,\Vert w\Vert. \]

For matrices, the comparison is done using the Hilbert-Schmidt norm.
See also
isApprox(), isMuchSmallerThan(const RealScalar&, RealScalar) const

◆ isMuchSmallerThan() [2/2]

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.
Note
The fuzzy compares are done multiplicatively. A vector $ v $ is considered to be much smaller than $ x $ within precision $ p $ if

\[ \Vert v \Vert \leqslant p\,\vert x\vert. \]

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.

See also
isApprox(), isMuchSmallerThan(const DenseBase<OtherDerived>&, RealScalar) const

◆ isOnes()

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
See also
class CwiseNullaryOp, Ones()

◆ isZero()

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
See also
class CwiseNullaryOp, Zero()

◆ lazyAssign()

template<typename Derived>
template<typename OtherDerived>
EIGEN_DEPRECATED Derived & Eigen::DenseBase< Derived >::lazyAssign ( const DenseBase< OtherDerived > & other)

◆ leftCols() [1/2]

template<typename Derived>
template<int N>
NColsBlockXpr< N >::Type Eigen::DenseBase< Derived >::leftCols ( Index n = N)
inline
Returns
a block consisting of the left columns of *this.
Template Parameters
Nthe number of columns in the block as specified at compile-time
Parameters
nthe number of columns in the block as specified at run-time

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:

Here is the array a:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is a.leftCols<2>():
 1804289383 -1550966999
 -465790871 -1122281286
 -189735855 -1364114958
  719885386  2044897763
Now the array a is:
          0           0  1365180540   336465782
          0           0   304089172 -1868760786
          0           0    35005211    -2309581
          0           0 -1852781081  1101513929
See also
block(Index,Index,NRowsType,NColsType), class Block

◆ leftCols() [2/2]

template<typename Derived>
template<typename NColsType>
NColsBlockXpr<... >::Type Eigen::DenseBase< Derived >::leftCols ( NColsType n)
inline
Returns
a block consisting of the left columns of *this.
Parameters
nthe number of columns in the block
Template Parameters
NColsTypethe type of the value handling the number of columns in the block, typically Index.

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:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is a.leftCols(2):
 1804289383 -1550966999
 -465790871 -1122281286
 -189735855 -1364114958
  719885386  2044897763
Now the array a is:
          0           0  1365180540   336465782
          0           0   304089172 -1868760786
          0           0    35005211    -2309581
          0           0 -1852781081  1101513929

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.

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

◆ LinSpaced() [1/4]

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

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
See also
setLinSpaced(Index,const Scalar&,const Scalar&), CwiseNullaryOp Special version for fixed size types which does not require the size parameter.

◆ LinSpaced() [2/4]

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

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
See also
setLinSpaced(Index,const Scalar&,const Scalar&), CwiseNullaryOp

◆ LinSpaced() [3/4]

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

◆ LinSpaced() [4/4]

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

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(Sequential, 4, 7, 10).transpose() << endl;
cout << VectorXd::LinSpaced(Sequential, 5, 0.0, 1.0).transpose() << endl;

Output:

 7  8  9 10
   0 0.25  0.5 0.75    1
See also
LinSpaced(Index,const Scalar&, const Scalar&), setLinSpaced(Index,const Scalar&,const Scalar&)

◆ maxCoeff() [1/3]

template<typename Derived>
template<int NaNPropagation>
internal::traits< Derived >::Scalar Eigen::DenseBase< Derived >::maxCoeff ( ) const
inline
Returns
the maximum of all coefficients of *this. In case *this contains NaN, NaNPropagation determines the behavior: NaNPropagation == PropagateFast : undefined NaNPropagation == PropagateNaN : result is NaN NaNPropagation == PropagateNumbers : result is maximum of elements that are not NaN
Warning
the matrix must be not empty, otherwise an assertion is triggered.

◆ maxCoeff() [2/3]

template<typename Derived>
template<int NaNPropagation, 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.

In case *this contains NaN, NaNPropagation determines the behavior: NaNPropagation == PropagateFast : undefined NaNPropagation == PropagateNaN : result is NaN NaNPropagation == PropagateNumbers : result is maximum of elements that are not NaN

Warning
the matrix must be not empty, otherwise an assertion is triggered.
See also
DenseBase::maxCoeff(IndexType*,IndexType*), DenseBase::minCoeff(IndexType*,IndexType*), DenseBase::visitor(), DenseBase::maxCoeff()

◆ maxCoeff() [3/3]

template<typename Derived>
template<int NaNPropagation, typename IndexType>
internal::traits< Derived >::Scalar Eigen::DenseBase< Derived >::maxCoeff ( IndexType * rowId,
IndexType * colId ) const
Returns
the maximum of all coefficients of *this and puts in *row and *col its location.

In case *this contains NaN, NaNPropagation determines the behavior: NaNPropagation == PropagateFast : undefined NaNPropagation == PropagateNaN : result is NaN NaNPropagation == PropagateNumbers : result is maximum of elements that are not NaN

Warning
the matrix must be not empty, otherwise an assertion is triggered.
See also
DenseBase::minCoeff(IndexType*,IndexType*), DenseBase::visit(), DenseBase::maxCoeff()

◆ mean()

template<typename Derived>
internal::traits< Derived >::Scalar Eigen::DenseBase< Derived >::mean ( ) const
inline
Returns
the mean of all coefficients of *this
See also
trace(), prod(), sum()

◆ middleCols() [1/2]

template<typename Derived>
template<int N>
NColsBlockXpr< N >::Type Eigen::DenseBase< Derived >::middleCols ( Index startCol,
Index n = N )
inline
Returns
a block consisting of a range of columns of *this.
Template Parameters
Nthe number of columns in the block as specified at compile-time
Parameters
startColthe index of the first column in the block
nthe number of columns in the block as specified at run-time

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>
int main() {
int const N = 5;
Eigen::MatrixXi A(N, N);
A.setRandom();
std::cout << "A =\n" << A << '\n' << std::endl;
std::cout << "A(:,1..3) =\n" << A.middleCols<3>(1) << std::endl;
return 0;
}
Matrix< int, Dynamic, Dynamic > MatrixXi
Dynamic×Dynamic matrix of type int.
Definition Matrix.h:477

Output:

A =
 1804289383 -1122281286    35005211  1101513929 -1016307419
 -465790871 -1364114958 -1852781081  1315634022 -1287999227
 -189735855  2044897763   336465782  -778350579 -1539069864
  719885386  1365180540 -1868760786  1059961393  1734575198
-1550966999   304089172    -2309581   628175011   149798315

A(:,1..3) =
-1122281286    35005211  1101513929
-1364114958 -1852781081  1315634022
 2044897763   336465782  -778350579
 1365180540 -1868760786  1059961393
  304089172    -2309581   628175011
See also
block(Index,Index,NRowsType,NColsType), class Block

◆ middleCols() [2/2]

template<typename Derived>
template<typename NColsType>
NColsBlockXpr<... >::Type Eigen::DenseBase< Derived >::middleCols ( Index startCol,
NColsType numCols )
inline
Returns
a block consisting of a range of columns of *this.
Parameters
startColthe index of the first column in the block
numColsthe number of columns in the block
Template Parameters
NColsTypethe type of the value handling the number of columns in the block, typically Index.

Example:

#include <Eigen/Core>
#include <iostream>
int main() {
int const N = 5;
Eigen::MatrixXi A(N, N);
A.setRandom();
std::cout << "A =\n" << A << '\n' << std::endl;
std::cout << "A(1..3,:) =\n" << A.middleCols(1, 3) << std::endl;
return 0;
}

Output:

A =
 1804289383 -1122281286    35005211  1101513929 -1016307419
 -465790871 -1364114958 -1852781081  1315634022 -1287999227
 -189735855  2044897763   336465782  -778350579 -1539069864
  719885386  1365180540 -1868760786  1059961393  1734575198
-1550966999   304089172    -2309581   628175011   149798315

A(1..3,:) =
-1122281286    35005211  1101513929
-1364114958 -1852781081  1315634022
 2044897763   336465782  -778350579
 1365180540 -1868760786  1059961393
  304089172    -2309581   628175011

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.

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

◆ middleRows() [1/2]

template<typename Derived>
template<int N>
NRowsBlockXpr< N >::Type Eigen::DenseBase< Derived >::middleRows ( Index startRow,
Index n = N )
inline
Returns
a block consisting of a range of rows of *this.
Template Parameters
Nthe number of rows in the block as specified at compile-time
Parameters
startRowthe index of the first row in the block
nthe number of rows in the block as specified at run-time

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>
int main() {
int const N = 5;
Eigen::MatrixXi A(N, N);
A.setRandom();
std::cout << "A =\n" << A << '\n' << std::endl;
std::cout << "A(1..3,:) =\n" << A.middleRows<3>(1) << std::endl;
return 0;
}

Output:

A =
 1804289383 -1122281286    35005211  1101513929 -1016307419
 -465790871 -1364114958 -1852781081  1315634022 -1287999227
 -189735855  2044897763   336465782  -778350579 -1539069864
  719885386  1365180540 -1868760786  1059961393  1734575198
-1550966999   304089172    -2309581   628175011   149798315

A(1..3,:) =
 -465790871 -1364114958 -1852781081  1315634022 -1287999227
 -189735855  2044897763   336465782  -778350579 -1539069864
  719885386  1365180540 -1868760786  1059961393  1734575198
See also
block(Index,Index,NRowsType,NColsType), class Block

◆ middleRows() [2/2]

template<typename Derived>
template<typename NRowsType>
NRowsBlockXpr<... >::Type Eigen::DenseBase< Derived >::middleRows ( Index startRow,
NRowsType n )
inline
Returns
a block consisting of a range of rows of *this.
Parameters
startRowthe index of the first row in the block
nthe number of rows in the block
Template Parameters
NRowsTypethe type of the value handling the number of rows in the block, typically Index.

Example:

#include <Eigen/Core>
#include <iostream>
int main() {
int const N = 5;
Eigen::MatrixXi A(N, N);
A.setRandom();
std::cout << "A =\n" << A << '\n' << std::endl;
std::cout << "A(2..3,:) =\n" << A.middleRows(2, 2) << std::endl;
return 0;
}

Output:

A =
 1804289383 -1122281286    35005211  1101513929 -1016307419
 -465790871 -1364114958 -1852781081  1315634022 -1287999227
 -189735855  2044897763   336465782  -778350579 -1539069864
  719885386  1365180540 -1868760786  1059961393  1734575198
-1550966999   304089172    -2309581   628175011   149798315

A(2..3,:) =
 -189735855  2044897763   336465782  -778350579 -1539069864
  719885386  1365180540 -1868760786  1059961393  1734575198

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.

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

◆ minCoeff() [1/3]

template<typename Derived>
template<int NaNPropagation>
internal::traits< Derived >::Scalar Eigen::DenseBase< Derived >::minCoeff ( ) const
inline
Returns
the minimum of all coefficients of *this. In case *this contains NaN, NaNPropagation determines the behavior: NaNPropagation == PropagateFast : undefined NaNPropagation == PropagateNaN : result is NaN NaNPropagation == PropagateNumbers : result is minimum of elements that are not NaN
Warning
the matrix must be not empty, otherwise an assertion is triggered.

◆ minCoeff() [2/3]

template<typename Derived>
template<int NaNPropagation, 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.

In case *this contains NaN, NaNPropagation determines the behavior: NaNPropagation == PropagateFast : undefined NaNPropagation == PropagateNaN : result is NaN NaNPropagation == PropagateNumbers : result is maximum of elements that are not NaN

Warning
the matrix must be not empty, otherwise an assertion is triggered.
See also
DenseBase::minCoeff(IndexType*,IndexType*), DenseBase::maxCoeff(IndexType*,IndexType*), DenseBase::visit(), DenseBase::minCoeff()

◆ minCoeff() [3/3]

template<typename Derived>
template<int NaNPropagation, typename IndexType>
internal::traits< Derived >::Scalar Eigen::DenseBase< Derived >::minCoeff ( IndexType * rowId,
IndexType * colId ) const
Returns
the minimum of all coefficients of *this and puts in *row and *col its location.

In case *this contains NaN, NaNPropagation determines the behavior: NaNPropagation == PropagateFast : undefined NaNPropagation == PropagateNaN : result is NaN NaNPropagation == PropagateNumbers : result is maximum of elements that are not NaN

Warning
the matrix must be not empty, otherwise an assertion is triggered.
See also
DenseBase::minCoeff(Index*), DenseBase::maxCoeff(Index*,Index*), DenseBase::visit(), DenseBase::minCoeff()

◆ nestByValue()

template<typename Derived>
const NestByValue< Derived > Eigen::DenseBase< Derived >::nestByValue ( ) const
inline
Returns
an expression of the temporary version of *this.

◆ NullaryExpr() [1/3]

template<typename Derived>
template<typename CustomNullaryOp>
const CwiseNullaryOp< CustomNullaryOp, PlainObject > Eigen::DenseBase< Derived >::NullaryExpr ( const CustomNullaryOp & func)
inlinestatic
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.

See also
class CwiseNullaryOp

◆ NullaryExpr() [2/3]

template<typename Derived>
template<typename CustomNullaryOp>
const CwiseNullaryOp< CustomNullaryOp, PlainObject > Eigen::DenseBase< Derived >::NullaryExpr ( Index rows,
Index cols,
const CustomNullaryOp & func )
inlinestatic
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.

See also
class CwiseNullaryOp

◆ NullaryExpr() [3/3]

template<typename Derived>
template<typename CustomNullaryOp>
const CwiseNullaryOp< CustomNullaryOp, PlainObject > Eigen::DenseBase< Derived >::NullaryExpr ( Index size,
const CustomNullaryOp & func )
inlinestatic
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>
int main() {
std::default_random_engine generator;
std::poisson_distribution<int> distribution(4.1);
auto poisson = [&]() { return distribution(generator); };
Eigen::RowVectorXi v = Eigen::RowVectorXi::NullaryExpr(10, poisson);
std::cout << v << "\n";
}
Matrix< int, 1, Dynamic > RowVectorXi
1×Dynamic vector of type int.
Definition Matrix.h:477

Output:

2 3 1 4 3 4 4 3 2 3
See also
class CwiseNullaryOp

◆ Ones() [1/3]

template<typename Derived>
const DenseBase< Derived >::ConstantReturnType Eigen::DenseBase< Derived >::Ones ( )
inlinestatic
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
See also
Ones(Index), Ones(Index,Index), isOnes(), class Ones

◆ Ones() [2/3]

template<typename Derived>
const DenseBase< Derived >::ConstantReturnType Eigen::DenseBase< Derived >::Ones ( Index rows,
Index cols )
inlinestatic
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
See also
Ones(), Ones(Index), isOnes(), class Ones

◆ Ones() [3/3]

template<typename Derived>
const DenseBase< Derived >::ConstantReturnType Eigen::DenseBase< Derived >::Ones ( Index newSize)
inlinestatic
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
See also
Ones(), Ones(Index,Index), isOnes(), class Ones

◆ operator()() [1/2]

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.

◆ operator()() [2/2]

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\placeholders\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:

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.

See also
operator()(const Indices&), class Block, class IndexedView, DenseBase\block(Index,Index,Index,Index)

◆ operator-()

template<typename Derived>
const NegativeReturnType Eigen::DenseBase< Derived >::operator- ( ) const
inline
Returns
an expression of the opposite of *this

◆ operator<<() [1/2]

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

◆ operator<<() [2/2]

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

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:

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;
v1 << 14, 15;
m2 << v1.transpose(), 16, v1, m1.block(1, 1, 2, 2);
cout << m2 << endl;
Matrix< int, 3, 3 > Matrix3i
3×3 matrix of type int.
Definition Matrix.h:477
Matrix< int, 2, 1 > Vector2i
2×1 vector of type int.
Definition Matrix.h:477

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
Note
According the c++ standard, the argument expressions of this comma initializer are evaluated in arbitrary order.
See also
CommaInitializer::finished(), class CommaInitializer

◆ operator=() [1/3]

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

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

◆ operator=() [2/3]

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

Copies other into *this.

Returns
a reference to *this.

◆ operator=() [3/3]

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

Copies the generic expression other into *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

Returns
a reference to *this.

◆ outerSize()

template<typename Derived>
EIGEN_CONSTEXPR Index Eigen::DenseBase< Derived >::outerSize ( ) const
inline
Returns
the outer size.
Note
For a vector, this returns just 1. For a matrix (non-vector), this is the major dimension with respect to the storage order, i.e., the number of columns for a column-major matrix, and the number of rows for a row-major matrix.

◆ prod()

template<typename Derived>
internal::traits< Derived >::Scalar Eigen::DenseBase< Derived >::prod ( ) const
inline
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.696  0.334  0.445
 0.205  -0.47 -0.633
-0.415  0.928 0.0241
Here is the product of all the coefficients:
-5.85e-05
See also
sum(), mean(), trace()

◆ Random() [1/3]

template<typename Derived>
const DenseBase< Derived >::RandomReturnType Eigen::DenseBase< Derived >::Random ( )
inlinestatic
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:

   40311868 -1793716316
  665553156 -1025905432

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.

Warning
This function is not re-entrant.
See also
DenseBase::setRandom(), DenseBase::Random(Index,Index), DenseBase::Random(Index)

◆ Random() [2/3]

template<typename Derived>
const DenseBase< Derived >::RandomReturnType Eigen::DenseBase< Derived >::Random ( Index rows,
Index cols )
inlinestatic
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.

Warning
This function is not re-entrant.

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:

 1804289383  -189735855 -1550966999
 -465790871   719885386 -1122281286

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.

See also
DenseBase::setRandom(), DenseBase::Random(Index), DenseBase::Random()

◆ Random() [3/3]

template<typename Derived>
const DenseBase< Derived >::RandomReturnType Eigen::DenseBase< Derived >::Random ( Index size)
inlinestatic
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.

Warning
This function is not re-entrant.

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:

1804289383
-465790871

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.

See also
DenseBase::setRandom(), DenseBase::Random(Index,Index), DenseBase::Random()

◆ real() [1/2]

template<typename Derived>
NonConstRealReturnType Eigen::DenseBase< Derived >::real ( )
inline
Returns
a non const expression of the real part of *this.
See also
imag()

◆ real() [2/2]

template<typename Derived>
RealReturnType Eigen::DenseBase< Derived >::real ( ) const
inline
Returns
a read-only expression of the real part of *this.
See also
imag()

◆ redux()

template<typename Derived>
template<typename Func>
internal::traits< Derived >::Scalar Eigen::DenseBase< Derived >::redux ( const Func & func) const
inline
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.

Warning
the matrix must be not empty, otherwise an assertion is triggered.
See also
DenseBase::sum(), DenseBase::minCoeff(), DenseBase::maxCoeff(), MatrixBase::colwise(), MatrixBase::rowwise()

◆ replicate() [1/2]

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:
 1804289383  -189735855 -1550966999
 -465790871   719885386 -1122281286
m.replicate<3,2>() = ...
 1804289383  -189735855 -1550966999  1804289383  -189735855 -1550966999
 -465790871   719885386 -1122281286  -465790871   719885386 -1122281286
 1804289383  -189735855 -1550966999  1804289383  -189735855 -1550966999
 -465790871   719885386 -1122281286  -465790871   719885386 -1122281286
 1804289383  -189735855 -1550966999  1804289383  -189735855 -1550966999
 -465790871   719885386 -1122281286  -465790871   719885386 -1122281286
See also
VectorwiseOp::replicate(), DenseBase::replicate(Index,Index), class Replicate

◆ replicate() [2/2]

template<typename Derived>
const Replicate< Derived, Dynamic, Dynamic > Eigen::DenseBase< Derived >::replicate ( Index rowFactor,
Index colFactor ) const
inline
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;
Matrix< int, 3, 1 > Vector3i
3×1 vector of type int.
Definition Matrix.h:477

Output:

Here is the vector v:
1804289383
-465790871
-189735855
v.replicate(2,5) = ...
1804289383 1804289383 1804289383 1804289383 1804289383
-465790871 -465790871 -465790871 -465790871 -465790871
-189735855 -189735855 -189735855 -189735855 -189735855
1804289383 1804289383 1804289383 1804289383 1804289383
-465790871 -465790871 -465790871 -465790871 -465790871
-189735855 -189735855 -189735855 -189735855 -189735855
See also
VectorwiseOp::replicate(), DenseBase::replicate<int,int>(), class Replicate

◆ reshaped() [1/2]

template<typename Derived>
template<int Order = ColMajor>
Reshaped< Derived,... > Eigen::DenseBase< Derived >::reshaped ( )
inline
Returns
an expression of *this with columns (or rows) stacked to a linear column vector
Template Parameters
Orderspecifies 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.

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;
TransposeReturnType transpose()
Definition Transpose.h:162
@ RowMajor
Definition Constants.h:320

Output:

Here is the matrix m:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.reshaped().transpose():
 1804289383  -465790871  -189735855   719885386 -1550966999 -1122281286 -1364114958  2044897763  1365180540   304089172    35005211 -1852781081   336465782 -1868760786    -2309581  1101513929
Here is m.reshaped<RowMajor>().transpose():  
 1804289383 -1550966999  1365180540   336465782  -465790871 -1122281286   304089172 -1868760786  -189735855 -1364114958    35005211    -2309581   719885386  2044897763 -1852781081  1101513929

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

See also
reshaped(NRowsType,NColsType), class Reshaped

◆ reshaped() [2/2]

template<typename Derived>
template<int Order = ColMajor, typename NRowsType, typename NColsType>
Reshaped< Derived,... > Eigen::DenseBase< Derived >::reshaped ( NRowsType nRows,
NColsType nCols )
inline
Returns
an expression of *this with reshaped sizes.
Parameters
nRowsthe number of rows in the reshaped expression, specified at either run-time or compile-time, or AutoSize
nColsthe number of columns in the reshaped expression, specified at either run-time or compile-time, or AutoSize
Template Parameters
Orderspecifies 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.
NRowsTypethe type of the value handling the number of rows, typically Index.
NColsTypethe type of the value handling the number of columns, typically Index.

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:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.reshaped(2, 8):
 1804289383  -189735855 -1550966999 -1364114958  1365180540    35005211   336465782    -2309581
 -465790871   719885386 -1122281286  2044897763   304089172 -1852781081 -1868760786  1101513929

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:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.reshaped(fix<2>,fix<8>):
 1804289383  -189735855 -1550966999 -1364114958  1365180540    35005211   336465782    -2309581
 -465790871   719885386 -1122281286  2044897763   304089172 -1852781081 -1868760786  1101513929

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:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.reshaped(2, AutoSize):
 1804289383  -189735855 -1550966999 -1364114958  1365180540    35005211   336465782    -2309581
 -465790871   719885386 -1122281286  2044897763   304089172 -1852781081 -1868760786  1101513929
Here is m.reshaped<RowMajor>(AutoSize, fix<8>):
 1804289383 -1550966999  1365180540   336465782  -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581   719885386  2044897763 -1852781081  1101513929

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.

See also
class Reshaped, fix, fix<N>(int)

◆ resize() [1/2]

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

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.

◆ resize() [2/2]

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

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.

◆ reverse() [1/2]

template<typename Derived>
DenseBase< Derived >::ReverseReturnType Eigen::DenseBase< Derived >::reverse ( )
inline
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:
 1804289383   719885386 -1364114958   304089172
 -465790871 -1550966999  2044897763    35005211
 -189735855 -1122281286  1365180540 -1852781081
Here is the reverse of m:
-1852781081  1365180540 -1122281286  -189735855
   35005211  2044897763 -1550966999  -465790871
  304089172 -1364114958   719885386  1804289383
Here is the coefficient (1,0) in the reverse of m:
35005211
Let us overwrite this coefficient with the value 4.
Now the matrix m is:
 1804289383   719885386 -1364114958   304089172
 -465790871 -1550966999  2044897763           4
 -189735855 -1122281286  1365180540 -1852781081

◆ reverse() [2/2]

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

This is the const version of reverse().

◆ reverseInPlace()

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

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.)
See also
VectorwiseOp::reverseInPlace(), reverse()

◆ rightCols() [1/2]

template<typename Derived>
template<int N>
NColsBlockXpr< N >::Type Eigen::DenseBase< Derived >::rightCols ( Index n = N)
inline
Returns
a block consisting of the right columns of *this.
Template Parameters
Nthe number of columns in the block as specified at compile-time
Parameters
nthe number of columns in the block as specified at run-time

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:

Here is the array a:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is a.rightCols<2>():
 1365180540   336465782
  304089172 -1868760786
   35005211    -2309581
-1852781081  1101513929
Now the array a is:
 1804289383 -1550966999           0           0
 -465790871 -1122281286           0           0
 -189735855 -1364114958           0           0
  719885386  2044897763           0           0
See also
block(Index,Index,NRowsType,NColsType), class Block

◆ rightCols() [2/2]

template<typename Derived>
template<typename NColsType>
NColsBlockXpr<... >::Type Eigen::DenseBase< Derived >::rightCols ( NColsType n)
inline
Returns
a block consisting of the right columns of *this.
Parameters
nthe number of columns in the block
Template Parameters
NColsTypethe type of the value handling the number of columns in the block, typically Index.

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:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is a.rightCols(2):
 1365180540   336465782
  304089172 -1868760786
   35005211    -2309581
-1852781081  1101513929
Now the array a is:
 1804289383 -1550966999           0           0
 -465790871 -1122281286           0           0
 -189735855 -1364114958           0           0
  719885386  2044897763           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.

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

◆ row()

template<typename Derived>
RowXpr Eigen::DenseBase< Derived >::row ( Index i)
inline
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:

1 0 0
4 5 6
0 0 1
See also
col(), class Block

◆ rowwise() [1/2]

template<typename Derived>
DenseBase< Derived >::RowwiseReturnType Eigen::DenseBase< Derived >::rowwise ( )
inline
Returns
a writable VectorwiseOp wrapper of *this providing additional partial reduction operations
See also
colwise(), class VectorwiseOp, Reductions, visitors and broadcasting

◆ rowwise() [2/2]

template<typename Derived>
ConstRowwiseReturnType Eigen::DenseBase< Derived >::rowwise ( ) const
inline
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.696  0.334  0.445
 0.205  -0.47 -0.633
-0.415  0.928 0.0241
Here is the sum of each row:
  1.48
-0.897
 0.537
Here is the maximum absolute value of each row:
0.696
0.633
0.928
See also
colwise(), class VectorwiseOp, Reductions, visitors and broadcasting

◆ segment() [1/2]

template<typename Derived>
template<int N>
FixedSegmentReturnType< N >::Type Eigen::DenseBase< Derived >::segment ( Index start,
Index n = N )
inline
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.

Template Parameters
Nthe number of coefficients in the segment as specified at compile-time
Parameters
startthe index of the first element in the segment
nthe number of coefficients in the segment as specified at compile-time

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:
1804289383 -465790871 -189735855  719885386
Here is v.segment<2>(1):
-465790871 -189735855
Now the vector v is:
1804289383 -465790871          0          0
See also
segment(Index,NType), class Block

◆ segment() [2/2]

template<typename Derived>
template<typename NType>
FixedSegmentReturnType<... >::Type Eigen::DenseBase< Derived >::segment ( Index start,
NType n )
inline
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.

Parameters
startthe first coefficient in the segment
nthe number of coefficients in the segment
Template Parameters
NTypethe type of the value handling the number of coefficients in the segment, typically Index.

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:
1804289383 -465790871 -189735855  719885386
Here is v.segment(1, 2):
-465790871 -189735855
Now the vector v is:
1804289383          0          0  719885386

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.

Note
Even in the case that the returned expression has dynamic size, in the case when it is applied to a fixed-size vector, it inherits a fixed maximal size, which means that evaluating it does not cause a dynamic memory allocation.
See also
block(Index,Index,NRowsType,NColsType), fix<N>, fix<N>(int), class Block

◆ select() [1/3]

template<typename Derived>
template<typename ThenDerived, typename ElseDerived>
CwiseTernaryOp< internal::scalar_boolean_select_op< typename DenseBase< ThenDerived >::Scalar, typename DenseBase< ElseDerived >::Scalar, typename DenseBase< Derived >::Scalar >, ThenDerived, ElseDerived, Derived > Eigen::DenseBase< Derived >::select ( const DenseBase< ThenDerived > & thenMatrix,
const DenseBase< ElseDerived > & elseMatrix ) const
inline
Returns
a matrix where each coefficient (i,j) is equal to thenMatrix(i,j) if *this(i,j) != Scalar(0), 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
See also
DenseBase::bitwiseSelect(const DenseBase<ThenDerived>&, const DenseBase<ElseDerived>&)

◆ select() [2/3]

template<typename Derived>
template<typename ThenDerived>
CwiseTernaryOp< internal::scalar_boolean_select_op< typename DenseBase< ThenDerived >::Scalar, typename DenseBase< ThenDerived >::Scalar, typename DenseBase< Derived >::Scalar >, ThenDerived, typename DenseBase< ThenDerived >::ConstantReturnType, Derived > Eigen::DenseBase< Derived >::select ( const DenseBase< ThenDerived > & thenMatrix,
const typename DenseBase< ThenDerived >::Scalar & elseScalar ) const
inline

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

See also
DenseBase::booleanSelect(const DenseBase<ThenDerived>&, const DenseBase<ElseDerived>&) const, class Select

◆ select() [3/3]

template<typename Derived>
template<typename ElseDerived>
CwiseTernaryOp< internal::scalar_boolean_select_op< typename DenseBase< ElseDerived >::Scalar, typename DenseBase< ElseDerived >::Scalar, typename DenseBase< Derived >::Scalar >, typename DenseBase< ElseDerived >::ConstantReturnType, ElseDerived, Derived > Eigen::DenseBase< Derived >::select ( const typename DenseBase< ElseDerived >::Scalar & thenScalar,
const DenseBase< ElseDerived > & elseMatrix ) const
inline

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

See also
DenseBase::booleanSelect(const DenseBase<ThenDerived>&, const DenseBase<ElseDerived>&) const, class Select

◆ setConstant()

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

Sets all coefficients in this expression to value val.

See also
fill(), setConstant(Index,const Scalar&), setConstant(Index,Index,const Scalar&), setZero(), setOnes(), Constant(), class CwiseNullaryOp, setZero(), setOnes()

◆ setLinSpaced() [1/2]

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

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 .

See also
LinSpaced(Index,const Scalar&,const Scalar&), setLinSpaced(Index, const Scalar&, const Scalar&), CwiseNullaryOp

◆ setLinSpaced() [2/2]

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

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:

v.setLinSpaced(5, 0.5f, 1.5f);
cout << v << endl;
Matrix< float, Dynamic, 1 > VectorXf
Dynamic×1 vector of type float.
Definition Matrix.h:478

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 .

See also
LinSpaced(Index,const Scalar&,const Scalar&), CwiseNullaryOp

◆ setOnes()

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

Sets all coefficients in this expression to one.

Example:

Matrix4i m = Matrix4i::Random();
m.row(1).setOnes();
cout << m << endl;
Derived & setOnes(Index size)
Definition CwiseNullaryOp.h:708

Output:

 1804289383 -1550966999  1365180540   336465782
          1           1           1           1
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
See also
class CwiseNullaryOp, Ones()

◆ setRandom()

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

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.

Warning
This function is not re-entrant.

Example:

Matrix4i m = Matrix4i::Zero();
m.col(1).setRandom();
cout << m << endl;
Derived & setRandom(Index size)
Definition Random.h:147

Output:

         0 1804289383          0          0
         0 -465790871          0          0
         0 -189735855          0          0
         0  719885386          0          0
See also
class CwiseNullaryOp, setRandom(Index), setRandom(Index,Index)

◆ setZero()

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

Sets all coefficients in this expression to zero.

Example:

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

Output:

 1804289383 -1550966999  1365180540   336465782
          0           0           0           0
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
See also
class CwiseNullaryOp, Zero()

◆ subVector() [1/2]

template<typename Derived>
template<DirectionType Direction>
std::conditional_t< Direction==Vertical, ColXpr, RowXpr > Eigen::DenseBase< Derived >::subVector ( Index i)
inline
Returns
the i-th subvector (column or vector) according to the Direction
See also
subVectors()

◆ subVector() [2/2]

template<typename Derived>
template<DirectionType Direction>
std::conditional_t< Direction==Vertical, ConstColXpr, ConstRowXpr > Eigen::DenseBase< Derived >::subVector ( Index i) const
inline

This is the const version of subVector(Index)

◆ subVectors()

template<typename Derived>
template<DirectionType Direction>
EIGEN_CONSTEXPR Index Eigen::DenseBase< Derived >::subVectors ( ) const
inline
Returns
the number of subvectors (rows or columns) in the direction Direction
See also
subVector(Index)

◆ sum()

template<typename Derived>
internal::traits< Derived >::Scalar Eigen::DenseBase< Derived >::sum ( ) const
inline
Returns
the sum of all coefficients of *this

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

See also
trace(), prod(), mean()

◆ swap() [1/2]

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

swaps *this with the expression other.

◆ swap() [2/2]

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

swaps *this with the matrix or array other.

◆ tail() [1/2]

template<typename Derived>
template<int N>
FixedSegmentReturnType< N >::Type Eigen::DenseBase< Derived >::tail ( Index n = N)
inline
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.

Template Parameters
Nthe number of coefficients in the segment as specified at compile-time
Parameters
nthe number of coefficients in the segment as specified at run-time

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:
1804289383 -465790871 -189735855  719885386
Here is v.tail(2):
-189735855  719885386
Now the vector v is:
1804289383 -465790871          0          0
See also
tail(NType), class Block

◆ tail() [2/2]

template<typename Derived>
template<typename NType>
FixedSegmentReturnType<... >::Type Eigen::DenseBase< Derived >::tail ( NType n)
inline
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.

Parameters
nthe number of coefficients in the segment
Template Parameters
NTypethe type of the value handling the number of coefficients in the segment, typically Index.

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:
1804289383 -465790871 -189735855  719885386
Here is v.tail(2):
-189735855  719885386
Now the vector v is:
1804289383 -465790871          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.

Note
Even in the case that the returned expression has dynamic size, in the case when it is applied to a fixed-size vector, it inherits a fixed maximal size, which means that evaluating it does not cause a dynamic memory allocation.
See also
class Block, block(Index,Index)

◆ topLeftCorner() [1/3]

template<typename Derived>
template<int CRows, int CCols>
FixedBlockXpr< CRows, CCols >::Type Eigen::DenseBase< Derived >::topLeftCorner ( )
inline
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:

Here is the matrix m:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.topLeftCorner<2,2>():
 1804289383 -1550966999
 -465790871 -1122281286
Now the matrix m is:
          0           0  1365180540   336465782
          0           0   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
See also
block(Index,Index,NRowsType,NColsType), class Block

◆ topLeftCorner() [2/3]

template<typename Derived>
template<int CRows, int CCols>
FixedBlockXpr< CRows, CCols >::Type Eigen::DenseBase< Derived >::topLeftCorner ( Index cRows,
Index cCols )
inline
Returns
an expression of a top-left corner of *this.
Template Parameters
CRowsnumber of rows in corner as specified at compile-time
CColsnumber of columns in corner as specified at compile-time
Parameters
cRowsnumber of rows in corner as specified at run-time
cColsnumber of columns in corner as specified at run-time

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:

Here is the matrix m:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.topLeftCorner<2,Dynamic>(2,2):
 1804289383 -1550966999
 -465790871 -1122281286
Now the matrix m is:
          0           0  1365180540   336465782
          0           0   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
See also
class Block

◆ topLeftCorner() [3/3]

template<typename Derived>
template<typename NRowsType, typename NColsType>
FixedBlockXpr<...,... >::Type Eigen::DenseBase< Derived >::topLeftCorner ( NRowsType cRows,
NColsType cCols )
inline
Returns
an expression of a top-left corner of *this with either dynamic or fixed sizes.
Parameters
cRowsthe number of rows in the corner
cColsthe number of columns in the corner
Template Parameters
NRowsTypethe type of the value handling the number of rows in the block, typically Index.
NColsTypethe type of the value handling the number of columns in the block, typically Index.

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:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.topLeftCorner(2, 2):
 1804289383 -1550966999
 -465790871 -1122281286
Now the matrix m is:
          0           0  1365180540   336465782
          0           0   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929

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.

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

◆ topRightCorner() [1/3]

template<typename Derived>
template<int CRows, int CCols>
FixedBlockXpr< CRows, CCols >::Type Eigen::DenseBase< Derived >::topRightCorner ( )
inline
Returns
an expression of a fixed-size top-right corner of *this.
Template Parameters
CRowsthe number of rows in the corner
CColsthe number of columns in the corner

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:

Here is the matrix m:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.topRightCorner<2,2>():
 1365180540   336465782
  304089172 -1868760786
Now the matrix m is:
 1804289383 -1550966999           0           0
 -465790871 -1122281286           0           0
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
See also
class Block, block<int,int>(Index,Index)

◆ topRightCorner() [2/3]

template<typename Derived>
template<int CRows, int CCols>
FixedBlockXpr< CRows, CCols >::Type Eigen::DenseBase< Derived >::topRightCorner ( Index cRows,
Index cCols )
inline
Returns
an expression of a top-right corner of *this.
Template Parameters
CRowsnumber of rows in corner as specified at compile-time
CColsnumber of columns in corner as specified at compile-time
Parameters
cRowsnumber of rows in corner as specified at run-time
cColsnumber of columns in corner as specified at run-time

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:

Here is the matrix m:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.topRightCorner<2,Dynamic>(2,2):
 1365180540   336465782
  304089172 -1868760786
Now the matrix m is:
 1804289383 -1550966999           0           0
 -465790871 -1122281286           0           0
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
See also
class Block

◆ topRightCorner() [3/3]

template<typename Derived>
template<typename NRowsType, typename NColsType>
FixedBlockXpr<...,... >::Type Eigen::DenseBase< Derived >::topRightCorner ( NRowsType cRows,
NColsType cCols )
inline
Returns
a expression of a top-right corner of *this with either dynamic or fixed sizes.
Parameters
cRowsthe number of rows in the corner
cColsthe number of columns in the corner
Template Parameters
NRowsTypethe type of the value handling the number of rows in the block, typically Index.
NColsTypethe type of the value handling the number of columns in the block, typically Index.

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:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is m.topRightCorner(2, 2):
 1365180540   336465782
  304089172 -1868760786
Now the matrix m is:
 1804289383 -1550966999           0           0
 -465790871 -1122281286           0           0
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929

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.

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

◆ topRows() [1/2]

template<typename Derived>
template<int N>
NRowsBlockXpr< N >::Type Eigen::DenseBase< Derived >::topRows ( Index n = N)
inline
Returns
a block consisting of the top rows of *this.
Template Parameters
Nthe number of rows in the block as specified at compile-time
Parameters
nthe number of rows in the block as specified at run-time

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:

Here is the array a:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is a.topRows<2>():
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
Now the array a is:
          0           0           0           0
          0           0           0           0
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
See also
block(Index,Index,NRowsType,NColsType), class Block

◆ topRows() [2/2]

template<typename Derived>
template<typename NRowsType>
NRowsBlockXpr<... >::Type Eigen::DenseBase< Derived >::topRows ( NRowsType n)
inline
Returns
a block consisting of the top rows of *this.
Parameters
nthe number of rows in the block
Template Parameters
NRowsTypethe type of the value handling the number of rows in the block, typically Index.

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:
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929
Here is a.topRows(2):
 1804289383 -1550966999  1365180540   336465782
 -465790871 -1122281286   304089172 -1868760786
Now the array a is:
          0           0           0           0
          0           0           0           0
 -189735855 -1364114958    35005211    -2309581
  719885386  2044897763 -1852781081  1101513929

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.

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

◆ transpose() [1/2]

template<typename Derived>
DenseBase< Derived >::TransposeReturnType Eigen::DenseBase< Derived >::transpose ( )
inline
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;
Matrix< int, 2, 2 > Matrix2i
2×2 matrix of type int.
Definition Matrix.h:477

Output:

Here is the matrix m:
1804289383 -189735855
-465790871  719885386
Here is the transpose of m:
1804289383 -465790871
-189735855  719885386
Here is the coefficient (1,0) in the transpose of m:
-189735855
Let us overwrite this coefficient with the value 0.
Now the matrix m is:
1804289383          0
-465790871  719885386
Warning
If you want to replace a matrix by its own transpose, do NOT do this:
m = m.transpose(); // bug!!! caused by aliasing effect
Instead, use the transposeInPlace() method:
m.transposeInPlace();
which gives Eigen good opportunities for optimization, or alternatively you can also do:
m = m.transpose().eval();
See also
transposeInPlace(), adjoint()

◆ transpose() [2/2]

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

This is the const version of transpose().

Make sure you read the warning for transpose() !

See also
transposeInPlace(), adjoint()

◆ transposeInPlace()

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

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().

Note
if the matrix is not square, then *this must be a resizable matrix. This excludes (non-square) fixed-size matrices, block-expressions and maps.
See also
transpose(), adjoint(), adjointInPlace()

◆ unaryExpr()

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

Apply a unary operator coefficient-wise.

Parameters
[in]funcFunctor implementing the unary operator
Template Parameters
CustomUnaryOpType of func
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>
// define function to be applied coefficient-wise
double ramp(double x) {
if (x > 0)
return x;
else
return 0;
}
int main(int, char**) {
Eigen::Matrix4d m1 = Eigen::Matrix4d::Random();
std::cout << m1 << std::endl << "becomes: " << std::endl << m1.unaryExpr(std::ptr_fun(ramp)) << std::endl;
return 0;
}
Matrix< double, 4, 4 > Matrix4d
4×4 matrix of type double.
Definition Matrix.h:479

Output:

   0.696    -0.47   0.0241    0.134
   0.205    0.928   0.0723    -0.16
  -0.415    0.445    0.432 -0.00986
   0.334   -0.633   -0.046   -0.498
becomes: 
 0.696      0 0.0241  0.134
 0.205  0.928 0.0723      0
     0  0.445  0.432      0
 0.334      0      0      0

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

Example:

#include <Eigen/Core>
#include <iostream>
// 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**) {
Eigen::Matrix4d m1 = Eigen::Matrix4d::Random();
std::cout << m1 << std::endl
<< "becomes: " << std::endl
<< m1.unaryExpr(CwiseClampOp<double>(-0.5, 0.5)) << std::endl;
return 0;
}
internal::traits< Derived >::Scalar Scalar
Definition DenseBase.h:62

Output:

   0.696    -0.47   0.0241    0.134
   0.205    0.928   0.0723    -0.16
  -0.415    0.445    0.432 -0.00986
   0.334   -0.633   -0.046   -0.498
becomes: 
     0.5    -0.47   0.0241    0.134
   0.205      0.5   0.0723    -0.16
  -0.415    0.445    0.432 -0.00986
   0.334     -0.5   -0.046   -0.498
See also
unaryViewExpr, binaryExpr, class CwiseUnaryOp

◆ unaryViewExpr() [1/2]

template<typename Derived>
template<typename CustomViewOp>
CwiseUnaryView< CustomViewOp, Derived > Eigen::DenseBase< Derived >::unaryViewExpr ( const CustomViewOp & func = CustomViewOp())
inline
Returns
a non-const expression of a custom coefficient-wise unary view func of *this

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

See also
unaryExpr, binaryExpr class CwiseUnaryOp

◆ unaryViewExpr() [2/2]

template<typename Derived>
template<typename CustomViewOp>
const CwiseUnaryView< CustomViewOp, const Derived > Eigen::DenseBase< Derived >::unaryViewExpr ( const CustomViewOp & func = CustomViewOp()) const
inline
Returns
a const 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>
// 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**) {
Eigen::Matrix4d m1 = Eigen::Matrix4d::Random();
std::cout << m1 << std::endl
<< "becomes: " << std::endl
<< m1.unaryExpr(CwiseClampOp<double>(-0.5, 0.5)) << std::endl;
return 0;
}

Output:

   0.696    -0.47   0.0241    0.134
   0.205    0.928   0.0723    -0.16
  -0.415    0.445    0.432 -0.00986
   0.334   -0.633   -0.046   -0.498
becomes: 
     0.5    -0.47   0.0241    0.134
   0.205      0.5   0.0723    -0.16
  -0.415    0.445    0.432 -0.00986
   0.334     -0.5   -0.046   -0.498
See also
unaryExpr, binaryExpr class CwiseUnaryOp

◆ value()

template<typename Derived>
CoeffReturnType Eigen::DenseBase< Derived >::value ( ) const
inline
Returns
the unique coefficient of a 1x1 expression

◆ visit()

template<typename Derived>
template<typename Visitor>
void Eigen::DenseBase< Derived >::visit ( Visitor & visitor) 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);
};
CoeffReturnType value() const
Definition DenseBase.h:481
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
Note
compared to one or two for loops, visitors offer automatic unrolling for small fixed size matrix.
if the matrix is empty, then the visitor is left unchanged.
See also
minCoeff(Index*,Index*), maxCoeff(Index*,Index*), DenseBase::redux()

◆ Zero() [1/3]

template<typename Derived>
const DenseBase< Derived >::ZeroReturnType Eigen::DenseBase< Derived >::Zero ( )
inlinestatic
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
See also
Zero(Index), Zero(Index,Index)

◆ Zero() [2/3]

template<typename Derived>
const DenseBase< Derived >::ZeroReturnType Eigen::DenseBase< Derived >::Zero ( Index rows,
Index cols )
inlinestatic
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
See also
Zero(), Zero(Index)

◆ Zero() [3/3]

template<typename Derived>
const DenseBase< Derived >::ZeroReturnType Eigen::DenseBase< Derived >::Zero ( Index size)
inlinestatic
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
See also
Zero(), Zero(Index,Index)

Friends And Related Symbol Documentation

◆ operator<<()

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

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.

See also
DenseBase::format()

The documentation for this class was generated from the following files: