Eigen  3.4.90 (git rev a4098ac676528a83cfb73d4d26ce1b42ec05f47c)

Detailed Description

This is the main module of Eigen providing dense matrix and vector support (both fixed and dynamic size) with all the features corresponding to a BLAS library and much more...

#include <Eigen/Core>

Modules

 Global array typedefs
 
 Global matrix typedefs
 
 Flags
 
 Enumerations
 

Namespaces

namespace  Eigen::indexing
 
namespace  Eigen::symbolic
 

Classes

class  Eigen::aligned_allocator< T >
 STL compatible allocator to use with types requiring a non-standard alignment. More...
 
class  Eigen::ArithmeticSequence< FirstType, SizeType, IncrType >
 
class  Eigen::Array< Scalar_, Rows_, Cols_, Options_, MaxRows_, MaxCols_ >
 General-purpose arrays with easy API for coefficient-wise operations. More...
 
class  Eigen::ArrayBase< Derived >
 Base class for all 1D and 2D array, and related expressions. More...
 
class  Eigen::ArrayWrapper< ExpressionType >
 Expression of a mathematical vector or matrix as an array object. More...
 
class  Eigen::symbolic::BaseExpr< Derived >
 
class  Eigen::Block< XprType, BlockRows, BlockCols, InnerPanel >
 Expression of a fixed-size or dynamic-size block. More...
 
class  Eigen::CommaInitializer< XprType >
 Helper class used by the comma initializer operator. More...
 
class  Eigen::CwiseBinaryOp< BinaryOp, LhsType, RhsType >
 Generic expression where a coefficient-wise binary operator is applied to two expressions. More...
 
class  Eigen::CwiseNullaryOp< NullaryOp, PlainObjectType >
 Generic expression of a matrix where all coefficients are defined by a functor. More...
 
class  Eigen::CwiseTernaryOp< TernaryOp, Arg1Type, Arg2Type, Arg3Type >
 Generic expression where a coefficient-wise ternary operator is applied to two expressions. More...
 
class  Eigen::CwiseUnaryOp< UnaryOp, XprType >
 Generic expression where a coefficient-wise unary operator is applied to an expression. More...
 
class  Eigen::CwiseUnaryView< ViewOp, MatrixType >
 Generic lvalue expression of a coefficient-wise unary operator of a matrix or a vector. More...
 
class  Eigen::DenseBase< Derived >
 Base class for all dense matrices, vectors, and arrays. More...
 
class  Eigen::DenseCoeffsBase< Derived, DirectAccessors >
 Base class providing direct read-only coefficient access to matrices and arrays. More...
 
class  Eigen::DenseCoeffsBase< Derived, DirectWriteAccessors >
 Base class providing direct read/write coefficient access to matrices and arrays. More...
 
class  Eigen::DenseCoeffsBase< Derived, ReadOnlyAccessors >
 Base class providing read-only coefficient access to matrices and arrays. More...
 
class  Eigen::DenseCoeffsBase< Derived, WriteAccessors >
 Base class providing read/write coefficient access to matrices and arrays. More...
 
class  Eigen::Diagonal< MatrixType, DiagIndex_ >
 Expression of a diagonal/subdiagonal/superdiagonal in a matrix. More...
 
class  Eigen::DiagonalMatrix< Scalar_, SizeAtCompileTime, MaxSizeAtCompileTime >
 Represents a diagonal matrix with its storage. More...
 
class  Eigen::DiagonalWrapper< DiagonalVectorType_ >
 Expression of a diagonal matrix. More...
 
class  Eigen::EigenBase< Derived >
 
class  Eigen::ForceAlignedAccess< ExpressionType >
 Enforce aligned packet loads and stores regardless of what is requested. More...
 
class  Eigen::IndexedView< XprType, RowIndices, ColIndices >
 Expression of a non-sequential sub-matrix defined by arbitrary sequences of row and column indices. More...
 
class  Eigen::IOFormat
 Stores a set of parameters controlling the way matrices are printed. More...
 
class  Eigen::Map< PlainObjectType, MapOptions, StrideType >
 A matrix or vector expression mapping an existing array of data. More...
 
