![]() |
Eigen
3.4.90 (git rev 9589cc4e7fd8e4538bedef80dd36c7738977a8be)
|
#include <Eigen/src/Core/MatrixBase.h>
Base class for all dense matrices, vectors, and expressions.
This class is the base that is inherited by all matrix, vector, and related expression types. Most of the Eigen API is contained in this class, and its base classes. Other important classes for the Eigen API are Matrix, and VectorwiseOp.
Note that some methods are defined in other modules such as the LU module LU module for all functions related to matrix inversions.
Derived | is the derived type, e.g. a matrix type, or an expression, etc. |
When writing a function taking Eigen objects as argument, if you want your function to take as argument any matrix, vector, or expression, just let it take a MatrixBase argument. As an example, here is a function printFirstRow which, given a matrix, vector, or expression x, prints the first row of x.
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_MATRIXBASE_PLUGIN
.
Public Member Functions | |
const MatrixFunctionReturnValue< Derived > | acosh () const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise inverse hyperbolic cosine use ArrayBase::acosh . | |
const AdjointReturnType | adjoint () const |
void | adjointInPlace () |
template<typename EssentialPart> | |
void | applyHouseholderOnTheLeft (const EssentialPart &essential, const Scalar &tau, Scalar *workspace) |
template<typename EssentialPart> | |
void | applyHouseholderOnTheRight (const EssentialPart &essential, const Scalar &tau, Scalar *workspace) |
template<typename OtherDerived> | |
void | applyOnTheLeft (const EigenBase< OtherDerived > &other) |
template<typename OtherScalar> | |
void | applyOnTheLeft (Index p, Index q, const JacobiRotation< OtherScalar > &j) |
template<typename OtherDerived> | |
void | applyOnTheRight (const EigenBase< OtherDerived > &other) |
template<typename OtherScalar> | |
void | applyOnTheRight (Index p, Index q, const JacobiRotation< OtherScalar > &j) |
ArrayWrapper< Derived > | array () |
const ArrayWrapper< const Derived > | array () const |
const DiagonalWrapper< const Derived > | asDiagonal () const |
const MatrixFunctionReturnValue< Derived > | asinh () const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise inverse hyperbolic sine use ArrayBase::asinh . | |
const SkewSymmetricWrapper< const Derived > | asSkewSymmetric () const |
const MatrixFunctionReturnValue< Derived > | atanh () const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise inverse hyperbolic cosine use ArrayBase::atanh . | |
template<int Options> | |
BDCSVD< typename MatrixBase< Derived >::PlainObject, Options > | bdcSvd () const |
template<int Options> | |
BDCSVD< typename MatrixBase< Derived >::PlainObject, Options > | bdcSvd (unsigned int computationOptions) const |
template<typename CustomBinaryOp, typename OtherDerived> | |
const CwiseBinaryOp< CustomBinaryOp, const Derived, const OtherDerived > | binaryExpr (const Eigen::MatrixBase< OtherDerived > &other, const CustomBinaryOp &func=CustomBinaryOp()) const |
RealScalar | blueNorm () const |
Matrix< Scalar, 3, 1 > | canonicalEulerAngles (Index a0, Index a1, Index a2) const |
template<typename PermutationIndexType> | |
const ColPivHouseholderQR< typename MatrixBase< Derived >::PlainObject, PermutationIndexType > | colPivHouseholderQr () const |
template<typename PermutationIndex> | |
const CompleteOrthogonalDecomposition< typename MatrixBase< Derived >::PlainObject, PermutationIndex > | completeOrthogonalDecomposition () const |
template<typename ResultType> | |
void | computeInverseAndDetWithCheck (ResultType &inverse, typename ResultType::Scalar &determinant, bool &invertible, const RealScalar &absDeterminantThreshold=NumTraits< Scalar >::dummy_precision()) const |
template<typename ResultType> | |
void | computeInverseWithCheck (ResultType &inverse, bool &invertible, const RealScalar &absDeterminantThreshold=NumTraits< Scalar >::dummy_precision()) const |
const MatrixFunctionReturnValue< Derived > | cos () const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise cosine use ArrayBase::cos . | |
const MatrixFunctionReturnValue< Derived > | cosh () const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise hyperbolic cosine use ArrayBase::cosh . | |
template<typename OtherDerived> | |
internal::cross_impl< Derived, OtherDerived >::return_type | cross (const MatrixBase< OtherDerived > &other) const |
template<typename OtherDerived> | |
PlainObject | cross3 (const MatrixBase< OtherDerived > &other) const |
const CwiseAbsReturnType | cwiseAbs () const |
const CwiseAbs2ReturnType | cwiseAbs2 () const |
const CwiseArgReturnType | cwiseArg () const |
const CwiseCbrtReturnType | cwiseCbrt () const |
template<typename OtherDerived> | |
const CwiseBinaryEqualReturnType< OtherDerived > | cwiseEqual (const Eigen::MatrixBase< OtherDerived > &other) const |
const CwiseScalarEqualReturnType | cwiseEqual (const Scalar &s) const |
template<typename OtherDerived> | |
const CwiseBinaryGreaterReturnType< OtherDerived > | cwiseGreater (const Eigen::MatrixBase< OtherDerived > &other) const |
const CwiseScalarGreaterReturnType | cwiseGreater (const Scalar &s) const |
template<typename OtherDerived> | |
const CwiseBinaryGreaterOrEqualReturnType< OtherDerived > | cwiseGreaterOrEqual (const Eigen::MatrixBase< OtherDerived > &other) const |
const CwiseScalarGreaterOrEqualReturnType | cwiseGreaterOrEqual (const Scalar &s) const |
const CwiseInverseReturnType | cwiseInverse () const |
template<typename OtherDerived> | |
const CwiseBinaryLessReturnType< OtherDerived > | cwiseLess (const Eigen::MatrixBase< OtherDerived > &other) const |
const CwiseScalarLessReturnType | cwiseLess (const Scalar &s) const |
template<typename OtherDerived> | |
const CwiseBinaryLessOrEqualReturnType< OtherDerived > | cwiseLessOrEqual (const Eigen::MatrixBase< OtherDerived > &other) const |
const CwiseScalarLessOrEqualReturnType | cwiseLessOrEqual (const Scalar &s) const |
template<int NaNPropagation = PropagateFast, typename OtherDerived> | |
const CwiseBinaryOp< internal::scalar_max_op< Scalar, Scalar, NaNPropagation >, const Derived, const OtherDerived > | cwiseMax (const Eigen::MatrixBase< OtherDerived > &other) const |
template<int NaNPropagation = PropagateFast> | |
const CwiseBinaryOp< internal::scalar_max_op< Scalar, Scalar, NaNPropagation >, const Derived, const ConstantReturnType > | cwiseMax (const Scalar &other) const |
template<int NaNPropagation = PropagateFast, typename OtherDerived> | |
const CwiseBinaryOp< internal::scalar_min_op< Scalar, Scalar, NaNPropagation >, const Derived, const OtherDerived > | cwiseMin (const Eigen::MatrixBase< OtherDerived > &other) const |
template<int NaNPropagation = PropagateFast> | |
const CwiseBinaryOp< internal::scalar_min_op< Scalar, Scalar, NaNPropagation >, const Derived, const ConstantReturnType > | cwiseMin (const Scalar &other) const |
template<typename OtherDerived> | |
const CwiseBinaryNotEqualReturnType< OtherDerived > | cwiseNotEqual (const Eigen::MatrixBase< OtherDerived > &other) const |
const CwiseScalarNotEqualReturnType | cwiseNotEqual (const Scalar &s) const |
template<typename OtherDerived> | |
const CwiseBinaryOp< internal::scalar_product_op< Derived ::Scalar, OtherDerived ::Scalar >, const Derived, const OtherDerived > | cwiseProduct (const Eigen::MatrixBase< OtherDerived > &other) const |
template<typename OtherDerived> | |
const CwiseBinaryOp< internal::scalar_quotient_op< Scalar >, const Derived, const OtherDerived > | cwiseQuotient (const Eigen::MatrixBase< OtherDerived > &other) const |
const CwiseSignReturnType | cwiseSign () const |
const CwiseSqrtReturnType | cwiseSqrt () const |
const CwiseSquareReturnType | cwiseSquare () const |
Scalar | determinant () const |
template<int Index_> | |
Diagonal< Derived, Index_ > | diagonal () |
DiagonalReturnType | diagonal () |
template<int Index_> | |
const Diagonal< const Derived, Index_ > | diagonal () const |
const ConstDiagonalReturnType | diagonal () const |
Diagonal< Derived, DynamicIndex > | diagonal (Index index) |
const Diagonal< const Derived, DynamicIndex > | diagonal (Index index) const |
Index | diagonalSize () const |
template<typename OtherDerived> | |
ScalarBinaryOpTraits< typenameinternal::traits< Derived >::Scalar, typenameinternal::traits< OtherDerived >::Scalar >::ReturnType | dot (const MatrixBase< OtherDerived > &other) const |
EigenvaluesReturnType | eigenvalues () const |
Computes the eigenvalues of a matrix. | |
EIGEN_DEPRECATED Matrix< Scalar, 3, 1 > | eulerAngles (Index a0, Index a1, Index a2) const |
const MatrixExponentialReturnValue< Derived > | exp () const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise exponential use ArrayBase::exp . | |
Derived & | forceAlignedAccess () |
const Derived & | forceAlignedAccess () const |
template<bool Enable> | |
std::conditional_t< Enable, ForceAlignedAccess< Derived >, Derived & > | forceAlignedAccessIf () |
template<bool Enable> | |
add_const_on_value_type_t< std::conditional_t< Enable, ForceAlignedAccess< Derived >, Derived & > > | forceAlignedAccessIf () const |
template<typename PermutationIndex> | |
const FullPivHouseholderQR< typename MatrixBase< Derived >::PlainObject, PermutationIndex > | fullPivHouseholderQr () const |
template<typename PermutationIndex> | |
const FullPivLU< typename MatrixBase< Derived >::PlainObject, PermutationIndex > | fullPivLu () const |
const HNormalizedReturnType | hnormalized () const |
homogeneous normalization | |
HomogeneousReturnType | homogeneous () const |
const HouseholderQR< PlainObject > | householderQr () const |
RealScalar | hypotNorm () const |
const Inverse< Derived > | inverse () const |
bool | isDiagonal (const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const |
bool | isIdentity (const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const |
bool | isLowerTriangular (const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const |
template<typename OtherDerived> | |
bool | isOrthogonal (const MatrixBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const |
bool | isSkewSymmetric (const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const |
bool | isUnitary (const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const |
bool | isUpperTriangular (const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const |
template<int Options> | |
JacobiSVD< typename MatrixBase< Derived >::PlainObject, Options > | jacobiSvd () const |
template<typename OtherDerived> | |
const Product< Derived, OtherDerived, LazyProduct > | lazyProduct (const MatrixBase< OtherDerived > &other) const |
const LDLT< PlainObject > | ldlt () const |
const LLT< PlainObject > | llt () const |
const MatrixLogarithmReturnValue< Derived > | log () const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise logarithm use ArrayBase::log . | |
template<int p> | |
RealScalar | lpNorm () const |
template<typename PermutationIndex> | |
const PartialPivLU< typename MatrixBase< Derived >::PlainObject, PermutationIndex > | lu () const |
template<typename EssentialPart> | |
void | makeHouseholder (EssentialPart &essential, Scalar &tau, RealScalar &beta) const |
void | makeHouseholderInPlace (Scalar &tau, RealScalar &beta) |
const MatrixFunctionReturnValue< Derived > | matrixFunction (StemFunction f) const |
Helper function for the unsupported MatrixFunctions module. | |
NoAlias< Derived, Eigen::MatrixBase > | noalias () |
RealScalar | norm () const |
void | normalize () |
const PlainObject | normalized () const |
template<typename OtherDerived> | |
bool | operator!= (const MatrixBase< OtherDerived > &other) const |
template<typename OtherDerived> | |
const CwiseBinaryOp< internal::scalar_bitwise_and_op< Scalar >, const Derived, const OtherDerived > | operator& (const Eigen::MatrixBase< OtherDerived > &other) const |
template<typename OtherDerived> | |
const CwiseBinaryOp< internal::scalar_boolean_and_op< Scalar >, const Derived, const OtherDerived > | operator&& (const Eigen::MatrixBase< OtherDerived > &other) const |
template<typename DiagonalDerived> | |
const Product< Derived, DiagonalDerived, LazyProduct > | operator* (const DiagonalBase< DiagonalDerived > &diagonal) const |
template<typename OtherDerived> | |
const Product< Derived, OtherDerived > | operator* (const MatrixBase< OtherDerived > &other) const |
template<typename SkewDerived> | |
const Product< Derived, SkewDerived, LazyProduct > | operator* (const SkewSymmetricBase< SkewDerived > &skew) const |
template<typename OtherDerived> | |
Derived & | operator*= (const EigenBase< OtherDerived > &other) |
template<typename OtherDerived> | |
const CwiseBinaryOp< sum< Scalar >, const Derived, const OtherDerived > | operator+ (const Eigen::MatrixBase< OtherDerived > &other) const |
template<typename OtherDerived> | |
Derived & | operator+= (const MatrixBase< OtherDerived > &other) |
template<typename OtherDerived> | |
const CwiseBinaryOp< difference< Scalar >, const Derived, const OtherDerived > | operator- (const Eigen::MatrixBase< OtherDerived > &other) const |
template<typename OtherDerived> | |
Derived & | operator-= (const MatrixBase< OtherDerived > &other) |
Derived & | operator= (const MatrixBase &other) |
template<typename OtherDerived> | |
bool | operator== (const MatrixBase< OtherDerived > &other) const |
template<typename OtherDerived> | |
const CwiseBinaryOp< internal::scalar_bitwise_xor_op< Scalar >, const Derived, const OtherDerived > | operator^ (const Eigen::MatrixBase< OtherDerived > &other) const |
RealScalar | operatorNorm () const |
Computes the L2 operator norm. | |
template<typename OtherDerived> | |
const CwiseBinaryOp< internal::scalar_bitwise_or_op< Scalar >, const Derived, const OtherDerived > | operator| (const Eigen::MatrixBase< OtherDerived > &other) const |
template<typename OtherDerived> | |
const CwiseBinaryOp< internal::scalar_boolean_or_op< Scalar >, const Derived, const OtherDerived > | operator|| (const Eigen::MatrixBase< OtherDerived > &other) const |
template<typename PermutationIndex> | |
const PartialPivLU< typename MatrixBase< Derived >::PlainObject, PermutationIndex > | partialPivLu () const |
const MatrixPowerReturnValue< Derived > | pow (const RealScalar &p) const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise power to p use ArrayBase::pow . | |
const MatrixComplexPowerReturnValue< Derived > | pow (const std::complex< RealScalar > &p) const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise power to p use ArrayBase::pow . | |
template<unsigned int UpLo> | |
MatrixBase< Derived >::template SelfAdjointViewReturnType< UpLo >::Type | selfadjointView () |
template<unsigned int UpLo> | |
MatrixBase< Derived >::template ConstSelfAdjointViewReturnType< UpLo >::Type | selfadjointView () const |
Derived & | setIdentity () |
Derived & | setIdentity (Index rows, Index cols) |
Resizes to the given size, and writes the identity expression (not necessarily square) into *this. | |
Derived & | setUnit (Index i) |
Set the coefficients of *this to the i-th unit (basis) vector. | |
Derived & | setUnit (Index newSize, Index i) |
Resizes to the given newSize, and writes the i-th unit (basis) vector into *this. | |
const MatrixFunctionReturnValue< Derived > | sin () const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise sine use ArrayBase::sin . | |
const MatrixFunctionReturnValue< Derived > | sinh () const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise hyperbolic sine use ArrayBase::sinh . | |
const SparseView< Derived > | sparseView (const Scalar &m_reference=Scalar(0), const typename NumTraits< Scalar >::Real &m_epsilon=NumTraits< Scalar >::dummy_precision()) const |
const MatrixSquareRootReturnValue< Derived > | sqrt () const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise square root use ArrayBase::sqrt . | |
RealScalar | squaredNorm () const |
RealScalar | stableNorm () const |
void | stableNormalize () |
const PlainObject | stableNormalized () const |
Scalar | trace () const |
template<unsigned int Mode> | |
MatrixBase< Derived >::template TriangularViewReturnType< Mode >::Type | triangularView () |
template<unsigned int Mode> | |
MatrixBase< Derived >::template ConstTriangularViewReturnType< Mode >::Type | triangularView () const |
PlainObject | unitOrthogonal (void) const |
![]() | |
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, Dynamic > | replicate (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 |
![]() | |
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 IdentityReturnType | Identity () |
static const IdentityReturnType | Identity (Index rows, Index cols) |
static const BasisReturnType | Unit (Index i) |
static const BasisReturnType | Unit (Index size, Index i) |
static const BasisReturnType | UnitW () |
static const BasisReturnType | UnitX () |
static const BasisReturnType | UnitY () |
static const BasisReturnType | UnitZ () |
![]() | |
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, PlainObject > | NullaryExpr (const CustomNullaryOp &func) |
template<typename CustomNullaryOp> | |
static const CwiseNullaryOp< CustomNullaryOp, PlainObject > | NullaryExpr (Index rows, Index cols, const CustomNullaryOp &func) |
template<typename CustomNullaryOp> | |
static const CwiseNullaryOp< CustomNullaryOp, PlainObject > | NullaryExpr (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) |
Additional Inherited Members | |
![]() | |
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 >::MaxColsAtCompileTime > | PlainArray |
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 > | PlainMatrix |
typedef std::conditional_t< internal::is_same< typename internal::traits< Derived >::XprKind, MatrixXpr >::value, PlainMatrix, PlainArray > | PlainObject |
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 |
![]() | |
constexpr | DenseBase ()=default |
![]() | |
template<typename Derived> | |
std::ostream & | operator<< (std::ostream &s, const DenseBase< Derived > &m) |
const MatrixFunctionReturnValue< Derived > Eigen::MatrixBase< Derived >::acosh | ( | ) | const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise inverse hyperbolic cosine use ArrayBase::acosh .
*this
.
|
inline |
Example:
Output:
Here is the 2x2 complex matrix m: (0.924,-0.824) (0.146,-0.237) (-0.199,-0.0532) (0.334,0.634) Here is the adjoint of m: (0.924,0.824) (-0.199,0.0532) (0.146,0.237) (0.334,-0.634)
|
inline |
This is the "in place" version of adjoint(): it replaces *this
by its own transpose. Thus, doing
has the same effect on m as doing
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 adjoint. If you just need the adjoint of a matrix, use adjoint().
*this
must be a resizable matrix. This excludes (non-square) fixed-size matrices, block-expressions and maps.void Eigen::MatrixBase< Derived >::applyHouseholderOnTheLeft | ( | const EssentialPart & | essential, |
const Scalar & | tau, | ||
Scalar * | workspace ) |
Apply the elementary reflector H given by
On input:
essential | the essential part of the vector v |
tau | the scaling factor of the Householder transformation |
workspace | a pointer to working space with at least this->cols() entries |
void Eigen::MatrixBase< Derived >::applyHouseholderOnTheRight | ( | const EssentialPart & | essential, |
const Scalar & | tau, | ||
Scalar * | workspace ) |
Apply the elementary reflector H given by
On input:
essential | the essential part of the vector v |
tau | the scaling factor of the Householder transformation |
workspace | a pointer to working space with at least this->rows() entries |
|
inline |
replaces *this
by other * *this
.
Example:
Output:
At start, A = -0.824 -0.199 0.634 0.924 -0.237 0.334 -0.0532 0.146 -0.779 After applyOnTheLeft, A = 0.924 -0.237 0.334 -0.0532 0.146 -0.779 -0.824 -0.199 0.634
|
inline |
This is defined in the Jacobi module.
Applies the rotation in the plane j to the rows p and q of *this
, i.e., it computes B = J * B, with
|
inline |
replaces *this
by *this
* other. It is equivalent to MatrixBase::operator*=().
Example:
Output:
At start, A = -0.824 -0.199 0.634 0.924 -0.237 0.334 -0.0532 0.146 -0.779 After A *= B, A = 0.634 -0.824 -0.199 0.334 0.924 -0.237 -0.779 -0.0532 0.146 After applyOnTheRight, A = -0.199 0.634 -0.824 -0.237 0.334 0.924 0.146 -0.779 -0.0532
|
inline |
This is defined in the Jacobi module.
Applies the rotation in the plane j to the columns p and q of *this
, i.e., it computes B = B * J with
|
inline |
|
inline |
|
inline |
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:
Output:
2 0 0 0 5 0 0 0 6
const MatrixFunctionReturnValue< Derived > Eigen::MatrixBase< Derived >::asinh | ( | ) | const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise inverse hyperbolic sine use ArrayBase::asinh .
*this
.
|
inline |
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.
const MatrixFunctionReturnValue< Derived > Eigen::MatrixBase< Derived >::atanh | ( | ) | const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise inverse hyperbolic cosine use ArrayBase::atanh .
*this
. BDCSVD< typename MatrixBase< Derived >::PlainObject, Options > Eigen::MatrixBase< Derived >::bdcSvd | ( | ) | const |
This is defined in the SVD module.
*this
computed by Divide & Conquer algorithmBDCSVD< typename MatrixBase< Derived >::PlainObject, Options > Eigen::MatrixBase< Derived >::bdcSvd | ( | unsigned int | computationOptions | ) | const |
This is defined in the SVD module.
*this
computed by Divide & Conquer algorithm
|
inline |
The template parameter CustomBinaryOp is the type of the functor of the custom operator (see class CwiseBinaryOp for an example)
Here is an example illustrating the use of custom functors:
Output:
(0.696,-0.727) (-0.47,-0.216) (0.0241,-0.778) (0.134,0.0072) (0.205,0.74) (0.928,0.852) (0.0723,-0.758) (-0.16,0.22) (-0.415,-0.757) (0.445,0.835) (0.432,-0.711) (-0.00986,0.859) (0.334,0.741) (-0.633,-0.834) (-0.046,-0.467) (-0.498,-0.316)
|
inline |
*this
using the Blue's algorithm. A Portable Fortran Program to Find the Euclidean Norm of a Vector, ACM TOMS, Vol 4, Issue 1, 1978.For architecture/scalar types without vectorization, this version is much faster than stableNorm(). Otherwise the stableNorm() is faster.
const ColPivHouseholderQR< typename MatrixBase< Derived >::PlainObject, PermutationIndexType > Eigen::MatrixBase< Derived >::colPivHouseholderQr | ( | ) | const |
*this
.const CompleteOrthogonalDecomposition< typename MatrixBase< Derived >::PlainObject, PermutationIndex > Eigen::MatrixBase< Derived >::completeOrthogonalDecomposition | ( | ) | const |
*this
.
|
inline |
This is defined in the LU module.
Computation of matrix inverse and determinant, with invertibility check.
This is only for fixed-size square matrices of size up to 4x4.
Notice that it will trigger a copy of input matrix when trying to do the inverse in place.
inverse | Reference to the matrix in which to store the inverse. |
determinant | Reference to the variable in which to store the determinant. |
invertible | Reference to the bool variable in which to store whether the matrix is invertible. |
absDeterminantThreshold | Optional parameter controlling the invertibility check. The matrix will be declared invertible if the absolute value of its determinant is greater than this threshold. |
Example:
Output:
Here is the matrix m: 0.696 0.334 0.445 0.205 -0.47 -0.633 -0.415 0.928 0.0241 Its determinant is 0.485 It is invertible, and its inverse is: 1.19 0.835 -0.00499 0.531 0.416 1.1 -0.00911 -1.62 -0.816
|
inline |
This is defined in the LU module.
Computation of matrix inverse, with invertibility check.
This is only for fixed-size square matrices of size up to 4x4.
Notice that it will trigger a copy of input matrix when trying to do the inverse in place.
inverse | Reference to the matrix in which to store the inverse. |
invertible | Reference to the bool variable in which to store whether the matrix is invertible. |
absDeterminantThreshold | Optional parameter controlling the invertibility check. The matrix will be declared invertible if the absolute value of its determinant is greater than this threshold. |
Example:
Output:
Here is the matrix m: 0.696 0.334 0.445 0.205 -0.47 -0.633 -0.415 0.928 0.0241 It is invertible, and its inverse is: 1.19 0.835 -0.00499 0.531 0.416 1.1 -0.00911 -1.62 -0.816
const MatrixFunctionReturnValue< Derived > Eigen::MatrixBase< Derived >::cos | ( | ) | const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise cosine use ArrayBase::cos .
*this
. const MatrixFunctionReturnValue< Derived > Eigen::MatrixBase< Derived >::cosh | ( | ) | const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise hyperbolic cosine use ArrayBase::cosh .
*this
.
|
inline |
*this
Example:
Output:
2 4 6 5 1 0
|
inline |
*this
Example:
Output:
4 16 36 25 1 0
|
inline |
*this
Example:
Output:
(0.924,-0.824) (0.146,-0.237) (0.633,-0.779) (-0.199,-0.0532) (0.334,0.634) (0.982,-0.573) -0.728 -1.02 -0.889 -2.88 1.09 -0.528
|
inline |
|
inline |
Example:
Output:
Comparing m with identity matrix: 1 1 0 1 Number of coefficients that are equal: 3
|
inline |
*this
and a scalar s
|
inline |
|
inline |
*this
and a scalar s
|
inline |
|
inline |
*this
and a scalar s
|
inline |
Example:
Output:
0.5 2 1 0.333 4 1
|
inline |
|
inline |
*this
and a scalar s
|
inline |
|
inline |
*this
and a scalar s
|
inline |
Example:
Output:
4 3 4
|
inline |
|
inline |
Example:
Output:
2 2 3
|
inline |
|
inline |
Example:
Output:
Comparing m with identity matrix: 0 0 1 0 Number of coefficients that are not equal: 1
|
inline |
*this
and a scalar s
|
inline |
Example:
Output:
a: 1804289383 719885386 -1364114958 -465790871 -1550966999 2044897763 -189735855 -1122281286 1365180540 b: 304089172 336465782 1101513929 35005211 -1868760786 1315634022 -1852781081 -2309581 -778350579 c: -1900151348 1274064924 1451757314 1178627603 -1395361186 1323254130 982755863 -343432946 1494074956
|
inline |
Example:
Output:
0.5 1.5 1.33
|
inline |
Example:
Output:
1 -1 1 -1 1 0
|
inline |
Example:
Output:
1 1.41 2
|
inline |
|
inline |
This is defined in the LU module.
|
inline |
*this
*this
is not required to be square.
The template parameter DiagIndex represent a super diagonal if DiagIndex > 0 and a sub diagonal otherwise. DiagIndex == 0 is equivalent to the main diagonal.
Example:
Output:
Here is the matrix m: 1804289383 -1550966999 1365180540 336465782 -465790871 -1122281286 304089172 -1868760786 -189735855 -1364114958 35005211 -2309581 719885386 2044897763 -1852781081 1101513929 Here are the coefficients on the 1st super-diagonal and 2nd sub-diagonal of m: -1550966999 304089172 -2309581 -189735855 2044897763
|
inline |
*this
*this
is not required to be square.
Example:
Output:
Here is the matrix m: 1804289383 719885386 -1364114958 -465790871 -1550966999 2044897763 -189735855 -1122281286 1365180540 Here are the coefficients on the main diagonal of m: 1804289383 -1550966999 1365180540
|
inline |
This is the const version of diagonal<int>().
|
inline |
This is the const version of diagonal().
|
inline |
*this
*this
is not required to be square.
The template parameter DiagIndex represent a super diagonal if DiagIndex > 0 and a sub diagonal otherwise. DiagIndex == 0 is equivalent to the main diagonal.
Example:
Output:
Here is the matrix m: 1804289383 -1550966999 1365180540 336465782 -465790871 -1122281286 304089172 -1868760786 -189735855 -1364114958 35005211 -2309581 719885386 2044897763 -1852781081 1101513929 Here are the coefficients on the 1st super-diagonal and 2nd sub-diagonal of m: -1550966999 304089172 -2309581 -189735855 2044897763
|
inline |
This is the const version of diagonal(Index).
|
inline |
|
inline |
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.
|
inline |
Computes the eigenvalues of a matrix.
This is defined in the Eigenvalues module.
This function computes the eigenvalues with the help of the EigenSolver class (for real matrices) or the ComplexEigenSolver class (for complex matrices).
The eigenvalues are repeated according to their algebraic multiplicity, so there are as many eigenvalues as rows in the matrix.
The SelfAdjointView class provides a better algorithm for selfadjoint matrices.
Example:
Output:
The eigenvalues of the 3x3 matrix of ones are: (-5.31e-17,0) (3,0) (0,0)
const MatrixExponentialReturnValue< Derived > Eigen::MatrixBase< Derived >::exp | ( | ) | const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise exponential use ArrayBase::exp .
*this
.
|
inline |
|
inline |
|
inline |
|
inline |
const FullPivHouseholderQR< typename MatrixBase< Derived >::PlainObject, PermutationIndex > Eigen::MatrixBase< Derived >::fullPivHouseholderQr | ( | ) | const |
*this
.
|
inline |
This is defined in the LU module.
*this
.
|
inline |
*this
.
|
inline |
*this
avoiding underflow and overflow. This version use a concatenation of hypot() calls, and it is very slow.
|
inlinestatic |
This variant is only for fixed-size MatrixBase types. For dynamic-size types, you need to use the variant taking size arguments.
Example:
Output:
1 0 0 0 0 1 0 0 0 0 1 0
|
inlinestatic |
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 Identity() should be used instead.
Example:
Output:
1 0 0 0 1 0 0 0 1 0 0 0
|
inline |
This is defined in the LU module.
For small fixed sizes up to 4x4, this method uses cofactors. In the general case, this method uses class PartialPivLU.
Here is the matrix m: 0.696 0.334 0.445 0.205 -0.47 -0.633 -0.415 0.928 0.0241 Its inverse is: 1.19 0.835 -0.00499 0.531 0.416 1.1 -0.00911 -1.62 -0.816
bool Eigen::MatrixBase< Derived >::isDiagonal | ( | const RealScalar & | prec = NumTraits<Scalar>::dummy_precision() | ) | const |
Example:
Output:
Here's the matrix m: 1e+04 0 1 0 1e+04 0 0 0 1e+04 m.isDiagonal() returns: 0 m.isDiagonal(1e-3) returns: 1
bool Eigen::MatrixBase< Derived >::isIdentity | ( | const RealScalar & | prec = NumTraits<Scalar>::dummy_precision() | ) | const |
Example:
Output:
Here's the matrix m: 1 0 0.0001 0 1 0 0 0 1 m.isIdentity() returns: 0 m.isIdentity(1e-3) returns: 1
bool Eigen::MatrixBase< Derived >::isLowerTriangular | ( | const RealScalar & | prec = NumTraits<Scalar>::dummy_precision() | ) | const |
bool Eigen::MatrixBase< Derived >::isOrthogonal | ( | const MatrixBase< OtherDerived > & | other, |
const RealScalar & | prec = NumTraits<Scalar>::dummy_precision() ) const |
Example:
Output:
Here's the vector v: 1 0 0 Here's the vector w: 0.0001 0 1 v.isOrthogonal(w) returns: 0 v.isOrthogonal(w,1e-3) returns: 1
bool Eigen::MatrixBase< Derived >::isSkewSymmetric | ( | const RealScalar & | prec = NumTraits<Scalar>::dummy_precision() | ) | const |
bool Eigen::MatrixBase< Derived >::isUnitary | ( | const RealScalar & | prec = NumTraits<Scalar>::dummy_precision() | ) | const |
m.isUnitary()
returns true if and only if the columns (equivalently, the rows) of m form an orthonormal basis.Example:
Output:
Here's the matrix m: 1 0 0.0001 0 1 0 0 0 1 m.isUnitary() returns: 0 m.isUnitary(1e-3) returns: 1
bool Eigen::MatrixBase< Derived >::isUpperTriangular | ( | const RealScalar & | prec = NumTraits<Scalar>::dummy_precision() | ) | const |
JacobiSVD< typename MatrixBase< Derived >::PlainObject, Options > Eigen::MatrixBase< Derived >::jacobiSvd | ( | ) | const |
This is defined in the SVD module.
*this
computed by two-sided Jacobi transformations.
|
inline |
*this
and other without implicit evaluation.The returned product will behave like any other expressions: the coefficients of the product will be computed once at a time as requested. This might be useful in some extremely rare cases when only a small and no coherent fraction of the result's coefficients have to be computed.
|
inline |
This is defined in the Cholesky module.
*this
|
inline |
This is defined in the Cholesky module.
*this
const MatrixLogarithmReturnValue< Derived > Eigen::MatrixBase< Derived >::log | ( | ) | const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise logarithm use ArrayBase::log .
*this
. MatrixBase< Derived >::RealScalar Eigen::MatrixBase< Derived >::lpNorm | ( | ) | const |
*this
, that is, returns the p-th root of the sum of the p-th powers of the absolute values of the coefficients of *this
. If p is the special value Eigen::Infinity, this function returns the *this
.In all cases, if *this
is empty, then the value 0 is returned.
*this
is a matrix, then its coefficients are interpreted as a 1D vector. Nonetheless, you can easily compute the 1-norm and
|
inline |
This is defined in the LU module.
Synonym of partialPivLu().
*this
.void Eigen::MatrixBase< Derived >::makeHouseholder | ( | EssentialPart & | essential, |
Scalar & | tau, | ||
RealScalar & | beta ) const |
Computes the elementary reflector H such that:
On output:
essential | the essential part of the vector v |
tau | the scaling factor of the Householder transformation |
beta | the result of H * *this |
void Eigen::MatrixBase< Derived >::makeHouseholderInPlace | ( | Scalar & | tau, |
RealScalar & | beta ) |
Computes the elementary reflector H such that:
The essential part of the vector v
is stored in *this.
On output:
tau | the scaling factor of the Householder transformation |
beta | the result of H * *this |
NoAlias< Derived, MatrixBase > Eigen::MatrixBase< Derived >::noalias | ( | ) |
*this
with an operator= assuming no aliasing between *this
and the source expression.More precisely, noalias() allows to bypass the EvalBeforeAssignBit flag. Currently, even though several expressions may alias, only product expressions have this flag. Therefore, noalias() is only useful when the source expression contains a matrix product.
Here are some examples where noalias is useful:
On the other hand the following example will lead to a wrong result:
because the result matrix A is also an operand of the matrix product. Therefore, there is no alternative than evaluating A * B in a temporary, that is the default behavior when you write:
|
inline |
*this
, and for matrices the Frobenius norm. In both cases, it consists in the square root of the sum of the square of all the matrix entries. For vectors, this is also equals to the square root of the dot product of *this
with itself.
|
inline |
Normalizes the vector, i.e. divides it by its own norm.
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
is left unchanged.
|
inline |
*this
by its own norm.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.
|
inline |
*this
and other are not exactly equal to each other.
|
inline |
*this
and other
|
inline |
*this
scaled by the scalar factor scalar T | is the scalar type of scalar. It must be compatible with the scalar type of the given expression. |
*this
divided by the scalar value scalar T | is the scalar type of scalar. It must be compatible with the scalar type of the given expression. |
*this
and other Example:
Output:
0 0 0
|
inline |
*this
by the diagonal matrix diagonal.
|
inline |
*this
and other.
|
inline |
*this
by the skew symmetric matrix skew.
|
inline |
replaces *this
by *this
* other.
*this
Example:
Output:
At start, A = -0.824 -0.199 0.634 0.924 -0.237 0.334 -0.0532 0.146 -0.779 After A *= B, A = 0.634 -0.824 -0.199 0.334 0.924 -0.237 -0.779 -0.0532 0.146 After applyOnTheRight, A = -0.199 0.634 -0.824 -0.237 0.334 0.924 0.146 -0.779 -0.0532
const CwiseBinaryOp< sum< Scalar >, const Derived, const OtherDerived > Eigen::MatrixBase< Derived >::operator+ | ( | const Eigen::MatrixBase< OtherDerived > & | other | ) | const |
*this
and other
|
inline |
replaces *this
by *this
+ other.
*this
const CwiseBinaryOp< difference< Scalar >, const Derived, const OtherDerived > Eigen::MatrixBase< Derived >::operator- | ( | const Eigen::MatrixBase< OtherDerived > & | other | ) | const |
*this
and other
|
inline |
replaces *this
by *this
- other.
*this
|
inline |
Special case of the template operator=, in order to prevent the compiler from generating a default operator= (issue hit with g++ 4.1)
|
inline |
*this
and other are all exactly equal.
|
inline |
|
inline |
Computes the L2 operator norm.
This is defined in the Eigenvalues module.
This function computes the L2 operator norm of a matrix, which is also known as the spectral norm. The norm of a matrix
where the maximum is over all vectors and the norm on the right is the Euclidean vector norm. The norm equals the largest singular value, which is the square root of the largest eigenvalue of the positive semi-definite matrix
The current implementation uses the eigenvalues of
Example:
Output:
The operator norm of the 3x3 matrix of ones is 3
|
inline |
*this
and other
|
inline |
*this
and other Example:
Output:
1 0 1
|
inline |
This is defined in the LU module.
*this
.const MatrixPowerReturnValue< Derived > Eigen::MatrixBase< Derived >::pow | ( | const RealScalar & | p | ) | const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise power to p
use ArrayBase::pow .
p
of *this
. const MatrixComplexPowerReturnValue< Derived > Eigen::MatrixBase< Derived >::pow | ( | const std::complex< RealScalar > & | p | ) | const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise power to p
use ArrayBase::pow .
p
of *this
. MatrixBase< Derived >::template SelfAdjointViewReturnType< UpLo >::Type Eigen::MatrixBase< Derived >::selfadjointView | ( | ) |
The parameter UpLo can be either Upper
or Lower
Example:
Output:
Here is the matrix m: 1804289383 719885386 -1364114958 -465790871 -1550966999 2044897763 -189735855 -1122281286 1365180540 Here is the symmetric matrix extracted from the upper part of m: 1804289383 719885386 -1364114958 719885386 -1550966999 2044897763 -1364114958 2044897763 1365180540 Here is the symmetric matrix extracted from the lower part of m: 1804289383 -465790871 -189735855 -465790871 -1550966999 -1122281286 -189735855 -1122281286 1365180540
MatrixBase< Derived >::template ConstSelfAdjointViewReturnType< UpLo >::Type Eigen::MatrixBase< Derived >::selfadjointView | ( | ) | const |
This is the const version of MatrixBase::selfadjointView()
|
inline |
Writes the identity expression (not necessarily square) into *this.
Example:
Output:
0 0 0 0 1 0 0 0 0 1 0 0 0 0 1 0
|
inline |
Resizes to the given size, and writes the identity expression (not necessarily square) into *this.
rows | the new number of rows |
cols | the new number of columns |
Example:
Output:
1 0 0 0 1 0 0 0 1
|
inline |
Set the coefficients of *this
to the i-th unit (basis) vector.
i | index of the unique coefficient to be set to 1 |
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.
|
inline |
Resizes to the given newSize, and writes the i-th unit (basis) vector into *this.
newSize | the new size of the vector |
i | index of the unique coefficient to be set to 1 |
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.
const MatrixFunctionReturnValue< Derived > Eigen::MatrixBase< Derived >::sin | ( | ) | const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise sine use ArrayBase::sin .
*this
. const MatrixFunctionReturnValue< Derived > Eigen::MatrixBase< Derived >::sinh | ( | ) | const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise hyperbolic sine use ArrayBase::sinh .
*this
. const MatrixSquareRootReturnValue< Derived > Eigen::MatrixBase< Derived >::sqrt | ( | ) | const |
This function requires the <a * href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module. To compute the * coefficient-wise square root use ArrayBase::sqrt .
*this
.
|
inline |
|
inline |
*this
avoiding underflow and overflow. This version use a blockwise two passes algorithm: 1 - find the absolute largest coefficient s
2 - compute For architecture/scalar types supporting vectorization, this version is faster than blueNorm(). Otherwise the blueNorm() is much faster.
|
inline |
Normalizes the vector while avoid underflow and overflow
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 method is analogue to the normalize() method, but it reduces the risk of underflow and overflow when computing the norm.
*this
is left unchanged.
|
inline |
*this
by its own norm while avoiding underflow and overflow.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 method is analogue to the normalized() method, but it reduces the risk of underflow and overflow when computing the norm.
|
inline |
*this
, i.e. the sum of the coefficients on the main diagonal.*this
can be any matrix, not necessarily square.
MatrixBase< Derived >::template TriangularViewReturnType< Mode >::Type Eigen::MatrixBase< Derived >::triangularView | ( | ) |
The parameter Mode can have the following values: Upper
, StrictlyUpper
, UnitUpper
, Lower
, StrictlyLower
, UnitLower
.
Example:
Output:
Here is the matrix m: 1804289383 719885386 -1364114958 -465790871 -1550966999 2044897763 -189735855 -1122281286 1365180540 Here is the upper-triangular matrix extracted from m: 1804289383 719885386 -1364114958 0 -1550966999 2044897763 0 0 1365180540 Here is the strictly-upper-triangular matrix extracted from m: 0 719885386 -1364114958 0 0 2044897763 0 0 0 Here is the unit-lower-triangular matrix extracted from m: 1 0 0 -465790871 1 0 -189735855 -1122281286 1
MatrixBase< Derived >::template ConstTriangularViewReturnType< Mode >::Type Eigen::MatrixBase< Derived >::triangularView | ( | ) | const |
This is the const version of MatrixBase::triangularView()
|
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.
This variant is for fixed-size vector only.
|
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.
|
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.
|
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.
|
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.
|
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.