Eigen::MatrixBase< Derived > Class Template Reference

Base class for all dense matrices, vectors, and expressions. More...

+ Inheritance diagram for Eigen::MatrixBase< Derived >:

Classes

struct  ConstSelfAdjointViewReturnType
 
struct  ConstTriangularViewReturnType
 
struct  SelfAdjointViewReturnType
 
struct  TriangularViewReturnType
 

Public Types

enum  { HomogeneousReturnTypeDirection }
 
enum  { SizeMinusOne }
 
typedef Diagonal< const Derived > ConstDiagonalReturnType
 
typedef Block< const Derived, internal::traits< Derived >::ColsAtCompileTime==1 ? SizeMinusOne :1, internal::traits< Derived >::ColsAtCompileTime==1 ? 1 :SizeMinusOneConstStartMinusOne
 
typedef Diagonal< Derived > DiagonalReturnType
 
typedef Homogeneous< Derived, HomogeneousReturnTypeDirectionHomogeneousReturnType
 
typedef Base::PlainObject PlainObject
 
typedef internal::stem_function< Scalar >::type StemFunction
 
- Public Types inherited from Eigen::DenseBase< Derived >
enum  {
  RowsAtCompileTime ,
  ColsAtCompileTime ,
  SizeAtCompileTime ,
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime ,
  MaxSizeAtCompileTime ,
  IsVectorAtCompileTime ,
  NumDimensions ,
  Flags ,
  IsRowMajor ,
  InnerSizeAtCompileTime ,
  InnerStrideAtCompileTime ,
  OuterStrideAtCompileTime
}
 
enum  { IsPlainObjectBase }
 
typedef DenseCoeffsBase< Derived, internal::accessors_level< Derived >::valueBase
 
typedef Base::CoeffReturnType CoeffReturnType
 
typedef VectorwiseOp< Derived, VerticalColwiseReturnType
 
typedef random_access_iterator_type const_iterator
 
typedef const VectorwiseOp< const Derived, VerticalConstColwiseReturnType
 
typedef const Reverse< const Derived, BothDirectionsConstReverseReturnType
 
typedef const VectorwiseOp< const Derived, HorizontalConstRowwiseReturnType
 
typedef Transpose< const Derived > ConstTransposeReturnType
 
typedef internal::add_const_on_value_type_t< typename internal::eval< Derived >::type > EvalReturnType
 
typedef random_access_iterator_type iterator
 
typedef internal::find_best_packet< Scalar, SizeAtCompileTime >::type PacketScalar
 
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. More...
 
typedef CwiseNullaryOp< internal::scalar_random_op< Scalar >, PlainObjectRandomReturnType
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef Reverse< Derived, BothDirectionsReverseReturnType
 
typedef VectorwiseOp< Derived, HorizontalRowwiseReturnType
 
typedef internal::traits< Derived >::Scalar Scalar
 
typedef internal::traits< Derived >::StorageIndex StorageIndex
 The type used to store indices. More...
 
typedef internal::traits< Derived >::StorageKind StorageKind
 
typedef Transpose< Derived > TransposeReturnType
 
typedef Scalar value_type
 
- Public Types inherited from Eigen::DenseCoeffsBase< Derived, DirectWriteAccessors >
typedef DenseCoeffsBase< Derived, WriteAccessorsBase
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef internal::traits< Derived >::Scalar Scalar
 
- Public Types inherited from Eigen::DenseCoeffsBase< Derived, WriteAccessors >
typedef DenseCoeffsBase< Derived, ReadOnlyAccessorsBase
 
typedef internal::packet_traits< Scalar >::type PacketScalar
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef internal::traits< Derived >::Scalar Scalar
 
typedef internal::traits< Derived >::StorageKind StorageKind
 
- Public Types inherited from Eigen::DenseCoeffsBase< Derived, ReadOnlyAccessors >
typedef EigenBase< Derived > Base
 
typedef std::conditional_t< bool(internal::traits< Derived >::Flags &LvalueBit), const Scalar &, std::conditional_t< internal::is_arithmetic< Scalar >::value, Scalar, const Scalar > > CoeffReturnType
 
typedef internal::add_const_on_value_type_if_arithmetic< typename internal::packet_traits< Scalar >::type >::type PacketReturnType
 
typedef internal::packet_traits< Scalar >::type PacketScalar
 
typedef internal::traits< Derived >::Scalar Scalar
 
typedef internal::traits< Derived >::StorageKind StorageKind
 
- Public Types inherited from Eigen::EigenBase< Derived >
typedef Eigen::Index Index
 The interface type of indices. More...
 
typedef internal::traits< Derived >::StorageKind StorageKind
 

Public Member Functions

const MatrixFunctionReturnValue< Derived > acosh () const
 This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise inverse hyperbolic cosine use ArrayBase::acosh . More...
 
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 unsupported MatrixFunctions module. To compute the coefficient-wise inverse hyperbolic sine use ArrayBase::asinh . More...
 
const PermutationWrapper< const Derived > asPermutation () const
 
const SkewSymmetricWrapper< const Derived > asSkewSymmetric () const
 
const MatrixFunctionReturnValue< Derived > atanh () const
 This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise inverse hyperbolic cosine use ArrayBase::atanh . More...
 
template<int Options = 0>
BDCSVD< PlainObject, Options > bdcSvd () const
 
template<int Options>
BDCSVD< typename MatrixBase< Derived >::PlainObject, Options > bdcSvd () const
 
template<int Options = 0>
EIGEN_DEPRECATED BDCSVD< PlainObject, Options > bdcSvd (unsigned int computationOptions) const
 
template<int Options>
BDCSVD< typename MatrixBase< Derived >::PlainObject, Options > bdcSvd (unsigned int computationOptions) const
 
RealScalar blueNorm () const
 
Matrix< Scalar, 3, 1 > canonicalEulerAngles (Index a0, Index a1, Index a2) const
 
template<typename PermutationIndex = DefaultPermutationIndex>
const ColPivHouseholderQR< PlainObject, PermutationIndex > colPivHouseholderQr () const
 
template<typename PermutationIndexType >
const ColPivHouseholderQR< typename MatrixBase< Derived >::PlainObject, PermutationIndexType > colPivHouseholderQr () const
 
template<typename PermutationIndex = DefaultPermutationIndex>
const CompleteOrthogonalDecomposition< PlainObject, PermutationIndex > completeOrthogonalDecomposition () 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 unsupported MatrixFunctions module. To compute the coefficient-wise cosine use ArrayBase::cos . More...
 
const MatrixFunctionReturnValue< Derived > cosh () const
 This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise hyperbolic cosine use ArrayBase::cosh . More...
 
template<typename OtherDerived >
internal::cross_impl< Derived, OtherDerived >::return_type cross (const MatrixBase< OtherDerived > &other) const
 
template<typename OtherDerived >
std::conditional_t< SizeAtCompileTime==2, Scalar, PlainObjectcross (const MatrixBase< OtherDerived > &other) const
 
template<typename OtherDerived >
PlainObject cross3 (const MatrixBase< OtherDerived > &other) const
 
template<typename OtherDerived >
const SparseMatrixBase< OtherDerived >::template CwiseProductDenseReturnType< Derived >::Type cwiseProduct (const SparseMatrixBase< OtherDerived > &other) const
 
Scalar determinant () const
 
DiagonalReturnType diagonal ()
 
template<int Index>
Diagonal< Derived, Indexdiagonal ()
 
const ConstDiagonalReturnType diagonal () const
 
template<int Index>
const Diagonal< const Derived, Indexdiagonal () const
 
Diagonal< Derived, DynamicIndexdiagonal (Index index)
 
const Diagonal< const Derived, DynamicIndexdiagonal (Index index) const
 
Index diagonalSize () const
 
template<typename OtherDerived >
ScalarBinaryOpTraits< typename internal::traits< Derived >::Scalar, typename internal::traits< OtherDerived >::Scalar >::ReturnType dot (const MatrixBase< OtherDerived > &other) const
 
typedef EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE (ConstStartMinusOne, Scalar, quotient) HNormalizedReturnType
 
EigenvaluesReturnType eigenvalues () const
 Computes the eigenvalues of a matrix. More...
 
EIGEN_DEPRECATED Matrix< Scalar, 3, 1 > eulerAngles (Index a0, Index a1, Index a2) const
 
const MatrixExponentialReturnValue< Derived > exp () const
 This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise exponential use ArrayBase::exp . More...
 
Derived & forceAlignedAccess ()
 
const Derived & forceAlignedAccess () const
 
template<bool Enable>
std::conditional_t< Enable, ForceAlignedAccess< Derived >, Derived & > forceAlignedAccessIf ()
 
template<bool Enable>
Derived & forceAlignedAccessIf ()
 
template<bool Enable>
add_const_on_value_type_t< std::conditional_t< Enable, ForceAlignedAccess< Derived >, Derived & > > forceAlignedAccessIf () const
 
template<bool Enable>
const Derived & forceAlignedAccessIf () const
 
template<typename PermutationIndex = DefaultPermutationIndex>
const FullPivHouseholderQR< PlainObject, PermutationIndex > fullPivHouseholderQr () const
 
template<typename PermutationIndex >
const FullPivHouseholderQR< typename MatrixBase< Derived >::PlainObject, PermutationIndex > fullPivHouseholderQr () const
 
template<typename PermutationIndex = DefaultPermutationIndex>
const FullPivLU< PlainObject, PermutationIndex > fullPivLu () const
 
template<typename PermutationIndex >
const FullPivLU< typename MatrixBase< Derived >::PlainObject, PermutationIndex > fullPivLu () const
 
const HNormalizedReturnType hnormalized () const
 homogeneous normalization More...
 
HomogeneousReturnType homogeneous () const
 
const HouseholderQR< PlainObjecthouseholderQr () 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 = 0>
JacobiSVD< PlainObject, Options > jacobiSvd () const
 
template<int Options>
JacobiSVD< typename MatrixBase< Derived >::PlainObject, Options > jacobiSvd () const
 
template<int Options = 0>
EIGEN_DEPRECATED JacobiSVD< PlainObject, Options > jacobiSvd (unsigned int computationOptions) const
 
template<int Options>
JacobiSVD< typename MatrixBase< Derived >::PlainObject, Options > jacobiSvd (unsigned int computationOptions) const
 
template<typename OtherDerived >
const Product< Derived, OtherDerived, LazyProductlazyProduct (const MatrixBase< OtherDerived > &other) const
 
const LDLT< PlainObjectldlt () const
 
const LLT< PlainObjectllt () const
 
const MatrixLogarithmReturnValue< Derived > log () const
 This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise logarithm use ArrayBase::log . More...
 
template<int p>
RealScalar lpNorm () const
 
template<typename PermutationIndex = DefaultPermutationIndex>
const PartialPivLU< PlainObject, PermutationIndex > lu () 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)
 
MatrixBase< Derived > & matrix ()
 
const MatrixBase< Derived > & matrix () const
 
const MatrixFunctionReturnValue< Derived > matrixFunction (StemFunction f) const
 Helper function for the unsupported MatrixFunctions module. More...
 
NoAlias< Derived, Eigen::MatrixBasenoalias ()
 
RealScalar norm () const
 
void normalize ()
 
const PlainObject normalized () const
 
template<typename OtherDerived >
bool operator!= (const MatrixBase< OtherDerived > &other) const
 