class  Eigen::MapBase< Derived, ReadOnlyAccessors >
 Base class for dense Map and Block expression with direct access. More...
 
class  Eigen::MapBase< Derived, WriteAccessors >
 Base class for non-const dense Map and Block expression with direct access. More...
 
class  Eigen::Matrix< Scalar_, Rows_, Cols_, Options_, MaxRows_, MaxCols_ >
 The matrix class, also used for vectors and row-vectors. More...
 
class  Eigen::MatrixBase< Derived >
 Base class for all dense matrices, vectors, and expressions. More...
 
class  Eigen::MatrixWrapper< ExpressionType >
 Expression of an array as a mathematical vector or matrix. More...
 
class  Eigen::NestByValue< ExpressionType >
 Expression which must be nested by value. More...
 
class  Eigen::NoAlias< ExpressionType, StorageBase >
 Pseudo expression providing an operator = assuming no aliasing. More...
 
class  Eigen::NumTraits< T >
 Holds information about the various numeric (i.e. scalar) types allowed by Eigen. More...
 
class  Eigen::PartialReduxExpr< MatrixType, MemberOp, Direction >
 Generic expression of a partially reduxed matrix. More...
 
class  Eigen::PermutationBase< Derived >
 Base class for permutations. More...
 
class  Eigen::PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, StorageIndex_ >
 Permutation matrix. More...
 
class  Eigen::PermutationWrapper< IndicesType_ >
 Class to view a vector of integers as a permutation matrix. More...
 
class  Eigen::PlainObjectBase< Derived >
 Dense storage base class for matrices and arrays. More...
 
class  Eigen::Product< Lhs_, Rhs_, Option >
 Expression of the product of two arbitrary matrices or vectors. More...
 
class  Eigen::Ref< PlainObjectType, Options, StrideType >
 A matrix or vector expression mapping an existing expression. More...
 
class  Eigen::Replicate< MatrixType, RowFactor, ColFactor >
 Expression of the multiple replication of a matrix or vector. More...
 
class  Eigen::Reshaped< XprType, Rows, Cols, Order >
 Expression of a fixed-size or dynamic-size reshape. More...
 
class  Eigen::Reverse< MatrixType, Direction >
 Expression of the reverse of a vector or matrix. More...
 
class  Eigen::ScalarBinaryOpTraits< ScalarA, ScalarB, BinaryOp >
 Determines whether the given binary operation of two numeric types is allowed and what the scalar return type is. More...
 
class  Eigen::Select< ConditionMatrixType, ThenMatrixType, ElseMatrixType >
 Expression of a coefficient wise version of the C++ ternary operator ?: More...
 
class  Eigen::SelfAdjointView< MatrixType_, UpLo >
 Expression of a selfadjoint matrix from a triangular part of a dense matrix. More...
 
class  Eigen::Solve< Decomposition, RhsType >
 Pseudo expression representing a solving operation. More...
 
class  Eigen::Stride< OuterStrideAtCompileTime_, InnerStrideAtCompileTime_ >
 Holds strides information for Map. More...
 
class  Eigen::Transpose< MatrixType >
 Expression of the transpose of a matrix. More...
 
class  Eigen::Transpositions< SizeAtCompileTime, MaxSizeAtCompileTime, StorageIndex_ >
 Represents a sequence of transpositions (row/column interchange) More...
 
class  Eigen::TriangularBase< Derived >
 Base class for triangular part in a matrix. More...
 
class  Eigen::TriangularView< MatrixType_, Mode_ >
 Expression of a triangular part in a matrix. More...
 
class  Eigen::TriangularViewImpl< MatrixType_, Mode_, Dense >
 Base class for a triangular part in a dense matrix. More...
 
class  Eigen::VectorBlock< VectorType, Size >
 Expression of a fixed-size or dynamic-size sub-vector. More...
 
class  Eigen::VectorwiseOp< ExpressionType, Direction >
 Pseudo expression providing broadcasting and partial reduction operations. More...
 
class  Eigen::WithFormat< ExpressionType >
 Pseudo expression providing matrix output with given format. More...
 

Functions

template<int N>
static const auto Eigen::fix ()
 
template<int N>
static const auto Eigen::fix (int val)
 

Variables