template<typename DiagonalDerived >
const Product< Derived, DiagonalDerived, LazyProductoperator* (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, LazyProductoperator* (const SkewSymmetricBase< SkewDerived > &skew) const
 
template<typename OtherDerived >
Derived & operator*= (const EigenBase< OtherDerived > &other)
 
template<typename OtherDerived >
Derived & operator+= (const MatrixBase< OtherDerived > &other)
 
template<typename OtherDerived >
Derived & operator-= (const MatrixBase< OtherDerived > &other)
 
template<typename OtherDerived >
Derived & operator= (const DenseBase< OtherDerived > &other)
 
template<typename OtherDerived >
Derived & operator= (const EigenBase< OtherDerived > &other)
 
Derived & operator= (const MatrixBase &other)
 
template<typename OtherDerived >
Derived & operator= (const ReturnByValue< OtherDerived > &other)
 
template<typename OtherDerived >
bool operator== (const MatrixBase< OtherDerived > &other) const
 
RealScalar operatorNorm () const
 Computes the L2 operator norm. More...
 
template<typename PermutationIndex = DefaultPermutationIndex>
const PartialPivLU< PlainObject, PermutationIndex > partialPivLu () 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 unsupported MatrixFunctions module. To compute the coefficient-wise power to p use ArrayBase::pow . More...
 
const MatrixComplexPowerReturnValue< Derived > pow (const std::complex< RealScalar > &p) const
 This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise power to p use ArrayBase::pow . More...
 
template<unsigned int UpLo>
SelfAdjointViewReturnType< UpLo >::Type selfadjointView ()
 
template<unsigned int UpLo>
MatrixBase< Derived >::template SelfAdjointViewReturnType< UpLo >::Type selfadjointView ()
 
template<unsigned int UpLo>
ConstSelfAdjointViewReturnType< UpLo >::Type selfadjointView () const
 
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. More...
 
Derived & setUnit (Index i)
 Set the coefficients of *this to the i-th unit (basis) vector. More...
 
Derived & setUnit (Index newSize, Index i)
 Resizes to the given newSize, and writes the i-th unit (basis) vector into *this. More...
 
const MatrixFunctionReturnValue< Derived > sin () const
 This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise sine use ArrayBase::sin . More...
 
const MatrixFunctionReturnValue< Derived > sinh () const
 This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise hyperbolic sine use ArrayBase::sinh . More...
 
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 unsupported MatrixFunctions module. To compute the coefficient-wise square root use ArrayBase::sqrt . More...
 
RealScalar squaredNorm () const
 
RealScalar stableNorm () const
 
void stableNormalize ()
 
const PlainObject stableNormalized () const
 
Scalar trace () const
 
template<unsigned int Mode>
TriangularViewReturnType< Mode >::Type triangularView ()
 
template<unsigned int Mode>
MatrixBase< Derived >::template TriangularViewReturnType< Mode >::Type triangularView ()
 
template<unsigned int Mode>
ConstTriangularViewReturnType< Mode >::Type triangularView () const
 
template<unsigned int Mode>
MatrixBase< Derived >::template ConstTriangularViewReturnType< Mode >::Type triangularView () const
 
PlainObject unitOrthogonal (void) const
 
- Public Member Functions inherited from Eigen::DenseBase< Derived >
bool all () const
 
bool allFinite () const
 
bool any () const
 
iterator begin ()
 
const_iterator begin () const
 
const_iterator cbegin () const
 
const_iterator cend () const
 
ColwiseReturnType colwise ()
 
ConstColwiseReturnType colwise () const
 
Index count () const
 
iterator end ()
 
const_iterator end () const
 
EvalReturnType eval () const
 
template<typename Dest >
void evalTo (Dest &) const
 
void fill (const Scalar &value)
 
template<unsigned int Added, unsigned int Removed>
EIGEN_DEPRECATED const Derived & flagged () const
 
ForceAlignedAccess< Derived > forceAlignedAccess ()
 
const ForceAlignedAccess< Derived > forceAlignedAccess () const
 
template<bool Enable>
std::conditional_t< Enable, ForceAlignedAccess< Derived >, Derived & > forceAlignedAccessIf ()
 
template<bool Enable>
const std::conditional_t< Enable, ForceAlignedAccess< Derived >, Derived & > forceAlignedAccessIf () const
 
const WithFormat< Derived > format (const IOFormat &fmt) const
 
bool hasNaN () const
 
EIGEN_CONSTEXPR Index innerSize () 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
 
bool isMuchSmallerThan (const RealScalar &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 >
Derived & lazyAssign (const DenseBase< OtherDerived > &other)
 
template<typename OtherDerived >
EIGEN_DEPRECATED Derived & lazyAssign (const DenseBase< OtherDerived > &other)
 
template<int p>
RealScalar lpNorm () const
 
template<int NaNPropagation>
internal::traits< Derived >::Scalar maxCoeff () const
 
internal::traits< Derived >::Scalar maxCoeff () const
 
template<int NaNPropagation, typename IndexType >
internal::traits< Derived >::Scalar maxCoeff (IndexType *index) const
 
template<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
 
template<typename IndexType >
internal::traits< Derived >::Scalar maxCoeff (IndexType *row, IndexType *col) const
 
Scalar mean () const
 
template<int NaNPropagation>
internal::traits< Derived >::Scalar minCoeff () const
 
internal::traits< Derived >::Scalar minCoeff () const
 
template<int NaNPropagation, typename IndexType >
internal::traits< Derived >::Scalar minCoeff (IndexType *index) const
 
template<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
 
template<typename IndexType >
internal::traits< Derived >::Scalar minCoeff (IndexType *row, IndexType *col) const
 
const NestByValue< Derived > nestByValue () const
 
Derived & operator*= (const Scalar &other)
 
template<typename OtherDerived >
Derived & operator+= (const EigenBase< OtherDerived > &other)
 
template<typename OtherDerived >
Derived & operator-= (const EigenBase< OtherDerived > &other)
 
Derived & operator/= (const Scalar &other)
 
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. More...
 
template<typename OtherDerived >
Derived & operator= (const ReturnByValue< OtherDerived > &func)
 
EIGEN_CONSTEXPR Index outerSize () const
 
Scalar prod () const
 
template<typename BinaryOp >
Scalar redux (const BinaryOp &func) 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
 
void resize (Index newSize)
 
void resize (Index rows, Index cols)
 
ReverseReturnType reverse ()
 
ConstReverseReturnType reverse () const
 
void reverseInPlace ()
 
RowwiseReturnType rowwise ()
 
ConstRowwiseReturnType rowwise () const
 
template<typename ThenDerived , typename ElseDerived >
CwiseTernaryOp< internal::scalar_boolean_select_op< typename DenseBase< ThenDerived >::Scalar, typename DenseBase< ElseDerived >::Scalar, Scalar >, ThenDerived, ElseDerived, Derived > select (const DenseBase< ThenDerived > &thenMatrix, const DenseBase< ElseDerived > &elseMatrix) const
 
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, Scalar >, ThenDerived, typename DenseBase< ThenDerived >::ConstantReturnType, Derived > select (const DenseBase< ThenDerived > &thenMatrix, const typename DenseBase< ThenDerived >::Scalar &elseScalar) 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, Scalar >, typename DenseBase< ElseDerived >::ConstantReturnType, ElseDerived, Derived > select (const typename DenseBase< ElseDerived >::Scalar &thenScalar, const DenseBase< ElseDerived > &elseMatrix) 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 & setEqualSpaced (const Scalar &low, const Scalar &step)
 
Derived & setEqualSpaced (Index size, const Scalar &low, const Scalar &step)
 
Derived & setLinSpaced (const Scalar &low, const Scalar &high)
 Sets a linearly spaced vector. More...
 
Derived & setLinSpaced (Index size, const Scalar &low, const Scalar &high)
 Sets a linearly spaced vector. More...
 
Derived & setOnes ()
 
Derived & setRandom ()
 
Derived & setZero ()
 
Scalar sum () const
 
template<typename OtherDerived >
void swap (const DenseBase< OtherDerived > &other)
 
template<typename OtherDerived >
void swap (PlainObjectBase< OtherDerived > &other)
 
Scalar trace () const
 
TransposeReturnType transpose ()
 
const ConstTransposeReturnType transpose () const
 
void transposeInPlace ()
 
CoeffReturnType value () const
 
template<typename Visitor >
void visit (Visitor &func) const
 
- Public Member Functions inherited from Eigen::DenseCoeffsBase< Derived, DirectWriteAccessors >
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index colStride () const EIGEN_NOEXCEPT
 
Derived & derived ()
 
const Derived & derived () const
 
EIGEN_CONSTEXPR Index innerStride () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index outerStride () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index rowStride () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index stride () const EIGEN_NOEXCEPT
 
- Public Member Functions inherited from Eigen::DenseCoeffsBase< Derived, WriteAccessors >
ScalarcoeffRef (Index index)
 
ScalarcoeffRef (Index row, Index col)
 
ScalarcoeffRefByOuterInner (Index outer, Index inner)
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & derived ()
 
const Derived & derived () const
 
Scalaroperator() (Index index)
 
Scalaroperator() (Index row, Index col)
 
Scalaroperator[] (Index index)
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
Scalarw ()
 
Scalarx ()
 
Scalary ()
 
Scalarz ()
 
- Public Member Functions inherited from Eigen::DenseCoeffsBase< Derived, ReadOnlyAccessors >
CoeffReturnType coeff (Index index) const
 
CoeffReturnType coeff (Index row, Index col) const
 
CoeffReturnType coeffByOuterInner (Index outer, Index inner) const
 
Index colIndexByOuterInner (Index outer, Index inner) const
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & derived ()
 
const Derived & derived () const
 
CoeffReturnType operator() (Index index) const
 
CoeffReturnType operator() (Index row, Index col) const
 
CoeffReturnType operator[] (Index index) const
 
template<int LoadMode>
PacketReturnType packet (Index index) const
 
template<int LoadMode>
PacketReturnType packet (Index row, Index col) const
 
template<int LoadMode>
PacketReturnType packetByOuterInner (Index outer, Index inner) const
 
Index rowIndexByOuterInner (Index outer, Index inner) const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
CoeffReturnType w () const
 
CoeffReturnType x () const
 
CoeffReturnType y () const
 
CoeffReturnType z () const
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
template<typename Dest >
void addTo (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheLeft (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheRight (Dest &dst) const
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & const_cast_derived () const
 
const Derived & const_derived () const
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Dest >
void evalTo (Dest &dst) const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
template<typename Dest >
void subTo (Dest &dst) const
 

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 Public Member Functions inherited from Eigen::DenseBase< Derived >
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 RandomAccessEqualSpacedReturnType EqualSpaced (const Scalar &low, const Scalar &step)
 
static const RandomAccessEqualSpacedReturnType EqualSpaced (Index size, const Scalar &low, const Scalar &step)
 
static const RandomAccessLinSpacedReturnType LinSpaced (const Scalar &low, const Scalar &high)
 Sets a linearly spaced vector. More...
 
static const RandomAccessLinSpacedReturnType LinSpaced (Index size, const Scalar &low, const Scalar &high)
 Sets a linearly spaced vector. More...
 
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 ConstantReturnType Zero ()
 
static const ConstantReturnType Zero (Index rows, Index cols)
 
static const ConstantReturnType Zero (Index size)
 

Protected Member Functions

template<typename OtherDerived >
Derived & operator+= (const ArrayBase< OtherDerived > &)
 
template<typename OtherDerived >
Derived & operator-= (const ArrayBase< OtherDerived > &)
 
- Protected Member Functions inherited from Eigen::DenseBase< Derived >
constexpr DenseBase ()
 
- Protected Member Functions inherited from Eigen::DenseCoeffsBase< Derived, ReadOnlyAccessors >
void coeffRef ()
 
void coeffRefByOuterInner ()
 
void colStride ()
 
void copyCoeff ()
 
void copyCoeffByOuterInner ()
 
void copyPacket ()
 
void copyPacketByOuterInner ()
 
void innerStride ()
 
void outerStride ()
 
void rowStride ()
 
void stride ()
 
void writePacket ()
 
void writePacketByOuterInner ()
 

Private Member Functions

template<typename OtherDerived >
 MatrixBase (const MatrixBase< OtherDerived > &)
 
 MatrixBase (int)
 
 MatrixBase (int, int)
 

Additional Inherited Members

Detailed Description

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

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.

Template Parameters
Derivedis 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.

template<typename Derived>
void printFirstRow(const Eigen::MatrixBase<Derived>& x)
{
cout << x.row(0) << endl;
}
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52

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.

See also
The class hierarchy

Definition at line 50 of file MatrixBase.h.

Member Typedef Documentation

◆ ConstDiagonalReturnType

template<typename Derived >
typedef Diagonal<const Derived> Eigen::MatrixBase< Derived >::ConstDiagonalReturnType

Definition at line 216 of file MatrixBase.h.

◆ ConstStartMinusOne

template<typename Derived >
typedef Block<const Derived, internal::traits<Derived>::ColsAtCompileTime==1 ? SizeMinusOne : 1, internal::traits<Derived>::ColsAtCompileTime==1 ? 1 : SizeMinusOne> Eigen::MatrixBase< Derived >::ConstStartMinusOne

Definition at line 420 of file MatrixBase.h.

◆ DiagonalReturnType

template<typename Derived >
typedef Diagonal<Derived> Eigen::MatrixBase< Derived >::DiagonalReturnType

Definition at line 212 of file MatrixBase.h.

◆ HomogeneousReturnType

template<typename Derived >
typedef Homogeneous<Derived, HomogeneousReturnTypeDirection> Eigen::MatrixBase< Derived >::HomogeneousReturnType

Definition at line 411 of file MatrixBase.h.

◆ PlainObject

template<typename Derived >
typedef Base::PlainObject Eigen::MatrixBase< Derived >::PlainObject

Definition at line 106 of file MatrixBase.h.

◆ StemFunction

template<typename Derived >
typedef internal::stem_function<Scalar>::type Eigen::MatrixBase< Derived >::StemFunction

Definition at line 464 of file MatrixBase.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename Derived >
anonymous enum
Enumerator
HomogeneousReturnTypeDirection 

Definition at line 409 of file MatrixBase.h.

@ HomogeneousReturnTypeDirection
Definition: MatrixBase.h:409
@ Horizontal
Definition: Constants.h:269
@ Vertical
Definition: Constants.h:266
const unsigned int RowMajorBit
Definition: Constants.h:68

◆ anonymous enum

template<typename Derived >
anonymous enum
Enumerator
SizeMinusOne 

Definition at line 415 of file MatrixBase.h.

415  {
417  };
const int Dynamic
Definition: Constants.h:24

Constructor & Destructor Documentation

◆ MatrixBase() [1/3]

template<typename Derived >
Eigen::MatrixBase< Derived >::MatrixBase ( int  )
explicitprivate

◆ MatrixBase() [2/3]

template<typename Derived >
Eigen::MatrixBase< Derived >::MatrixBase ( int  ,
int   
)
private

◆ MatrixBase() [3/3]

template<typename Derived >
template<typename OtherDerived >
Eigen::MatrixBase< Derived >::MatrixBase ( const MatrixBase< OtherDerived > &  )
explicitprivate

Member Function Documentation

◆ acosh()

template<typename Derived >
const MatrixFunctionReturnValue<Derived> Eigen::MatrixBase< Derived >::acosh ( ) const

This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise inverse hyperbolic cosine use ArrayBase::acosh .

Returns
an expression of the matrix inverse hyperbolic cosine of *this.

◆ adjoint()

template<typename Derived >
const MatrixBase< Derived >::AdjointReturnType Eigen::MatrixBase< Derived >::adjoint
inline
Returns
an expression of the adjoint (i.e. conjugate transpose) of *this.

Example:

cout << "Here is the 2x2 complex matrix m:" << endl << m << endl;
cout << "Here is the adjoint of m:" << endl << m.adjoint() << endl;
Matrix3f m
static const RandomReturnType Random()
Definition: Random.h:114
Matrix< std::complex< float >, 2, 2 > Matrix2cf
2×2 matrix of type std::complex<float>.
Definition: Matrix.h:503

Output:

Here is the 2x2 complex matrix m:
 (-0.211,0.68) (-0.605,0.823)
 (0.597,0.566)  (0.536,-0.33)
Here is the adjoint of m:
 (-0.211,-0.68)  (0.597,-0.566)
(-0.605,-0.823)    (0.536,0.33)
Warning
If you want to replace a matrix by its own adjoint, do NOT do this:
m = m.adjoint(); // bug!!! caused by aliasing effect
Instead, use the adjointInPlace() method:
m.adjointInPlace();
which gives Eigen good opportunities for optimization, or alternatively you can also do:
m = m.adjoint().eval();
See also
adjointInPlace(), transpose(), conjugate(), class Transpose, class internal::scalar_conjugate_op

Definition at line 223 of file Transpose.h.

224 {
225  return AdjointReturnType(this->transpose());
226 }
TransposeReturnType transpose()
Definition: Transpose.h:184

◆ adjointInPlace()

template<typename Derived >
void Eigen::MatrixBase< Derived >::adjointInPlace
inline

"in place" adjoint implementation This is the "in place" version of adjoint(): it replaces *this by its own transpose. Thus, doing

m.adjointInPlace();

has the same effect on m as doing

m = m.adjoint().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 adjoint. If you just need the adjoint of a matrix, use adjoint().

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(), transposeInPlace()

Definition at line 377 of file Transpose.h.

378 {
379  derived() = adjoint().eval();
380 }
const AdjointReturnType adjoint() const
Definition: Transpose.h:223

◆ applyHouseholderOnTheLeft()

template<typename Derived >
template<typename EssentialPart >
void Eigen::MatrixBase< Derived >::applyHouseholderOnTheLeft ( const EssentialPart &  essential,
const Scalar tau,
Scalar workspace 
)

Apply the elementary reflector H given by \( H = I - tau v v^*\) with \( v^T = [1 essential^T] \) from the left to a vector or matrix.

On input:

Parameters
essentialthe essential part of the vector v
tauthe scaling factor of the Householder transformation
workspacea pointer to working space with at least this->cols() entries
See also
MatrixBase::makeHouseholder(), MatrixBase::makeHouseholderInPlace(), MatrixBase::applyHouseholderOnTheRight()

Definition at line 118 of file Householder.h.

122 {
123  if(rows() == 1)
124  {
125  *this *= Scalar(1)-tau;
126  }
127  else if(!numext::is_exactly_zero(tau))
128  {
129  Map<typename internal::plain_row_type<PlainObject>::type> tmp(workspace,cols());
130  Block<Derived, EssentialPart::SizeAtCompileTime, Derived::ColsAtCompileTime> bottom(derived(), 1, 0, rows()-1, cols());
131  tmp.noalias() = essential.adjoint() * bottom;
132  tmp += this->row(0);
133  this->row(0) -= tau * tmp;
134  bottom.noalias() -= tau * essential * tmp;
135  }
136 }
RowXpr row(Index i)
This is the const version of row(). *‍/.
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:61
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: EigenBase.h:65
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:62
bool is_exactly_zero(const X &x)
Definition: Meta.h:475

◆ applyHouseholderOnTheRight()

template<typename Derived >
template<typename EssentialPart >
void Eigen::MatrixBase< Derived >::applyHouseholderOnTheRight ( const EssentialPart &  essential,
const Scalar tau,
Scalar workspace 
)

Apply the elementary reflector H given by \( H = I - tau v v^*\) with \( v^T = [1 essential^T] \) from the right to a vector or matrix.

On input:

Parameters
essentialthe essential part of the vector v
tauthe scaling factor of the Householder transformation
workspacea pointer to working space with at least this->rows() entries
See also
MatrixBase::makeHouseholder(), MatrixBase::makeHouseholderInPlace(), MatrixBase::applyHouseholderOnTheLeft()

Definition at line 156 of file Householder.h.

160 {
161  if(cols() == 1)
162  {
163  *this *= Scalar(1)-tau;
164  }
165  else if(!numext::is_exactly_zero(tau))
166  {
167  Map<typename internal::plain_col_type<PlainObject>::type> tmp(workspace,rows());
168  Block<Derived, Derived::RowsAtCompileTime, EssentialPart::SizeAtCompileTime> right(derived(), 0, 1, rows(), cols()-1);
169  tmp.noalias() = right * essential;
170  tmp += this->col(0);
171  this->col(0) -= tau * tmp;
172  right.noalias() -= tau * tmp * essential.adjoint();
173  }
174 }
ColXpr col(Index i)
This is the const version of col().

◆ applyOnTheLeft() [1/2]

template<typename Derived >
template<typename OtherDerived >
void Eigen::MatrixBase< Derived >::applyOnTheLeft ( const EigenBase< OtherDerived > &  other)
inline

replaces *this by other * *this.

Example:

B << 0,1,0,
0,0,1,
1,0,0;
cout << "At start, A = " << endl << A << endl;
A.applyOnTheLeft(B);
cout << "After applyOnTheLeft, A = " << endl << A << endl;
MatrixXcf A
MatrixXf B
Matrix< float, 3, 3 > Matrix3f
3×3 matrix of type float.
Definition: Matrix.h:501

Output:

At start, A = 
  0.68  0.597  -0.33
-0.211  0.823  0.536
 0.566 -0.605 -0.444
After applyOnTheLeft, A = 
-0.211  0.823  0.536
 0.566 -0.605 -0.444
  0.68  0.597  -0.33

Definition at line 544 of file MatrixBase.h.

545 {
546  other.derived().applyThisOnTheLeft(derived());
547 }

◆ applyOnTheLeft() [2/2]

template<typename Derived >
template<typename OtherScalar >
void Eigen::MatrixBase< Derived >::applyOnTheLeft ( Index  p,
Index  q,
const JacobiRotation< OtherScalar > &  j 
)
inline

This is defined in the Jacobi module.

#include <Eigen/Jacobi>

Applies the rotation in the plane j to the rows p and q of *this, i.e., it computes B = J * B, with \( B = \left ( \begin{array}{cc} \text{*this.row}(p) \\ \text{*this.row}(q) \end{array} \right ) \).

See also
class JacobiRotation, MatrixBase::applyOnTheRight(), internal::apply_rotation_in_the_plane()

Definition at line 297 of file Jacobi.h.

298 {
299  RowXpr x(this->row(p));
300  RowXpr y(this->row(q));
302 }
float * p
void apply_rotation_in_the_plane(DenseBase< VectorX > &xpr_x, DenseBase< VectorY > &xpr_y, const JacobiRotation< OtherScalar > &j)
Definition: Jacobi.h:455
std::ptrdiff_t j

◆ applyOnTheRight()

template<typename Derived >
template<typename OtherDerived >
void Eigen::MatrixBase< Derived >::applyOnTheRight ( const EigenBase< OtherDerived > &  other)
inline

replaces *this by *this * other. It is equivalent to MatrixBase::operator*=().

Example:

B << 0,1,0,
0,0,1,
1,0,0;
cout << "At start, A = " << endl << A << endl;
A *= B;
cout << "After A *= B, A = " << endl << A << endl;
A.applyOnTheRight(B); // equivalent to A *= B
cout << "After applyOnTheRight, A = " << endl << A << endl;

Output:

At start, A = 
  0.68  0.597  -0.33
-0.211  0.823  0.536
 0.566 -0.605 -0.444
After A *= B, A = 
 -0.33   0.68  0.597
 0.536 -0.211  0.823
-0.444  0.566 -0.605
After applyOnTheRight, A = 
 0.597  -0.33   0.68
 0.823  0.536 -0.211
-0.605 -0.444  0.566

Definition at line 532 of file MatrixBase.h.

533 {
534  other.derived().applyThisOnTheRight(derived());
535 }

◆ array() [1/2]

template<typename Derived >
ArrayWrapper<Derived> Eigen::MatrixBase< Derived >::array ( )
inline
Returns
an Array expression of this matrix
See also
ArrayBase::matrix()

Definition at line 324 of file MatrixBase.h.

324 { return ArrayWrapper<Derived>(derived()); }

◆ array() [2/2]

template<typename Derived >
const ArrayWrapper<const Derived> Eigen::MatrixBase< Derived >::array ( ) const
inline
Returns
a const Array expression of this matrix
See also
ArrayBase::matrix()

Definition at line 327 of file MatrixBase.h.

327 { return ArrayWrapper<const Derived>(derived()); }

◆ asDiagonal()

template<typename Derived >
const DiagonalWrapper< const Derived > Eigen::MatrixBase< Derived >::asDiagonal
inline
Returns
a pseudo-expression of a diagonal matrix with *this as vector of diagonal coefficients

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 << Matrix3i(Vector3i(2,5,6).asDiagonal()) << endl;
const DiagonalWrapper< const Derived > asDiagonal() const
Matrix< int, 3, 1 > Vector3i
3×1 vector of type int.
Definition: Matrix.h:500
Matrix< int, 3, 3 > Matrix3i
3×3 matrix of type int.
Definition: Matrix.h:500

Output:

2 0 0
0 5 0
0 0 6
See also
class DiagonalWrapper, class DiagonalMatrix, diagonal(), isDiagonal()

Definition at line 384 of file DiagonalMatrix.h.

385 {
386  return DiagonalWrapper<const Derived>(derived());
387 }

◆ asinh()

template<typename Derived >
const MatrixFunctionReturnValue<Derived> Eigen::MatrixBase< Derived >::asinh ( ) const

This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise inverse hyperbolic sine use ArrayBase::asinh .

Returns
an expression of the matrix inverse hyperbolic sine of *this.

◆ asPermutation()

template<typename Derived >
const PermutationWrapper< const Derived > Eigen::MatrixBase< Derived >::asPermutation

Definition at line 594 of file PermutationMatrix.h.

595 {
596  return derived();
597 }

◆ asSkewSymmetric()

template<typename Derived >
const SkewSymmetricWrapper< const Derived > Eigen::MatrixBase< Derived >::asSkewSymmetric
inline
Returns
a pseudo-expression of a skew symmetric matrix with *this as vector of coefficients

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
class SkewSymmetricWrapper, class SkewSymmetricMatrix3, vector(), isSkewSymmetric()

Definition at line 349 of file SkewSymmetricMatrix3.h.

350 {
351  return SkewSymmetricWrapper<const Derived>(derived());
352 }

◆ atanh()

template<typename Derived >
const MatrixFunctionReturnValue<Derived> Eigen::MatrixBase< Derived >::atanh ( ) const

This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise inverse hyperbolic cosine use ArrayBase::atanh .

Returns
an expression of the matrix inverse hyperbolic cosine of *this.

◆ bdcSvd() [1/4]

template<typename Derived >
template<int Options = 0>
BDCSVD<PlainObject, Options> Eigen::MatrixBase< Derived >::bdcSvd ( ) const
inline

◆ bdcSvd() [2/4]

template<typename Derived >
template<int Options>
BDCSVD<typename MatrixBase<Derived>::PlainObject, Options> Eigen::MatrixBase< Derived >::bdcSvd ( ) const

This is defined in the SVD module.

#include <Eigen/SVD>
Returns
the singular value decomposition of *this computed by Divide & Conquer algorithm
See also
class BDCSVD

Definition at line 1475 of file BDCSVD.h.

1475  {
1476  return BDCSVD<PlainObject, Options>(*this);
1477 }

◆ bdcSvd() [3/4]

template<typename Derived >
template<int Options = 0>
EIGEN_DEPRECATED BDCSVD<PlainObject, Options> Eigen::MatrixBase< Derived >::bdcSvd ( unsigned int  computationOptions) const
inline

◆ bdcSvd() [4/4]

template<typename Derived >
template<int Options>
BDCSVD<typename MatrixBase<Derived>::PlainObject, Options> Eigen::MatrixBase< Derived >::bdcSvd ( unsigned int  computationOptions) const

This is defined in the SVD module.

#include <Eigen/SVD>
Returns
the singular value decomposition of *this computed by Divide & Conquer algorithm
See also
class BDCSVD

Definition at line 1487 of file BDCSVD.h.

1488  {
1489  return BDCSVD<PlainObject, Options>(*this, computationOptions);
1490 }

◆ blueNorm()

template<typename Derived >
NumTraits< typename internal::traits< Derived >::Scalar >::Real Eigen::MatrixBase< Derived >::blueNorm
inline
Returns
the l2 norm of *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.

See also
norm(), stableNorm(), hypotNorm()

Definition at line 231 of file StableNorm.h.

232 {
233  return internal::blueNorm_impl(*this);
234 }
NumTraits< typename traits< Derived >::Scalar >::Real blueNorm_impl(const EigenBase< Derived > &_vec)
Definition: StableNorm.h:122

◆ colPivHouseholderQr() [1/2]

template<typename Derived >
template<typename PermutationIndex = DefaultPermutationIndex>
const ColPivHouseholderQR<PlainObject, PermutationIndex> Eigen::MatrixBase< Derived >::colPivHouseholderQr ( ) const
inline

◆ colPivHouseholderQr() [2/2]

template<typename Derived >
template<typename PermutationIndexType >
const ColPivHouseholderQR<typename MatrixBase<Derived>::PlainObject, PermutationIndexType> Eigen::MatrixBase< Derived >::colPivHouseholderQr ( ) const
Returns
the column-pivoting Householder QR decomposition of *this.
See also
class ColPivHouseholderQR

Definition at line 673 of file ColPivHouseholderQR.h.

674 {
675  return ColPivHouseholderQR<PlainObject, PermutationIndexType>(eval());
676 }
EvalReturnType eval() const
Definition: DenseBase.h:405

◆ completeOrthogonalDecomposition() [1/2]

template<typename Derived >
template<typename PermutationIndex = DefaultPermutationIndex>
const CompleteOrthogonalDecomposition<PlainObject, PermutationIndex> Eigen::MatrixBase< Derived >::completeOrthogonalDecomposition ( ) const
inline

◆ completeOrthogonalDecomposition() [2/2]

template<typename Derived >
template<typename PermutationIndex >
const CompleteOrthogonalDecomposition<typename MatrixBase<Derived>::PlainObject, PermutationIndex> Eigen::MatrixBase< Derived >::completeOrthogonalDecomposition ( ) const
Returns
the complete orthogonal decomposition of *this.
See also
class CompleteOrthogonalDecomposition

Definition at line 645 of file CompleteOrthogonalDecomposition.h.

645  {
646  return CompleteOrthogonalDecomposition<PlainObject>(eval());
647 }

◆ computeInverseAndDetWithCheck()

template<typename Derived >
template<typename ResultType >
void Eigen::MatrixBase< Derived >::computeInverseAndDetWithCheck ( ResultType &  inverse,
typename ResultType::Scalar &  determinant,
bool invertible,
const RealScalar absDeterminantThreshold = NumTraits<Scalar>::dummy_precision() 
) const
inline

This is defined in the LU module.

#include <Eigen/LU>

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.

Parameters
inverseReference to the matrix in which to store the inverse.
determinantReference to the variable in which to store the determinant.
invertibleReference to the bool variable in which to store whether the matrix is invertible.
absDeterminantThresholdOptional parameter controlling the invertibility check. The matrix will be declared invertible if the absolute value of its determinant is greater than this threshold.

Example:

cout << "Here is the matrix m:" << endl << m << endl;
bool invertible;
double determinant;
m.computeInverseAndDetWithCheck(inverse,determinant,invertible);
cout << "Its determinant is " << determinant << endl;
if(invertible) {
cout << "It is invertible, and its inverse is:" << endl << inverse << endl;
}
else {
cout << "It is not invertible." << endl;
}
const Inverse< Derived > inverse() const
Definition: InverseImpl.h:350
Scalar determinant() const
Definition: Determinant.h:110
Matrix< double, 3, 3 > Matrix3d
3×3 matrix of type double.
Definition: Matrix.h:502

Output:

Here is the matrix m:
  0.68  0.597  -0.33
-0.211  0.823  0.536
 0.566 -0.605 -0.444
Its determinant is 0.209
It is invertible, and its inverse is:
-0.199   2.23   2.84
  1.01 -0.555  -1.42
 -1.62   3.59   3.29
See also
inverse(), computeInverseWithCheck()

Definition at line 379 of file InverseImpl.h.

385 {
386  // i'd love to put some static assertions there, but SFINAE means that they have no effect...
387  eigen_assert(rows() == cols());
388  // for 2x2, it's worth giving a chance to avoid evaluating.
389  // for larger sizes, evaluating has negligible cost and limits code size.
390  typedef std::conditional_t<
391  RowsAtCompileTime == 2,
392  internal::remove_all_t<typename internal::nested_eval<Derived, 2>::type>,
394  > MatrixType;
395  internal::compute_inverse_and_det_with_check<MatrixType, ResultType>::run
396  (derived(), absDeterminantThreshold, inverse, determinant, invertible);
397 }
#define eigen_assert(x)
Definition: Macros.h:902
Matrix< float, 1, Dynamic > MatrixType
Base::PlainObject PlainObject
Definition: MatrixBase.h:106

◆ computeInverseWithCheck()

template<typename Derived >
template<typename ResultType >
void Eigen::MatrixBase< Derived >::computeInverseWithCheck ( ResultType &  inverse,
bool invertible,
const RealScalar absDeterminantThreshold = NumTraits<Scalar>::dummy_precision() 
) const
inline

This is defined in the LU module.

#include <Eigen/LU>

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.

Parameters
inverseReference to the matrix in which to store the inverse.
invertibleReference to the bool variable in which to store whether the matrix is invertible.
absDeterminantThresholdOptional parameter controlling the invertibility check. The matrix will be declared invertible if the absolute value of its determinant is greater than this threshold.

Example:

cout << "Here is the matrix m:" << endl << m << endl;
bool invertible;
m.computeInverseWithCheck(inverse,invertible);
if(invertible) {
cout << "It is invertible, and its inverse is:" << endl << inverse << endl;
}
else {
cout << "It is not invertible." << endl;
}

Output:

Here is the matrix m:
  0.68  0.597  -0.33
-0.211  0.823  0.536
 0.566 -0.605 -0.444
It is invertible, and its inverse is:
-0.199   2.23   2.84
  1.01 -0.555  -1.42
 -1.62   3.59   3.29
See also
inverse(), computeInverseAndDetWithCheck()

Definition at line 420 of file InverseImpl.h.

425 {
427  // i'd love to put some static assertions there, but SFINAE means that they have no effect...
428  eigen_assert(rows() == cols());
429  computeInverseAndDetWithCheck(inverse,determinant,invertible,absDeterminantThreshold);
430 }
void computeInverseAndDetWithCheck(ResultType &inverse, typename ResultType::Scalar &determinant, bool &invertible, const RealScalar &absDeterminantThreshold=NumTraits< Scalar >::dummy_precision()) const
Definition: InverseImpl.h:379

◆ cos()

template<typename Derived >
const MatrixFunctionReturnValue<Derived> Eigen::MatrixBase< Derived >::cos ( ) const

This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise cosine use ArrayBase::cos .

Returns
an expression of the matrix cosine of *this.

◆ cosh()

template<typename Derived >
const MatrixFunctionReturnValue<Derived> Eigen::MatrixBase< Derived >::cosh ( ) const

This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise hyperbolic cosine use ArrayBase::cosh .

Returns
an expression of the matrix hyperbolic cosine of *this.

◆ cross()

template<typename Derived >
template<typename OtherDerived >
internal::cross_impl<Derived, OtherDerived>::return_type Eigen::MatrixBase< Derived >::cross ( const MatrixBase< OtherDerived > &  other) const
inline

◆ cwiseProduct()

template<typename Derived >
template<typename OtherDerived >
const SparseMatrixBase<OtherDerived>::template CwiseProductDenseReturnType<Derived>::Type Eigen::MatrixBase< Derived >::cwiseProduct ( const SparseMatrixBase< OtherDerived > &  other) const
inline

Definition at line 457 of file MatrixBase.h.

458  {
459  return other.cwiseProduct(derived());
460  }

◆ determinant()

template<typename Derived >
internal::traits< Derived >::Scalar Eigen::MatrixBase< Derived >::determinant
inline

This is defined in the LU module.

#include <Eigen/LU>
Returns
the determinant of this matrix

Definition at line 110 of file Determinant.h.

111 {
112  eigen_assert(rows() == cols());
113  typedef typename internal::nested_eval<Derived,Base::RowsAtCompileTime>::type Nested;
114  return internal::determinant_impl<internal::remove_all_t<Nested>>::run(derived());
115 }

◆ diagonal() [1/6]

template<typename Derived >
Diagonal< Derived, Index_ > Eigen::MatrixBase< Derived >::diagonal
inline
Returns
an expression of the main diagonal of the matrix *this

*this is not required to be square.

Example:

cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here are the coefficients on the main diagonal of m:" << endl
<< m.diagonal() << endl;

Output:

Here is the matrix m:
 7  6 -3
-2  9  6
 6 -6 -5
Here are the coefficients on the main diagonal of m:
 7
 9
-5
See also
class Diagonal
Returns
an expression of the DiagIndex-th sub or super diagonal of the matrix *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:

cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here are the coefficients on the 1st super-diagonal and 2nd sub-diagonal of m:" << endl
<< m.diagonal<1>().transpose() << endl
<< m.diagonal<-2>().transpose() << endl;
Matrix< int, 4, 4 > Matrix4i
4×4 matrix of type int.
Definition: Matrix.h:500

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here are the coefficients on the 1st super-diagonal and 2nd sub-diagonal of m:
9 1 9
6 6
See also
MatrixBase::diagonal(), class Diagonal

Definition at line 189 of file Diagonal.h.

190 {
191  return DiagonalReturnType(derived());
192 }
Diagonal< Derived > DiagonalReturnType
Definition: MatrixBase.h:212

◆ diagonal() [2/6]

template<typename Derived >
template<int Index>
Diagonal<Derived, Index> Eigen::MatrixBase< Derived >::diagonal ( )

◆ diagonal() [3/6]

template<typename Derived >
const Diagonal< const Derived, Index_ > Eigen::MatrixBase< Derived >::diagonal
inline

This is the const version of diagonal().

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

Definition at line 198 of file Diagonal.h.

199 {
201 }
Diagonal< const Derived > ConstDiagonalReturnType
Definition: MatrixBase.h:216

◆ diagonal() [4/6]

template<typename Derived >
template<int Index>
const Diagonal<const Derived, Index> Eigen::MatrixBase< Derived >::diagonal ( ) const

◆ diagonal() [5/6]

template<typename Derived >
Diagonal< Derived, DynamicIndex > Eigen::MatrixBase< Derived >::diagonal ( Index  index)
inline
Returns
an expression of the DiagIndex-th sub or super diagonal of the matrix *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:

cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here are the coefficients on the 1st super-diagonal and 2nd sub-diagonal of m:" << endl
<< m.diagonal(1).transpose() << endl
<< m.diagonal(-2).transpose() << endl;

Output:

Here is the matrix m:
 7  9 -5 -3
-2 -6  1  0
 6 -3  0  9
 6  6  3  9
Here are the coefficients on the 1st super-diagonal and 2nd sub-diagonal of m:
9 1 9
6 6
See also
MatrixBase::diagonal(), class Diagonal

Definition at line 216 of file Diagonal.h.

217 {
218  return Diagonal<Derived, DynamicIndex>(derived(), index);
219 }

◆ diagonal() [6/6]

template<typename Derived >
const Diagonal< const Derived, DynamicIndex > Eigen::MatrixBase< Derived >::diagonal ( Index  index) const
inline

This is the const version of diagonal(Index).

Definition at line 224 of file Diagonal.h.

225 {
226  return Diagonal<const Derived, DynamicIndex>(derived(), index);
227 }

◆ diagonalSize()

template<typename Derived >
Index Eigen::MatrixBase< Derived >::diagonalSize ( ) const
inline
Returns
the size of the main diagonal, which is min(rows(),cols()).
See also
rows(), cols(), SizeAtCompileTime.

Definition at line 104 of file MatrixBase.h.

104 { return (numext::mini)(rows(),cols()); }
EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)

◆ dot()

template<typename Derived >
template<typename OtherDerived >
ScalarBinaryOpTraits< typename internal::traits< Derived >::Scalar, typename internal::traits< OtherDerived >::Scalar >::ReturnType Eigen::MatrixBase< Derived >::dot ( const MatrixBase< OtherDerived > &  other) const
inline
Returns
the dot product of *this with other.

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.

Note
If the scalar type is complex numbers, then this function returns the hermitian (sesquilinear) dot product, conjugate-linear in the first variable and linear in the second variable.
See also
squaredNorm(), norm()

Definition at line 69 of file Dot.h.

70 {
73  EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Derived,OtherDerived)
74 #if !(defined(EIGEN_NO_STATIC_ASSERT) && defined(EIGEN_NO_DEBUG))
76  Eigen::internal::scalar_conj_product_op<Scalar EIGEN_COMMA typename OtherDerived::Scalar>,
77  Scalar, typename OtherDerived::Scalar);
78 #endif
79 
80  eigen_assert(size() == other.size());
81 
82  return internal::dot_nocheck<Derived,OtherDerived>::run(*this, other);
83 }
#define EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(TYPE0, TYPE1)
Definition: StaticAssert.h:61
#define EIGEN_STATIC_ASSERT_VECTOR_ONLY(TYPE)
Definition: StaticAssert.h:36
#define EIGEN_CHECK_BINARY_COMPATIBILIY(BINOP, LHS, RHS)
Definition: XprHelper.h:900
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:69

◆ EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE()

template<typename Derived >
typedef Eigen::MatrixBase< Derived >::EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE ( ConstStartMinusOne  ,
Scalar  ,
quotient   
)

◆ eigenvalues()

template<typename Derived >
MatrixBase< Derived >::EigenvaluesReturnType Eigen::MatrixBase< Derived >::eigenvalues
inline

Computes the eigenvalues of a matrix.

Returns
Column vector containing the eigenvalues.

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:

VectorXcd eivals = ones.eigenvalues();
cout << "The eigenvalues of the 3x3 matrix of ones are:" << endl << eivals << endl;
VectorXcd eivals
static const ConstantReturnType Ones()
Matrix< std::complex< double >, Dynamic, 1 > VectorXcd
Dynamic×1 vector of type std::complex<double>.
Definition: Matrix.h:504
Matrix< double, Dynamic, Dynamic > MatrixXd
Dynamic×Dynamic matrix of type double.
Definition: Matrix.h:502

Output:

The eigenvalues of the 3x3 matrix of ones are:
(-5.31e-17,0)
        (3,0)
        (0,0)
See also
EigenSolver::eigenvalues(), ComplexEigenSolver::eigenvalues(), SelfAdjointView::eigenvalues()

Definition at line 69 of file MatrixBaseEigenvalues.h.

70 {
71  return internal::eigenvalues_selector<Derived, NumTraits<Scalar>::IsComplex>::run(derived());
72 }

◆ exp()

template<typename Derived >
const MatrixExponentialReturnValue<Derived> Eigen::MatrixBase< Derived >::exp ( ) const

This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise exponential use ArrayBase::exp .

Returns
an expression of the matrix exponential of *this.

◆ forceAlignedAccess() [1/2]

template<typename Derived >
ForceAlignedAccess< Derived > Eigen::MatrixBase< Derived >::forceAlignedAccess
inline
Returns
an expression of *this with forced aligned access
See also
forceAlignedAccessIf(), class ForceAlignedAccess

Definition at line 311 of file MatrixBase.h.

311 { return derived(); }

◆ forceAlignedAccess() [2/2]

template<typename Derived >
const ForceAlignedAccess< Derived > Eigen::MatrixBase< Derived >::forceAlignedAccess
inline
Returns
an expression of *this with forced aligned access
See also
forceAlignedAccessIf(),class ForceAlignedAccess

Definition at line 310 of file MatrixBase.h.

310 { return derived(); }

◆ forceAlignedAccessIf() [1/4]

template<typename Derived >
template<bool Enable>
std::conditional_t<Enable,ForceAlignedAccess<Derived>,Derived&> Eigen::MatrixBase< Derived >::forceAlignedAccessIf ( )
inline
Returns
an expression of *this with forced aligned access if Enable is true.
See also
forceAlignedAccess(), class ForceAlignedAccess

Definition at line 145 of file ForceAlignedAccess.h.

146 {
147  return derived(); // FIXME This should not work but apparently is never used
148 }

◆ forceAlignedAccessIf() [2/4]

template<typename Derived >
template<bool Enable>
Derived& Eigen::MatrixBase< Derived >::forceAlignedAccessIf ( )
inline

Definition at line 313 of file MatrixBase.h.

313 { return derived(); }

◆ forceAlignedAccessIf() [3/4]

template<typename Derived >
template<bool Enable>
add_const_on_value_type_t<std::conditional_t<Enable,ForceAlignedAccess<Derived>,Derived&> > Eigen::MatrixBase< Derived >::forceAlignedAccessIf ( ) const
inline
Returns
an expression of *this with forced aligned access if Enable is true.
See also
forceAlignedAccess(), class ForceAlignedAccess

Definition at line 134 of file ForceAlignedAccess.h.

135 {
136  return derived(); // FIXME This should not work but apparently is never used
137 }

◆ forceAlignedAccessIf() [4/4]

template<typename Derived >
template<bool Enable>
const Derived& Eigen::MatrixBase< Derived >::forceAlignedAccessIf ( ) const
inline

Definition at line 312 of file MatrixBase.h.

312 { return derived(); }

◆ fullPivHouseholderQr() [1/2]

template<typename Derived >
template<typename PermutationIndex = DefaultPermutationIndex>
const FullPivHouseholderQR<PlainObject, PermutationIndex> Eigen::MatrixBase< Derived >::fullPivHouseholderQr ( ) const
inline

◆ fullPivHouseholderQr() [2/2]

template<typename Derived >
template<typename PermutationIndex >
const FullPivHouseholderQR<typename MatrixBase<Derived>::PlainObject, PermutationIndex> Eigen::MatrixBase< Derived >::fullPivHouseholderQr ( ) const
Returns
the full-pivoting Householder QR decomposition of *this.
See also
class FullPivHouseholderQR

Definition at line 731 of file FullPivHouseholderQR.h.

732 {
733  return FullPivHouseholderQR<PlainObject, PermutationIndex>(eval());
734 }

◆ fullPivLu() [1/2]

template<typename Derived >
template<typename PermutationIndex = DefaultPermutationIndex>
const FullPivLU<PlainObject, PermutationIndex> Eigen::MatrixBase< Derived >::fullPivLu ( ) const
inline

◆ fullPivLu() [2/2]

template<typename Derived >
template<typename PermutationIndex >
const FullPivLU<typename MatrixBase<Derived>::PlainObject, PermutationIndex> Eigen::MatrixBase< Derived >::fullPivLu ( ) const
inline

This is defined in the LU module.

#include <Eigen/LU>
Returns
the full-pivoting LU decomposition of *this.
See also
class FullPivLU

Definition at line 870 of file FullPivLU.h.

871 {
872  return FullPivLU<PlainObject, PermutationIndex>(eval());
873 }

◆ householderQr()

template<typename Derived >
const HouseholderQR< typename MatrixBase< Derived >::PlainObject > Eigen::MatrixBase< Derived >::householderQr
inline
Returns
the Householder QR decomposition of *this.
See also
class HouseholderQR

Definition at line 527 of file HouseholderQR.h.

528 {
529  return HouseholderQR<PlainObject>(eval());
530 }

◆ hypotNorm()

template<typename Derived >
NumTraits< typename internal::traits< Derived >::Scalar >::Real Eigen::MatrixBase< Derived >::hypotNorm
inline
Returns
the l2 norm of *this avoiding undeflow and overflow. This version use a concatenation of hypot() calls, and it is very slow.
See also
norm(), stableNorm()

Definition at line 243 of file StableNorm.h.

244 {
245  if(size()==1)
246  return numext::abs(coeff(0,0));
247  else
248  return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
249 }
const CwiseAbsReturnType cwiseAbs() const
CoeffReturnType coeff(Index row, Index col) const
EIGEN_ALWAYS_INLINE std::enable_if_t< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typename NumTraits< T >::Real > abs(const T &x)

◆ Identity() [1/2]

template<typename Derived >
const MatrixBase< Derived >::IdentityReturnType Eigen::MatrixBase< Derived >::Identity
inlinestatic
Returns
an expression of the identity matrix (not necessarily square).

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

Example:

cout << Matrix<double, 3, 4>::Identity() << endl;

Output:

1 0 0 0
0 1 0 0
0 0 1 0
See also
Identity(Index,Index), setIdentity(), isIdentity()

Definition at line 828 of file CwiseNullaryOp.h.

829 {
831  return MatrixBase<Derived>::NullaryExpr(RowsAtCompileTime, ColsAtCompileTime, internal::scalar_identity_op<Scalar>());
832 }
#define EIGEN_STATIC_ASSERT_FIXED_SIZE(TYPE)
Definition: StaticAssert.h:41
static const CwiseNullaryOp< CustomNullaryOp, PlainObject > NullaryExpr(Index rows, Index cols, const CustomNullaryOp &func)

◆ Identity() [2/2]

template<typename Derived >
const MatrixBase< Derived >::IdentityReturnType Eigen::MatrixBase< Derived >::Identity ( Index  rows,
Index  cols 
)
inlinestatic
Returns
an expression of the identity matrix (not necessarily square).

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:

cout << MatrixXd::Identity(4, 3) << endl;
static const IdentityReturnType Identity()

Output:

1 0 0
0 1 0
0 0 1
0 0 0
See also
Identity(), setIdentity(), isIdentity()

Definition at line 811 of file CwiseNullaryOp.h.

812 {
813  return DenseBase<Derived>::NullaryExpr(rows, cols, internal::scalar_identity_op<Scalar>());
814 }

◆ inverse()

template<typename Derived >
const Inverse< Derived > Eigen::MatrixBase< Derived >::inverse
inline

This is defined in the LU module.

#include <Eigen/LU>
Returns
the matrix inverse of this matrix.

For small fixed sizes up to 4x4, this method uses cofactors. In the general case, this method uses class PartialPivLU.

Note
This matrix must be invertible, otherwise the result is undefined. If you need an invertibility check, do the following: Example:
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Its inverse is:" << endl << m.inverse() << endl;
Output:
Here is the matrix m:
  0.68  0.597  -0.33
-0.211  0.823  0.536
 0.566 -0.605 -0.444
Its inverse is:
-0.199   2.23   2.84
  1.01 -0.555  -1.42
 -1.62   3.59   3.29
See also
computeInverseAndDetWithCheck()

Definition at line 350 of file InverseImpl.h.

351 {
352  EIGEN_STATIC_ASSERT(!NumTraits<Scalar>::IsInteger,THIS_FUNCTION_IS_NOT_FOR_INTEGER_NUMERIC_TYPES)
353  eigen_assert(rows() == cols());
354  return Inverse<Derived>(derived());
355 }
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26

◆ isDiagonal()

template<typename Derived >
bool Eigen::MatrixBase< Derived >::isDiagonal ( const RealScalar prec = NumTraits<Scalar>::dummy_precision()) const
Returns
true if *this is approximately equal to a diagonal matrix, within the precision given by prec.

Example:

m(0,2) = 1;
cout << "Here's the matrix m:" << endl << m << endl;
cout << "m.isDiagonal() returns: " << m.isDiagonal() << endl;
cout << "m.isDiagonal(1e-3) returns: " << m.isDiagonal(1e-3) << endl;
Array< double, 1, 3 > e(1./3., 0.5, 2.)

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
See also
asDiagonal()

Definition at line 398 of file DiagonalMatrix.h.

399 {
400  if(cols() != rows()) return false;
401  RealScalar maxAbsOnDiagonal = static_cast<RealScalar>(-1);
402  for(Index j = 0; j < cols(); ++j)
403  {
404  RealScalar absOnDiagonal = numext::abs(coeff(j,j));
405  if(absOnDiagonal > maxAbsOnDiagonal) maxAbsOnDiagonal = absOnDiagonal;
406  }
407  for(Index j = 0; j < cols(); ++j)
408  for(Index i = 0; i < j; ++i)
409  {
410  if(!internal::isMuchSmallerThan(coeff(i, j), maxAbsOnDiagonal, prec)) return false;
411  if(!internal::isMuchSmallerThan(coeff(j, i), maxAbsOnDiagonal, prec)) return false;
412  }
413  return true;
414 }
NumTraits< Scalar >::Real RealScalar
Definition: DenseBase.h:68
bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41

◆ isIdentity()

template<typename Derived >
bool Eigen::MatrixBase< Derived >::isIdentity ( const RealScalar prec = NumTraits<Scalar>::dummy_precision()) const
Returns
true if *this is approximately equal to the identity matrix (not necessarily square), within the precision given by prec.

Example:

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

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
See also
class CwiseNullaryOp, Identity(), Identity(Index,Index), setIdentity()

Definition at line 844 of file CwiseNullaryOp.h.

846 {
847  typename internal::nested_eval<Derived,1>::type self(derived());
848  for(Index j = 0; j < cols(); ++j)
849  {
850  for(Index i = 0; i < rows(); ++i)
851  {
852  if(i == j)
853  {
854  if(!internal::isApprox(self.coeff(i, j), static_cast<Scalar>(1), prec))
855  return false;
856  }
857  else
858  {
859  if(!internal::isMuchSmallerThan(self.coeff(i, j), static_cast<RealScalar>(1), prec))
860  return false;
861  }
862  }
863  }
864  return true;
865 }
bool isApprox(const Scalar &x, const Scalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())

◆ isLowerTriangular()

template<typename Derived >
bool Eigen::MatrixBase< Derived >::isLowerTriangular ( const RealScalar prec = NumTraits<Scalar>::dummy_precision()) const
Returns
true if *this is approximately equal to a lower triangular matrix, within the precision given by prec.
See also
isUpperTriangular()

Definition at line 692 of file TriangularMatrix.h.

693 {
694  RealScalar maxAbsOnLowerPart = static_cast<RealScalar>(-1);
695  for(Index j = 0; j < cols(); ++j)
696  for(Index i = j; i < rows(); ++i)
697  {
698  RealScalar absValue = numext::abs(coeff(i,j));
699  if(absValue > maxAbsOnLowerPart) maxAbsOnLowerPart = absValue;
700  }
701  RealScalar threshold = maxAbsOnLowerPart * prec;
702  for(Index j = 1; j < cols(); ++j)
703  {
704  Index maxi = numext::mini(j, rows()-1);
705  for(Index i = 0; i < maxi; ++i)
706  if(numext::abs(coeff(i, j)) > threshold) return false;
707  }
708  return true;
709 }
EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)