static const Eigen::internal::all_t Eigen::placeholders::all
 
static const lastp1_t Eigen::placeholders::end
 
static const last_t Eigen::placeholders::last
 
static const auto Eigen::placeholders::lastp1
 

Function Documentation

◆ fix() [1/2]

template<int N>
Eigen::fix< N > ( )
static

This identifier permits to construct an object embedding a compile-time integer N.

Template Parameters
Nthe compile-time integer value

It is typically used in conjunction with the Eigen::seq and Eigen::seqN functions to pass compile-time values to them:

seqN(10,fix<4>,fix<-3>) // <=> [10 7 4 1]
ArithmeticSequence< typename internal::cleanup_index_type< FirstType >::type, typename internal::cleanup_index_type< SizeType >::type, typename internal::cleanup_seq_incr< IncrType >::type > seqN(FirstType first, SizeType size, IncrType incr)
Definition: ArithmeticSequence.h:99

See also the function fix(int) to pass both a compile-time and runtime value.

In c++14, it is implemented as:

template<int N> static const internal::FixedInt<N> fix{};
static const auto fix()

where internal::FixedInt<N> is an internal template class similar to std::integral_constant <int,N> Here, fix<N> is thus an object of type internal::FixedInt<N>.

In c++98/11, it is implemented as a function:

template<int N> inline internal::FixedInt<N> fix();

Here internal::FixedInt<N> is thus a pointer to function.

If for some reason you want a true object in c++98 then you can write:

fix<N>()

which is also valid in c++14.

See also
fix<N>(int), seq, seqN

◆ fix() [2/2]

template<int N>
Eigen::fix< N > ( int  val)
static

This function returns an object embedding both a compile-time integer N, and a fallback runtime value val.

Template Parameters
Nthe compile-time integer value
Parameters
valthe fallback runtime integer value

This function is a more general version of the fix identifier/function that can be used in template code where the compile-time value could turn out to actually mean "undefined at compile-time". For positive integers such as a size or a dimension, this case is identified by Eigen::Dynamic, whereas runtime signed integers (e.g., an increment/stride) are identified as Eigen::DynamicIndex. In such a case, the runtime value val will be used as a fallback.

A typical use case would be:

template<typename Derived> void foo(const MatrixBase<Derived> &mat) {
const int N = Derived::RowsAtCompileTime==Dynamic ? Dynamic : Derived::RowsAtCompileTime/2;
const int n = mat.rows()/2;
... mat( seqN(0,fix<N>(n) ) ...;
}
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:62
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
const int Dynamic
Definition: Constants.h:24

In this example, the function Eigen::seqN knows that the second argument is expected to be a size. If the passed compile-time value N equals Eigen::Dynamic, then the proxy object returned by fix will be dissmissed, and converted to an Eigen::Index of value n. Otherwise, the runtime-value n will be dissmissed, and the returned ArithmeticSequence will be of the exact same type as seqN(0,fix<N>) .

See also
fix, seqN, class ArithmeticSequence

Variable Documentation

◆ all

Eigen::placeholders::all
static

Can be used as a parameter to DenseBase::operator()(const RowIndices&, const ColIndices&) to index all rows or columns

◆ end

Eigen::placeholders::end
static
See also
lastp1

◆ last

Eigen::placeholders::last
static

Can be used as a parameter to Eigen::seq and Eigen::seqN functions to symbolically reference the last element/row/columns of the underlying vector or matrix once passed to DenseBase::operator()(const RowIndices&, const ColIndices&).

This symbolic placeholder supports standard arithmetic operations.

A typical usage example would be:

using namespace Eigen;
VectorXd v(n);
v(seq(2,last-2)).setOnes();
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182
static const last_t last
Definition: IndexedViewHelper.h:44
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
auto seq(FirstType f, LastType l, IncrType incr)
See also
end

◆ lastp1

Eigen::placeholders::lastp1
static

Can be used as a parameter to Eigen::seq and Eigen::seqN functions to symbolically reference the last+1 element/row/columns of the underlying vector or matrix once passed to DenseBase::operator()(const RowIndices&, const ColIndices&).

This symbolic placeholder supports standard arithmetic operations. It is essentially an alias to last+fix<1>.

See also
last