◆ isOrthogonal()

template<typename Derived >
template<typename OtherDerived >
bool Eigen::MatrixBase< Derived >::isOrthogonal ( const MatrixBase< OtherDerived > &  other,
const RealScalar prec = NumTraits<Scalar>::dummy_precision() 
) const
Returns
true if *this is approximately orthogonal to other, within the precision given by prec.

Example:

Vector3d v(1,0,0);
Vector3d w(1e-4,0,1);
cout << "Here's the vector v:" << endl << v << endl;
cout << "Here's the vector w:" << endl << w << endl;
cout << "v.isOrthogonal(w) returns: " << v.isOrthogonal(w) << endl;
cout << "v.isOrthogonal(w,1e-3) returns: " << v.isOrthogonal(w,1e-3) << endl;
Array< int, Dynamic, 1 > v
Matrix< double, 3, 1 > Vector3d
3×1 vector of type double.
Definition: Matrix.h:502

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

Definition at line 280 of file Dot.h.

282 {
283  typename internal::nested_eval<Derived,2>::type nested(derived());
284  typename internal::nested_eval<OtherDerived,2>::type otherNested(other.derived());
285  return numext::abs2(nested.dot(otherNested)) <= prec * prec * nested.squaredNorm() * otherNested.squaredNorm();
286 }
bool abs2(bool x)

◆ isSkewSymmetric()

template<typename Derived >
bool Eigen::MatrixBase< Derived >::isSkewSymmetric ( const RealScalar prec = NumTraits<Scalar>::dummy_precision()) const
Returns
true if *this is approximately equal to a skew symmetric matrix, within the precision given by prec.

Definition at line 358 of file SkewSymmetricMatrix3.h.

359 {
360  if(cols() != rows()) return false;
361  return (this->transpose() + *this).isZero(prec);
362 }

◆ isUnitary()

template<typename Derived >
bool Eigen::MatrixBase< Derived >::isUnitary ( const RealScalar prec = NumTraits<Scalar>::dummy_precision()) const
Returns
true if *this is approximately an unitary matrix, within the precision given by prec. In the case where the Scalar type is real numbers, a unitary matrix is an orthogonal matrix, whence the name.
Note
This can be used to check whether a family of vectors forms an orthonormal basis. Indeed, m.isUnitary() returns true if and only if the columns (equivalently, the rows) of m form an orthonormal basis.

Example:

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

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

Definition at line 300 of file Dot.h.

301 {
302  typename internal::nested_eval<Derived,1>::type self(derived());
303  for(Index i = 0; i < cols(); ++i)
304  {
305  if(!internal::isApprox(self.col(i).squaredNorm(), static_cast<RealScalar>(1), prec))
306  return false;
307  for(Index j = 0; j < i; ++j)
308  if(!internal::isMuchSmallerThan(self.col(i).dot(self.col(j)), static_cast<Scalar>(1), prec))
309  return false;
310  }
311  return true;
312 }
RealScalar squaredNorm() const
Definition: Dot.h:94
ScalarBinaryOpTraits< typename internal::traits< Derived >::Scalar, typename internal::traits< OtherDerived >::Scalar >::ReturnType dot(const MatrixBase< OtherDerived > &other) const
Definition: Dot.h:69

◆ isUpperTriangular()

template<typename Derived >
bool Eigen::MatrixBase< Derived >::isUpperTriangular ( const RealScalar prec = NumTraits<Scalar>::dummy_precision()) const
Returns
true if *this is approximately equal to an upper triangular matrix, within the precision given by prec.
See also
isLowerTriangular()

Definition at line 667 of file TriangularMatrix.h.

668 {
669  RealScalar maxAbsOnUpperPart = static_cast<RealScalar>(-1);
670  for(Index j = 0; j < cols(); ++j)
671  {
672  Index maxi = numext::mini(j, rows()-1);
673  for(Index i = 0; i <= maxi; ++i)
674  {
675  RealScalar absValue = numext::abs(coeff(i,j));
676  if(absValue > maxAbsOnUpperPart) maxAbsOnUpperPart = absValue;
677  }
678  }
679  RealScalar threshold = maxAbsOnUpperPart * prec;
680  for(Index j = 0; j < cols(); ++j)
681  for(Index i = j+1; i < rows(); ++i)
682  if(numext::abs(coeff(i, j)) > threshold) return false;
683  return true;
684 }

◆ jacobiSvd() [1/4]

template<typename Derived >
template<int Options = 0>
JacobiSVD<PlainObject, Options> Eigen::MatrixBase< Derived >::jacobiSvd ( ) const
inline

◆ jacobiSvd() [2/4]

template<typename Derived >
template<int Options>
JacobiSVD<typename MatrixBase<Derived>::PlainObject, Options> Eigen::MatrixBase< Derived >::jacobiSvd ( ) const

This is defined in the SVD module.

#include <Eigen/SVD>
Returns
the singular value decomposition of *this computed by two-sided Jacobi transformations.
See also
class JacobiSVD

Definition at line 830 of file JacobiSVD.h.

830  {
831  return JacobiSVD<PlainObject, Options>(*this);
832 }

◆ jacobiSvd() [3/4]

template<typename Derived >
template<int Options = 0>
EIGEN_DEPRECATED JacobiSVD<PlainObject, Options> Eigen::MatrixBase< Derived >::jacobiSvd ( unsigned int  computationOptions) const
inline

◆ jacobiSvd() [4/4]

template<typename Derived >
template<int Options>
JacobiSVD<typename MatrixBase<Derived>::PlainObject, Options> Eigen::MatrixBase< Derived >::jacobiSvd ( unsigned int  computationOptions) const

Definition at line 836 of file JacobiSVD.h.

837  {
838  return JacobiSVD<PlainObject, Options>(*this, computationOptions);
839 }

◆ lazyProduct()

template<typename Derived >
template<typename OtherDerived >
const Product< Derived, OtherDerived, LazyProduct > Eigen::MatrixBase< Derived >::lazyProduct ( const MatrixBase< OtherDerived > &  other) const
inline
Returns
an expression of the matrix product of *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.

Warning
This version of the matrix product can be much much slower. So use it only if you know what you are doing and that you measured a true speed improvement.
See also
operator*(const MatrixBase&)

Definition at line 444 of file GeneralProduct.h.

445 {
446  enum {
447  ProductIsValid = Derived::ColsAtCompileTime==Dynamic
448  || OtherDerived::RowsAtCompileTime==Dynamic
449  || int(Derived::ColsAtCompileTime)==int(OtherDerived::RowsAtCompileTime),
450  AreVectors = Derived::IsVectorAtCompileTime && OtherDerived::IsVectorAtCompileTime,
451  SameSizes = EIGEN_PREDICATE_SAME_MATRIX_SIZE(Derived,OtherDerived)
452  };
453  // note to the lost user:
454  // * for a dot product use: v1.dot(v2)
455  // * for a coeff-wise product use: v1.cwiseProduct(v2)
456  EIGEN_STATIC_ASSERT(ProductIsValid || !(AreVectors && SameSizes),
457  INVALID_VECTOR_VECTOR_PRODUCT__IF_YOU_WANTED_A_DOT_OR_COEFF_WISE_PRODUCT_YOU_MUST_USE_THE_EXPLICIT_FUNCTIONS)
458  EIGEN_STATIC_ASSERT(ProductIsValid || !(SameSizes && !AreVectors),
459  INVALID_MATRIX_PRODUCT__IF_YOU_WANTED_A_COEFF_WISE_PRODUCT_YOU_MUST_USE_THE_EXPLICIT_FUNCTION)
460  EIGEN_STATIC_ASSERT(ProductIsValid || SameSizes, INVALID_MATRIX_PRODUCT)
461 
462  return Product<Derived,OtherDerived,LazyProduct>(derived(), other.derived());
463 }
#define EIGEN_PREDICATE_SAME_MATRIX_SIZE(TYPE0, TYPE1)
Definition: StaticAssert.h:68

◆ ldlt()

template<typename Derived >
const LDLT< typename MatrixBase< Derived >::PlainObject > Eigen::MatrixBase< Derived >::ldlt
inline

This is defined in the Cholesky module.

Returns
the Cholesky decomposition with full pivoting without square root of *this
See also
SelfAdjointView::ldlt()

Definition at line 678 of file LDLT.h.

679 {
680  return LDLT<PlainObject>(derived());
681 }

◆ llt()

template<typename Derived >
const LLT< typename MatrixBase< Derived >::PlainObject > Eigen::MatrixBase< Derived >::llt
inline

This is defined in the Cholesky module.

#include <Eigen/Cholesky>
Returns
the LLT decomposition of *this
See also
SelfAdjointView::llt()

Definition at line 537 of file LLT.h.

538 {
539  return LLT<PlainObject>(derived());
540 }

◆ log()

template<typename Derived >
const MatrixLogarithmReturnValue<Derived> Eigen::MatrixBase< Derived >::log ( ) const

This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise logarithm use ArrayBase::log .

Returns
an expression of the matrix logarithm of *this.

◆ lpNorm()

template<typename Derived >
template<int p>
MatrixBase< Derived >::RealScalar Eigen::MatrixBase< Derived >::lpNorm
Returns
the coefficient-wise \( \ell^p \) norm of *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 \( \ell^\infty \) norm, that is the maximum of the absolute values of the coefficients of *this.

In all cases, if *this is empty, then the value 0 is returned.

Note
For matrices, this function does not compute the operator-norm. That is, if *this is a matrix, then its coefficients are interpreted as a 1D vector. Nonetheless, you can easily compute the 1-norm and \(\infty\)-norm matrix operator norms using partial reductions .
See also
norm()

Definition at line 265 of file Dot.h.

266 {
267  return internal::lpNorm_selector<Derived, p>::run(*this);
268 }

◆ lu() [1/2]

template<typename Derived >
template<typename PermutationIndex = DefaultPermutationIndex>
const PartialPivLU<PlainObject, PermutationIndex> Eigen::MatrixBase< Derived >::lu ( ) const
inline

◆ lu() [2/2]

template<typename Derived >
template<typename PermutationIndex >
const PartialPivLU<typename MatrixBase<Derived>::PlainObject, PermutationIndex> Eigen::MatrixBase< Derived >::lu ( ) const
inline

This is defined in the LU module.

#include <Eigen/LU>

Synonym of partialPivLu().

Returns
the partial-pivoting LU decomposition of *this.
See also
class PartialPivLU

Definition at line 616 of file PartialPivLU.h.

617 {
618  return PartialPivLU<PlainObject, PermutationIndex>(eval());
619 }

◆ makeHouseholder()

template<typename Derived >
template<typename EssentialPart >
void Eigen::MatrixBase< Derived >::makeHouseholder ( EssentialPart &  essential,
Scalar tau,
RealScalar beta 
) const

Computes the elementary reflector H such that: \( H *this = [ beta 0 ... 0]^T \) where the transformation H is: \( H = I - tau v v^*\) and the vector v is: \( v^T = [1 essential^T] \)

On output:

Parameters
essentialthe essential part of the vector v
tauthe scaling factor of the Householder transformation
betathe result of H * *this
See also
MatrixBase::makeHouseholderInPlace(), MatrixBase::applyHouseholderOnTheLeft(), MatrixBase::applyHouseholderOnTheRight()

Definition at line 69 of file Householder.h.

73 {
74  using numext::sqrt;
75  using numext::conj;
76 
77  EIGEN_STATIC_ASSERT_VECTOR_ONLY(EssentialPart)
78  VectorBlock<const Derived, EssentialPart::SizeAtCompileTime> tail(derived(), 1, size()-1);
79 
80  RealScalar tailSqNorm = size()==1 ? RealScalar(0) : tail.squaredNorm();
81  Scalar c0 = coeff(0);
83 
84  if(tailSqNorm <= tol && numext::abs2(numext::imag(c0))<=tol)
85  {
86  tau = RealScalar(0);
87  beta = numext::real(c0);
88  essential.setZero();
89  }
90  else
91  {
92  beta = sqrt(numext::abs2(c0) + tailSqNorm);
93  if (numext::real(c0)>=RealScalar(0))
94  beta = -beta;
95  essential = tail / (c0 - beta);
96  tau = conj((beta - c0) / beta);
97  }
98 }
FixedSegmentReturnType<... >::Type tail(NType n)
const ImagReturnType imag() const
RealReturnType real() const
const MatrixSquareRootReturnValue< Derived > sqrt() const
This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise square...
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
EIGEN_ALWAYS_INLINE float sqrt(const float &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)

◆ makeHouseholderInPlace()

template<typename Derived >
void Eigen::MatrixBase< Derived >::makeHouseholderInPlace ( Scalar tau,
RealScalar beta 
)

Computes the elementary reflector H such that: \( H *this = [ beta 0 ... 0]^T \) where the transformation H is: \( H = I - tau v v^*\) and the vector v is: \( v^T = [1 essential^T] \)

The essential part of the vector v is stored in *this.

On output:

Parameters
tauthe scaling factor of the Householder transformation
betathe result of H * *this
See also
MatrixBase::makeHouseholder(), MatrixBase::applyHouseholderOnTheLeft(), MatrixBase::applyHouseholderOnTheRight()

Definition at line 45 of file Householder.h.

46 {
47  VectorBlock<Derived, internal::decrement_size<Base::SizeAtCompileTime>::ret> essentialPart(derived(), 1, size()-1);
48  makeHouseholder(essentialPart, tau, beta);
49 }
void makeHouseholder(EssentialPart &essential, Scalar &tau, RealScalar &beta) const
Definition: Householder.h:69

◆ matrix() [1/2]

template<typename Derived >
MatrixBase<Derived>& Eigen::MatrixBase< Derived >::matrix ( )
inline

Definition at line 319 of file MatrixBase.h.

319 { return *this; }

◆ matrix() [2/2]

template<typename Derived >
const MatrixBase<Derived>& Eigen::MatrixBase< Derived >::matrix ( ) const
inline

Definition at line 320 of file MatrixBase.h.

320 { return *this; }

◆ matrixFunction()

template<typename Derived >
const MatrixFunctionReturnValue<Derived> Eigen::MatrixBase< Derived >::matrixFunction ( StemFunction  f) const

Helper function for the unsupported MatrixFunctions module.

◆ noalias()

template<typename Derived >
NoAlias< Derived, MatrixBase > Eigen::MatrixBase< Derived >::noalias
Returns
a pseudo expression of *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:

D.noalias() = A * B;
D.noalias() += A.transpose() * B;
D.noalias() -= 2 * A * B.adjoint();

On the other hand the following example will lead to a wrong result:

A.noalias() = A * B;

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:

A = A * B;
See also
class NoAlias

Definition at line 104 of file NoAlias.h.

105 {
106  return NoAlias<Derived, Eigen::MatrixBase >(derived());
107 }

◆ norm()

template<typename Derived >
NumTraits< typename internal::traits< Derived >::Scalar >::Real Eigen::MatrixBase< Derived >::norm
inline
Returns
, for vectors, the l2 norm of *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.
See also
lpNorm(), dot(), squaredNorm()

Definition at line 106 of file Dot.h.

107 {
108  return numext::sqrt(squaredNorm());
109 }

◆ normalize()

template<typename Derived >
void Eigen::MatrixBase< Derived >::normalize
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.

Warning
If the input vector is too small (i.e., this->norm()==0), then *this is left unchanged.
See also
norm(), normalized()

Definition at line 143 of file Dot.h.

144 {
146  // NOTE: after extensive benchmarking, this conditional does not impact performance, at least on recent x86 CPU
147  if(z>RealScalar(0))
148  derived() /= numext::sqrt(z);
149 }

◆ normalized()

template<typename Derived >
const MatrixBase< Derived >::PlainObject Eigen::MatrixBase< Derived >::normalized
inline
Returns
an expression of the quotient of *this by its own norm.
Warning
If the input vector is too small (i.e., this->norm()==0), then this function returns a copy of the input.

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
norm(), normalize()

Definition at line 122 of file Dot.h.

123 {
124  typedef typename internal::nested_eval<Derived,2>::type Nested_;
125  Nested_ n(derived());
126  RealScalar z = n.squaredNorm();
127  // NOTE: after extensive benchmarking, this conditional does not impact performance, at least on recent x86 CPU
128  if(z>RealScalar(0))
129  return n / numext::sqrt(z);
130  else
131  return n;
132 }
int n

◆ operator!=()

template<typename Derived >
template<typename OtherDerived >
bool Eigen::MatrixBase< Derived >::operator!= ( const MatrixBase< OtherDerived > &  other) const
inline
Returns
true if at least one pair of coefficients of *this and other are not exactly equal to each other.
Warning
When using floating point scalar values you probably should rather use a fuzzy comparison such as isApprox()
See also
isApprox(), operator==

Definition at line 303 of file MatrixBase.h.

304  { return cwiseNotEqual(other).any(); }
const CwiseBinaryNotEqualReturnType< OtherDerived > cwiseNotEqual(const EIGEN_CURRENT_STORAGE_BASE_CLASS< OtherDerived > &other) const

◆ operator*() [1/3]

template<typename Derived >
template<typename DiagonalDerived >
const Product< Derived, DiagonalDerived, LazyProduct > Eigen::MatrixBase< Derived >::operator* ( const DiagonalBase< DiagonalDerived > &  a_diagonal) const
inline
Returns
the diagonal matrix product of *this by the diagonal matrix diagonal.

Definition at line 23 of file DiagonalProduct.h.

24 {
25  return Product<Derived, DiagonalDerived, LazyProduct>(derived(),a_diagonal.derived());
26 }

◆ operator*() [2/3]

template<typename Derived >
template<typename OtherDerived >
const Product< Derived, OtherDerived > Eigen::MatrixBase< Derived >::operator* ( const MatrixBase< OtherDerived > &  other) const
inline

Implementation of matrix base methods

Returns
the matrix product of *this and other.
Note
If instead of the matrix product you want the coefficient-wise product, see Cwise::operator*().
See also
lazyProduct(), operator*=(const MatrixBase&), Cwise::operator*()

Definition at line 401 of file GeneralProduct.h.

402 {
403  // A note regarding the function declaration: In MSVC, this function will sometimes
404  // not be inlined since DenseStorage is an unwindable object for dynamic
405  // matrices and product types are holding a member to store the result.
406  // Thus it does not help tagging this function with EIGEN_STRONG_INLINE.
407  enum {
408  ProductIsValid = Derived::ColsAtCompileTime==Dynamic
409  || OtherDerived::RowsAtCompileTime==Dynamic
410  || int(Derived::ColsAtCompileTime)==int(OtherDerived::RowsAtCompileTime),
411  AreVectors = Derived::IsVectorAtCompileTime && OtherDerived::IsVectorAtCompileTime,
412  SameSizes = EIGEN_PREDICATE_SAME_MATRIX_SIZE(Derived,OtherDerived)
413  };
414  // note to the lost user:
415  // * for a dot product use: v1.dot(v2)
416  // * for a coeff-wise product use: v1.cwiseProduct(v2)
417  EIGEN_STATIC_ASSERT(ProductIsValid || !(AreVectors && SameSizes),
418  INVALID_VECTOR_VECTOR_PRODUCT__IF_YOU_WANTED_A_DOT_OR_COEFF_WISE_PRODUCT_YOU_MUST_USE_THE_EXPLICIT_FUNCTIONS)
419  EIGEN_STATIC_ASSERT(ProductIsValid || !(SameSizes && !AreVectors),
420  INVALID_MATRIX_PRODUCT__IF_YOU_WANTED_A_COEFF_WISE_PRODUCT_YOU_MUST_USE_THE_EXPLICIT_FUNCTION)
421  EIGEN_STATIC_ASSERT(ProductIsValid || SameSizes, INVALID_MATRIX_PRODUCT)
422 #ifdef EIGEN_DEBUG_PRODUCT
423  internal::product_type<Derived,OtherDerived>::debug();
424 #endif
425 
426  return Product<Derived, OtherDerived>(derived(), other.derived());
427 }

◆ operator*() [3/3]

template<typename Derived >
template<typename SkewDerived >
const Product< Derived, SkewDerived, LazyProduct > Eigen::MatrixBase< Derived >::operator* ( const SkewSymmetricBase< SkewDerived > &  skew) const
inline
Returns
the matrix product of *this by the skew symmetric matrix \skew.

Definition at line 369 of file SkewSymmetricMatrix3.h.

370 {
371  return Product<Derived, SkewDerived, LazyProduct>(derived(), skew.derived());
372 }

◆ operator*=()

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

Implementation of matrix base methods replaces *this by *this * other.

Returns
a reference to *this

Example:

B << 0,1,0,
0,0,1,
1,0,0;
cout << "At start, A = " << endl << A << endl;
A *= B;
cout << "After A *= B, A = " << endl << A << endl;
A.applyOnTheRight(B); // equivalent to A *= B
cout << "After applyOnTheRight, A = " << endl << A << endl;

Output:

At start, A = 
  0.68  0.597  -0.33
-0.211  0.823  0.536
 0.566 -0.605 -0.444
After A *= B, A = 
 -0.33   0.68  0.597
 0.536 -0.211  0.823
-0.444  0.566 -0.605
After applyOnTheRight, A = 
 0.597  -0.33   0.68
 0.823  0.536 -0.211
-0.605 -0.444  0.566

Definition at line 519 of file MatrixBase.h.

520 {
521  other.derived().applyThisOnTheRight(derived());
522  return derived();
523 }

◆ operator+=() [1/2]

template<typename Derived >
template<typename OtherDerived >
Derived& Eigen::MatrixBase< Derived >::operator+= ( const ArrayBase< OtherDerived > &  )
inlineprotected

Definition at line 497 of file MatrixBase.h.

498  {EIGEN_STATIC_ASSERT(std::ptrdiff_t(sizeof(typename OtherDerived::Scalar))==-1,YOU_CANNOT_MIX_ARRAYS_AND_MATRICES); return *this;}

◆ operator+=() [2/2]

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

replaces *this by *this + other.

Returns
a reference to *this

Definition at line 177 of file CwiseBinaryOp.h.

178 {
179  call_assignment(derived(), other.derived(), internal::add_assign_op<Scalar,typename OtherDerived::Scalar>());
180  return derived();
181 }
void call_assignment(Dst &dst, const Src &src)

◆ operator-=() [1/2]

template<typename Derived >
template<typename OtherDerived >
Derived& Eigen::MatrixBase< Derived >::operator-= ( const ArrayBase< OtherDerived > &  )
inlineprotected

Definition at line 500 of file MatrixBase.h.

501  {EIGEN_STATIC_ASSERT(std::ptrdiff_t(sizeof(typename OtherDerived::Scalar))==-1,YOU_CANNOT_MIX_ARRAYS_AND_MATRICES); return *this;}

◆ operator-=() [2/2]

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

replaces *this by *this - other.

Returns
a reference to *this

Definition at line 164 of file CwiseBinaryOp.h.

165 {
166  call_assignment(derived(), other.derived(), internal::sub_assign_op<Scalar,typename OtherDerived::Scalar>());
167  return derived();
168 }

◆ operator=() [1/4]

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

Definition at line 66 of file Assign.h.

67 {
68  internal::call_assignment(derived(), other.derived());
69  return derived();
70 }

◆ operator=() [2/4]

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

Definition at line 75 of file Assign.h.

76 {
77  internal::call_assignment(derived(), other.derived());
78  return derived();
79 }

◆ operator=() [3/4]

template<typename Derived >
Derived & Eigen::MatrixBase< Derived >::operator= ( const MatrixBase< 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)

Definition at line 57 of file Assign.h.

58 {
59  internal::call_assignment(derived(), other.derived());
60  return derived();
61 }

◆ operator=() [4/4]

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

Definition at line 84 of file Assign.h.

85 {
86  other.derived().evalTo(derived());
87  return derived();
88 }

◆ operator==()

template<typename Derived >
template<typename OtherDerived >
bool Eigen::MatrixBase< Derived >::operator== ( const MatrixBase< OtherDerived > &  other) const
inline
Returns
true if each coefficients of *this and other are all exactly equal.
Warning
When using floating point scalar values you probably should rather use a fuzzy comparison such as isApprox()
See also
isApprox(), operator!=

Definition at line 295 of file MatrixBase.h.

296  { return cwiseEqual(other).all(); }
const CwiseBinaryEqualReturnType< OtherDerived > cwiseEqual(const EIGEN_CURRENT_STORAGE_BASE_CLASS< OtherDerived > &other) const

◆ operatorNorm()

template<typename Derived >
MatrixBase< Derived >::RealScalar Eigen::MatrixBase< Derived >::operatorNorm
inline

Computes the L2 operator norm.

Returns
Operator norm of the matrix.

This is defined in the Eigenvalues module.

#include <Eigen/Eigenvalues>

This function computes the L2 operator norm of a matrix, which is also known as the spectral norm. The norm of a matrix \( A \) is defined to be

\[ \|A\|_2 = \max_x \frac{\|Ax\|_2}{\|x\|_2} \]

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 \( A^*A \).

The current implementation uses the eigenvalues of \( A^*A \), as computed by SelfAdjointView::eigenvalues(), to compute the operator norm of a matrix. The SelfAdjointView class provides a better algorithm for selfadjoint matrices.

Example:

cout << "The operator norm of the 3x3 matrix of ones is "
<< ones.operatorNorm() << endl;

Output:

The operator norm of the 3x3 matrix of ones is 3
See also
SelfAdjointView::eigenvalues(), SelfAdjointView::operatorNorm()

Definition at line 122 of file MatrixBaseEigenvalues.h.

123 {
124  using std::sqrt;
125  typename Derived::PlainObject m_eval(derived());
126  // FIXME if it is really guaranteed that the eigenvalues are already sorted,
127  // then we don't need to compute a maxCoeff() here, comparing the 1st and last ones is enough.
128  return sqrt((m_eval*m_eval.adjoint())
129  .eval()
130  .template selfadjointView<Lower>()
131  .eigenvalues()
132  .maxCoeff()
133  );
134 }
const SqrtReturnType sqrt() const

◆ partialPivLu() [1/2]

template<typename Derived >
template<typename PermutationIndex = DefaultPermutationIndex>
const PartialPivLU<PlainObject, PermutationIndex> Eigen::MatrixBase< Derived >::partialPivLu ( ) const
inline

◆ partialPivLu() [2/2]

template<typename Derived >
template<typename PermutationIndex >
const PartialPivLU<typename MatrixBase<Derived>::PlainObject, PermutationIndex> Eigen::MatrixBase< Derived >::partialPivLu ( ) const
inline

This is defined in the LU module.

#include <Eigen/LU>
Returns
the partial-pivoting LU decomposition of *this.
See also
class PartialPivLU

Definition at line 600 of file PartialPivLU.h.

601 {
602  return PartialPivLU<PlainObject, PermutationIndex>(eval());
603 }

◆ pow() [1/2]

template<typename Derived >
const MatrixPowerReturnValue<Derived> Eigen::MatrixBase< Derived >::pow ( const RealScalar p) const

This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise power to p use ArrayBase::pow .

Returns
an expression of the matrix power to p of *this.

◆ pow() [2/2]

template<typename Derived >
const MatrixComplexPowerReturnValue<Derived> Eigen::MatrixBase< Derived >::pow ( const std::complex< RealScalar > &  p) const

This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise power to p use ArrayBase::pow .

Returns
an expression of the matrix power to p of *this.

◆ selfadjointView() [1/4]

template<typename Derived >
template<unsigned int UpLo>
SelfAdjointViewReturnType<UpLo>::Type Eigen::MatrixBase< Derived >::selfadjointView ( )

◆ selfadjointView() [2/4]

template<typename Derived >
template<unsigned int UpLo>
MatrixBase<Derived>::template SelfAdjointViewReturnType<UpLo>::Type Eigen::MatrixBase< Derived >::selfadjointView ( )
Returns
an expression of a symmetric/self-adjoint view extracted from the upper or lower triangular part of the current matrix

The parameter UpLo can be either Upper or Lower

Example:

cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the symmetric matrix extracted from the upper part of m:" << endl
<< Matrix3i(m.selfadjointView<Upper>()) << endl;
cout << "Here is the symmetric matrix extracted from the lower part of m:" << endl
<< Matrix3i(m.selfadjointView<Lower>()) << endl;
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213

Output:

Here is the matrix m:
 7  6 -3
-2  9  6
 6 -6 -5
Here is the symmetric matrix extracted from the upper part of m:
 7  6 -3
 6  9  6
-3  6 -5
Here is the symmetric matrix extracted from the lower part of m:
 7 -2  6
-2  9 -6
 6 -6 -5
See also
class SelfAdjointView

Definition at line 358 of file SelfAdjointView.h.

359 {
361 }
SelfAdjointView< Derived, UpLo > Type
Definition: MatrixBase.h:243

◆ selfadjointView() [3/4]

template<typename Derived >
template<unsigned int UpLo>
ConstSelfAdjointViewReturnType<UpLo>::Type Eigen::MatrixBase< Derived >::selfadjointView ( ) const

◆ selfadjointView() [4/4]

template<typename Derived >
template<unsigned int UpLo>
MatrixBase<Derived>::template ConstSelfAdjointViewReturnType<UpLo>::Type Eigen::MatrixBase< Derived >::selfadjointView ( ) const

Implementation of MatrixBase methods This is the const version of MatrixBase::selfadjointView()

Definition at line 341 of file SelfAdjointView.h.

342 {
344 }
const SelfAdjointView< const Derived, UpLo > Type
Definition: MatrixBase.h:244

◆ setIdentity() [1/2]

template<typename Derived >
Derived & Eigen::MatrixBase< Derived >::setIdentity
inline

Writes the identity expression (not necessarily square) into *this.

Example:

m.block<3,3>(1,0).setIdentity();
cout << m << endl;
static const ConstantReturnType Zero()
Derived & setIdentity()

Output:

0 0 0 0
1 0 0 0
0 1 0 0
0 0 1 0
See also
class CwiseNullaryOp, Identity(), Identity(Index,Index), isIdentity()

Definition at line 902 of file CwiseNullaryOp.h.

903 {
904  return internal::setIdentity_impl<Derived>::run(derived());
905 }

◆ setIdentity() [2/2]

template<typename Derived >
Derived & Eigen::MatrixBase< Derived >::setIdentity ( Index  rows,
Index  cols 
)
inline

Resizes to the given size, and writes the identity expression (not necessarily square) into *this.

Parameters
rowsthe new number of rows
colsthe new number of columns

Example:

m.setIdentity(3, 3);
cout << m << endl;
Matrix< float, Dynamic, Dynamic > MatrixXf
Dynamic×Dynamic matrix of type float.
Definition: Matrix.h:501

Output:

1 0 0
0 1 0
0 0 1
See also
MatrixBase::setIdentity(), class CwiseNullaryOp, MatrixBase::Identity()

Definition at line 918 of file CwiseNullaryOp.h.

919 {
920  derived().resize(rows, cols);
921  return setIdentity();
922 }

◆ setUnit() [1/2]

template<typename Derived >
Derived & Eigen::MatrixBase< Derived >::setUnit ( Index  i)
inline

Set the coefficients of *this to the i-th unit (basis) vector.

Parameters
iindex 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.

See also
MatrixBase::setIdentity(), class CwiseNullaryOp, MatrixBase::Unit(Index,Index)

Definition at line 1001 of file CwiseNullaryOp.h.

1002 {
1004  eigen_assert(i<size());
1005  derived().setZero();
1006  derived().coeffRef(i) = Scalar(1);
1007  return derived();
1008 }

◆ setUnit() [2/2]

template<typename Derived >
Derived & Eigen::MatrixBase< Derived >::setUnit ( Index  newSize,
Index  i 
)
inline

Resizes to the given newSize, and writes the i-th unit (basis) vector into *this.

Parameters
newSizethe new size of the vector
iindex 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.

See also
MatrixBase::setIdentity(), class CwiseNullaryOp, MatrixBase::Unit(Index,Index)

Definition at line 1020 of file CwiseNullaryOp.h.

1021 {
1023  eigen_assert(i<newSize);
1024  derived().resize(newSize);
1025  return setUnit(i);
1026 }
Derived & setUnit(Index i)
Set the coefficients of *this to the i-th unit (basis) vector.

◆ sin()

template<typename Derived >
const MatrixFunctionReturnValue<Derived> Eigen::MatrixBase< Derived >::sin ( ) const

This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise sine use ArrayBase::sin .

Returns
an expression of the matrix sine of *this.

◆ sinh()

template<typename Derived >
const MatrixFunctionReturnValue<Derived> Eigen::MatrixBase< Derived >::sinh ( ) const

This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise hyperbolic sine use ArrayBase::sinh .

Returns
an expression of the matrix hyperbolic sine of *this.

◆ sqrt()

template<typename Derived >
const MatrixSquareRootReturnValue<Derived> Eigen::MatrixBase< Derived >::sqrt ( ) const

This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise square root use ArrayBase::sqrt .

Returns
an expression of the matrix square root of *this.

◆ squaredNorm()

template<typename Derived >
NumTraits< typename internal::traits< Derived >::Scalar >::Real Eigen::MatrixBase< Derived >::squaredNorm
inline
Returns
, for vectors, the squared l2 norm of *this, and for matrices the squared Frobenius norm. In both cases, it consists in the sum of the square of all the matrix entries. For vectors, this is also equals to the dot product of *this with itself.
See also
dot(), norm(), lpNorm()

Definition at line 94 of file Dot.h.

95 {
96  return numext::real((*this).cwiseAbs2().sum());
97 }

◆ stableNorm()

template<typename Derived >
NumTraits< typename internal::traits< Derived >::Scalar >::Real Eigen::MatrixBase< Derived >::stableNorm
inline
Returns
the l2 norm of *this avoiding underflow and overflow. This version use a blockwise two passes algorithm: 1 - find the absolute largest coefficient s 2 - compute \( s \Vert \frac{*this}{s} \Vert \) in a standard way

For architecture/scalar types supporting vectorization, this version is faster than blueNorm(). Otherwise the blueNorm() is much faster.

See also
norm(), blueNorm(), hypotNorm()

Definition at line 215 of file StableNorm.h.

216 {
218 }
VectorType::RealScalar stable_norm_impl(const VectorType &vec, std::enable_if_t< VectorType::IsVectorAtCompileTime > *=0)
Definition: StableNorm.h:84

◆ stableNormalize()

template<typename Derived >
void Eigen::MatrixBase< Derived >::stableNormalize
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.

Warning
If the input vector is too small (i.e., this->norm()==0), then *this is left unchanged.
See also
stableNorm(), stableNormalized(), normalize()

Definition at line 189 of file Dot.h.

190 {
191  RealScalar w = cwiseAbs().maxCoeff();
192  RealScalar z = (derived()/w).squaredNorm();
193  if(z>RealScalar(0))
194  derived() /= numext::sqrt(z)*w;
195 }

◆ stableNormalized()

template<typename Derived >
const MatrixBase< Derived >::PlainObject Eigen::MatrixBase< Derived >::stableNormalized
inline
Returns
an expression of the quotient of *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.

Warning
If the input vector is too small (i.e., this->norm()==0), then this function returns a copy of the input.
See also
stableNorm(), stableNormalize(), normalized()

Definition at line 165 of file Dot.h.

166 {
167  typedef typename internal::nested_eval<Derived,3>::type Nested_;
168  Nested_ n(derived());
169  RealScalar w = n.cwiseAbs().maxCoeff();
170  RealScalar z = (n/w).squaredNorm();
171  if(z>RealScalar(0))
172  return n / (numext::sqrt(z)*w);
173  else
174  return n;
175 }

◆ trace()

template<typename Derived >
internal::traits< Derived >::Scalar Eigen::MatrixBase< Derived >::trace
inline
Returns
the trace of *this, i.e. the sum of the coefficients on the main diagonal.

*this can be any matrix, not necessarily square.

See also
diagonal(), sum()

Definition at line 595 of file Redux.h.

596 {
597  return derived().diagonal().sum();
598 }

◆ triangularView() [1/4]

template<typename Derived >
template<unsigned int Mode>
TriangularViewReturnType<Mode>::Type Eigen::MatrixBase< Derived >::triangularView ( )

◆ triangularView() [2/4]

template<typename Derived >
template<unsigned int Mode>
MatrixBase<Derived>::template TriangularViewReturnType<Mode>::Type Eigen::MatrixBase< Derived >::triangularView ( )

Implementation of TriangularView methods Implementation of MatrixBase methods

Returns
an expression of a triangular view extracted from the current matrix

The parameter Mode can have the following values: Upper, StrictlyUpper, UnitUpper, Lower, StrictlyLower, UnitLower.

Example:

cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the upper-triangular matrix extracted from m:" << endl
<< Matrix3i(m.triangularView<Eigen::Upper>()) << endl;
cout << "Here is the strictly-upper-triangular matrix extracted from m:" << endl
<< Matrix3i(m.triangularView<Eigen::StrictlyUpper>()) << endl;
cout << "Here is the unit-lower-triangular matrix extracted from m:" << endl
<< Matrix3i(m.triangularView<Eigen::UnitLower>()) << endl;
// FIXME need to implement output for triangularViews (Bug 885)
@ StrictlyUpper
Definition: Constants.h:225
@ UnitLower
Definition: Constants.h:219

Output:

Here is the matrix m:
 7  6 -3
-2  9  6
 6 -6 -5
Here is the upper-triangular matrix extracted from m:
 7  6 -3
 0  9  6
 0  0 -5
Here is the strictly-upper-triangular matrix extracted from m:
 0  6 -3
 0  0  6
 0  0  0
Here is the unit-lower-triangular matrix extracted from m:
 1  0  0
-2  1  0
 6 -6  1
See also
class TriangularView

Definition at line 646 of file TriangularMatrix.h.

647 {
649 }
TriangularView< Derived, Mode > Type
Definition: MatrixBase.h:233

◆ triangularView() [3/4]

template<typename Derived >
template<unsigned int Mode>
ConstTriangularViewReturnType<Mode>::Type Eigen::MatrixBase< Derived >::triangularView ( ) const

◆ triangularView() [4/4]

template<typename Derived >
template<unsigned int Mode>
MatrixBase<Derived>::template ConstTriangularViewReturnType<Mode>::Type Eigen::MatrixBase< Derived >::triangularView ( ) const

This is the const version of MatrixBase::triangularView()

Definition at line 656 of file TriangularMatrix.h.

657 {
659 }
const TriangularView< const Derived, Mode > Type
Definition: MatrixBase.h:234

◆ Unit() [1/2]

template<typename Derived >
const MatrixBase< Derived >::BasisReturnType Eigen::MatrixBase< Derived >::Unit ( Index  i)
inlinestatic
Returns
an expression of the i-th unit (basis) vector.

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.

See also
MatrixBase::Unit(Index,Index), MatrixBase::UnitX(), MatrixBase::UnitY(), MatrixBase::UnitZ(), MatrixBase::UnitW()

Definition at line 946 of file CwiseNullaryOp.h.

947 {
949  return BasisReturnType(SquareMatrixType::Identity(),i);
950 }

◆ Unit() [2/2]

template<typename Derived >
const MatrixBase< Derived >::BasisReturnType Eigen::MatrixBase< Derived >::Unit ( Index  newSize,
Index  i 
)
inlinestatic
Returns
an expression of the i-th unit (basis) vector.

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
MatrixBase::Unit(Index), MatrixBase::UnitX(), MatrixBase::UnitY(), MatrixBase::UnitZ(), MatrixBase::UnitW()

Definition at line 931 of file CwiseNullaryOp.h.

932 {
934  return BasisReturnType(SquareMatrixType::Identity(newSize,newSize), i);
935 }

◆ UnitW()

template<typename Derived >
const MatrixBase< Derived >::BasisReturnType Eigen::MatrixBase< Derived >::UnitW
inlinestatic
Returns
an expression of the W axis unit vector (0,0,0,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.

See also
MatrixBase::Unit(Index,Index), MatrixBase::Unit(Index), MatrixBase::UnitY(), MatrixBase::UnitZ(), MatrixBase::UnitW()

Definition at line 989 of file CwiseNullaryOp.h.

990 { return Derived::Unit(3); }

◆ UnitX()

template<typename Derived >
const MatrixBase< Derived >::BasisReturnType Eigen::MatrixBase< Derived >::UnitX
inlinestatic
Returns
an expression of the X axis unit vector (1{,0}^*)

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
MatrixBase::Unit(Index,Index), MatrixBase::Unit(Index), MatrixBase::UnitY(), MatrixBase::UnitZ(), MatrixBase::UnitW()

Definition at line 959 of file CwiseNullaryOp.h.

960 { return Derived::Unit(0); }

◆ UnitY()

template<typename Derived >
const MatrixBase< Derived >::BasisReturnType Eigen::MatrixBase< Derived >::UnitY
inlinestatic
Returns
an expression of the Y axis unit vector (0,1{,0}^*)

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
MatrixBase::Unit(Index,Index), MatrixBase::Unit(Index), MatrixBase::UnitY(), MatrixBase::UnitZ(), MatrixBase::UnitW()

Definition at line 969 of file CwiseNullaryOp.h.

970 { return Derived::Unit(1); }

◆ UnitZ()

template<typename Derived >
const MatrixBase< Derived >::BasisReturnType Eigen::MatrixBase< Derived >::UnitZ
inlinestatic
Returns
an expression of the Z axis unit vector (0,0,1{,0}^*)

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
MatrixBase::Unit(Index,Index), MatrixBase::Unit(Index), MatrixBase::UnitY(), MatrixBase::UnitZ(), MatrixBase::UnitW()

Definition at line 979 of file CwiseNullaryOp.h.

980 { return Derived::Unit(2); }

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