template<typename Derived>
MatrixBase class
Base class for all dense matrices, vectors, and expressions.
Template parameters | |
---|---|
Derived | is the derived type, e.g. a matrix type, or an expression, etc. |
Contents
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.
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; }
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
.
Base classes
-
template<typename Derived>class DenseBase
- Base class for all dense matrices, vectors, and arrays.
Derived classes
-
template<typename _Scalar, int _Rows, int _Cols, int _Options, int _MaxRows, int _MaxCols>struct dense_xpr_base_dispatcher<Matrix<_Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols>>
Public types
- using CwiseAbs2ReturnType = CwiseUnaryOp<internal::scalar_abs2_op<Scalar>, const Derived>
- using CwiseAbsReturnType = CwiseUnaryOp<internal::scalar_abs_op<Scalar>, const Derived>
- using CwiseInverseReturnType = CwiseUnaryOp<internal::scalar_inverse_op<Scalar>, const Derived>
- using CwiseScalarEqualReturnType = CwiseBinaryOp<internal::scalar_cmp_op<Scalar, Scalar, internal::cmp_EQ>, const Derived, const ConstantReturnType>
- using CwiseSignReturnType = CwiseUnaryOp<internal::scalar_sign_op<Scalar>, const Derived>
- using CwiseSqrtReturnType = CwiseUnaryOp<internal::scalar_sqrt_op<Scalar>, const Derived>
Public static functions
- static auto Identity() -> const IdentityReturnType
- static auto Identity(Index rows, Index cols) -> const IdentityReturnType
- static auto Unit(Index size, Index i) -> const BasisReturnType
- static auto Unit(Index i) -> const BasisReturnType
- static auto UnitW() -> const BasisReturnType
- static auto UnitX() -> const BasisReturnType
- static auto UnitY() -> const BasisReturnType
- static auto UnitZ() -> const BasisReturnType
Public functions
- auto acosh() const -> const MatrixFunctionReturnValue<Derived>
- This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise inverse hyperbolic cosine use ArrayBase::
acosh . - auto adjoint() const -> const AdjointReturnType
- 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)
- auto array() -> ArrayWrapper<Derived>
- auto array() const -> const ArrayWrapper<const Derived>
- auto asDiagonal() const -> const DiagonalWrapper<const Derived>
- auto asinh() const -> const MatrixFunctionReturnValue<Derived>
- This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise inverse hyperbolic sine use ArrayBase::
asinh . - auto atanh() const -> const MatrixFunctionReturnValue<Derived>
- This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise inverse hyperbolic cosine use ArrayBase::
atanh . - auto bdcSvd(unsigned int computationOptions = 0) const -> BDCSVD<PlainObject>
-
template<typename CustomBinaryOp, typename OtherDerived>auto binaryExpr(const Eigen::
MatrixBase<OtherDerived>& other, const CustomBinaryOp& func = CustomBinaryOp()) const -> const CwiseBinaryOp<CustomBinaryOp, const Derived, const OtherDerived> - auto blueNorm() const -> RealScalar
- auto colPivHouseholderQr() const -> const ColPivHouseholderQR<PlainObject>
- auto completeOrthogonalDecomposition() const -> const CompleteOrthogonalDecomposition<PlainObject>
-
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
- auto cos() const -> const MatrixFunctionReturnValue<Derived>
- This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise cosine use ArrayBase::
cos . - auto cosh() const -> const MatrixFunctionReturnValue<Derived>
- This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise hyperbolic cosine use ArrayBase::
cosh . - auto cwiseAbs() const -> const CwiseAbsReturnType
- auto cwiseAbs2() const -> const CwiseAbs2ReturnType
-
template<typename OtherDerived>auto cwiseEqual(const Eigen::
MatrixBase<OtherDerived>& other) const -> const CwiseBinaryOp<std::equal_to<Scalar>, const Derived, const OtherDerived> - auto cwiseEqual(const Scalar& s) const -> const CwiseScalarEqualReturnType
- auto cwiseInverse() const -> const CwiseInverseReturnType
-
template<typename OtherDerived>auto cwiseMax(const Eigen::
MatrixBase<OtherDerived>& other) const -> const CwiseBinaryOp<internal::scalar_max_op<Scalar, Scalar>, const Derived, const OtherDerived> - auto cwiseMax(const Scalar& other) const -> const CwiseBinaryOp<internal::scalar_max_op<Scalar, Scalar>, const Derived, const ConstantReturnType>
-
template<typename OtherDerived>auto cwiseMin(const Eigen::
MatrixBase<OtherDerived>& other) const -> const CwiseBinaryOp<internal::scalar_min_op<Scalar, Scalar>, const Derived, const OtherDerived> - auto cwiseMin(const Scalar& other) const -> const CwiseBinaryOp<internal::scalar_min_op<Scalar, Scalar>, const Derived, const ConstantReturnType>
-
template<typename OtherDerived>auto cwiseNotEqual(const Eigen::
MatrixBase<OtherDerived>& other) const -> const CwiseBinaryOp<std::not_equal_to<Scalar>, const Derived, const OtherDerived> -
template<typename OtherDerived>auto cwiseProduct(const Eigen::
MatrixBase<OtherDerived>& other) const -> const CwiseBinaryOp<internal::scalar_product_op<Derived ::Scalar, OtherDerived ::Scalar>, const Derived, const OtherDerived> -
template<typename OtherDerived>auto cwiseQuotient(const Eigen::
MatrixBase<OtherDerived>& other) const -> const CwiseBinaryOp<internal::scalar_quotient_op<Scalar>, const Derived, const OtherDerived> - auto cwiseSign() const -> const CwiseSignReturnType
- auto cwiseSqrt() const -> const CwiseSqrtReturnType
- auto determinant() const -> Scalar
- auto diagonal() -> DiagonalReturnType
- auto diagonal() const -> ConstDiagonalReturnType
- auto diagonal(Index index) -> DiagonalDynamicIndexReturnType
- auto diagonal(Index index) const -> ConstDiagonalDynamicIndexReturnType
- auto diagonalSize() const -> Index
-
template<typename OtherDerived>auto dot(const MatrixBase<OtherDerived>& other) const -> ScalarBinaryOpTraits<typename internal::traits<Derived>::Scalar, typename internal::traits<OtherDerived>::Scalar>::ReturnType
- auto eigenvalues() const -> EigenvaluesReturnType
- Computes the eigenvalues of a matrix.
- auto exp() const -> const MatrixExponentialReturnValue<Derived>
- This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise exponential use ArrayBase::
exp . - auto forceAlignedAccess() const -> const Derived&
- auto forceAlignedAccess() -> Derived&
-
template<bool Enable>auto forceAlignedAccessIf() const -> internal::add_const_on_value_type<typename internal::conditional<Enable, ForceAlignedAccess<Derived>, Derived&>::type>::type
-
template<bool Enable>auto forceAlignedAccessIf() -> internal::conditional<Enable, ForceAlignedAccess<Derived>, Derived&>::type
- auto fullPivHouseholderQr() const -> const FullPivHouseholderQR<PlainObject>
- auto fullPivLu() const -> const FullPivLU<PlainObject>
- auto hnormalized() const -> const HNormalizedReturnType
- homogeneous normalization
- auto householderQr() const -> const HouseholderQR<PlainObject>
- auto hypotNorm() const -> RealScalar
- auto inverse() const -> const Inverse<Derived>
- auto isDiagonal(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
- auto isIdentity(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
- auto isLowerTriangular(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
-
template<typename OtherDerived>auto isOrthogonal(const MatrixBase<OtherDerived>& other, const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
- auto isUnitary(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
- auto isUpperTriangular(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const -> bool
- auto jacobiSvd(unsigned int computationOptions = 0) const -> JacobiSVD<PlainObject>
-
template<typename OtherDerived>auto lazyProduct(const MatrixBase<OtherDerived>& other) const -> const Product<Derived, OtherDerived, LazyProduct>
- auto ldlt() const -> const LDLT<PlainObject>
- auto llt() const -> const LLT<PlainObject>
- auto log() const -> const MatrixLogarithmReturnValue<Derived>
- This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise logarithm use ArrayBase::
log . -
template<int p>auto lpNorm() const -> RealScalar
- auto lu() const -> const PartialPivLU<PlainObject>
-
template<typename EssentialPart>void makeHouseholder(EssentialPart& essential, Scalar& tau, RealScalar& beta) const
- void makeHouseholderInPlace(Scalar& tau, RealScalar& beta)
- auto matrixFunction(StemFunction f) const -> const MatrixFunctionReturnValue<Derived>
- Helper function for the unsupported MatrixFunctions module.
-
auto noalias() -> NoAlias<Derived, Eigen::
MatrixBase> - auto norm() const -> RealScalar
- void normalize()
- auto normalized() const -> const PlainObject
-
template<typename OtherDerived>auto operator&&(const Eigen::
MatrixBase<OtherDerived>& other) const -> const CwiseBinaryOp<internal::scalar_boolean_and_op, const Derived, const OtherDerived> -
template<typename T>auto operator*(const T& scalar) const -> const CwiseBinaryOp<internal::scalar_product_op<Scalar, T>, Derived, Constant<T>>
-
template<typename OtherDerived>auto operator*(const MatrixBase<OtherDerived>& other) const -> const Product<Derived, OtherDerived>
-
template<typename DiagonalDerived>auto operator*(const DiagonalBase<DiagonalDerived>& diagonal) const -> const Product<Derived, DiagonalDerived, LazyProduct>
-
template<typename OtherDerived>auto operator*=(const EigenBase<OtherDerived>& other) -> Derived&
-
template<typename OtherDerived>auto operator!=(const MatrixBase<OtherDerived>& other) const -> bool
-
template<typename OtherDerived>auto operator+(const Eigen::
MatrixBase<OtherDerived>& other) const -> const CwiseBinaryOp<sum<Scalar>, const Derived, const OtherDerived> -
template<typename OtherDerived>auto operator+=(const MatrixBase<OtherDerived>& other) -> Derived&
-
template<typename OtherDerived>auto operator-(const Eigen::
MatrixBase<OtherDerived>& other) const -> const CwiseBinaryOp<difference<Scalar>, const Derived, const OtherDerived> -
template<typename OtherDerived>auto operator-=(const MatrixBase<OtherDerived>& other) -> Derived&
-
template<typename T>auto operator/(const T& scalar) const -> const CwiseBinaryOp<internal::scalar_quotient_op<Scalar, T>, Derived, Constant<T>>
- auto operator=(const MatrixBase& other) -> Derived&
-
template<typename OtherDerived>auto operator==(const MatrixBase<OtherDerived>& other) const -> bool
- auto operatorNorm() const -> RealScalar
- Computes the L2 operator norm.
-
template<typename OtherDerived>auto operator||(const Eigen::
MatrixBase<OtherDerived>& other) const -> const CwiseBinaryOp<internal::scalar_boolean_or_op, const Derived, const OtherDerived> - auto partialPivLu() const -> const PartialPivLU<PlainObject>
- auto pow(const RealScalar& p) const -> const MatrixPowerReturnValue<Derived>
- This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise power to
p
use ArrayBase::pow . - auto pow(const std::complex<RealScalar>& p) const -> const MatrixComplexPowerReturnValue<Derived>
- This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise power to
p
use ArrayBase::pow . -
template<unsigned int UpLo>auto selfadjointView() const -> MatrixBase<Derived>::template ConstSelfAdjointViewReturnType<UpLo>::Type
-
template<unsigned int UpLo>auto selfadjointView() -> MatrixBase<Derived>::template SelfAdjointViewReturnType<UpLo>::Type
- auto setIdentity() -> Derived&
- auto setIdentity(Index rows, Index cols) -> Derived&
- Resizes to the given size, and writes the identity expression (not necessarily square) into *this.
- auto setUnit(Index i) -> Derived&
- Set the coefficients of
*this
to the i-th unit (basis) vector. - auto setUnit(Index newSize, Index i) -> Derived&
- Resizes to the given newSize, and writes the i-th unit (basis) vector into *this.
- auto sin() const -> const MatrixFunctionReturnValue<Derived>
- This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise sine use ArrayBase::
sin . - auto sinh() const -> const MatrixFunctionReturnValue<Derived>
- This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise hyperbolic sine use ArrayBase::
sinh . - auto sqrt() const -> const MatrixSquareRootReturnValue<Derived>
- This function requires the unsupported MatrixFunctions module. To compute the coefficient-wise square root use ArrayBase::
sqrt . - auto squaredNorm() const -> RealScalar
- auto stableNorm() const -> RealScalar
- void stableNormalize()
- auto stableNormalized() const -> const PlainObject
- auto trace() const -> Scalar
-
template<unsigned int Mode>auto triangularView() -> MatrixBase<Derived>::template TriangularViewReturnType<Mode>::Type
-
template<unsigned int Mode>auto triangularView() const -> MatrixBase<Derived>::template ConstTriangularViewReturnType<Mode>::Type
Friends
-
template<typename T>auto operator*(const T& scalar, const StorageBaseType& expr) -> const CwiseBinaryOp<internal::scalar_product_op<T, Scalar>, Constant<T>, Derived>
Typedef documentation
template<typename Derived>
typedef CwiseUnaryOp<internal::scalar_abs2_op<Scalar>, const Derived> Eigen:: MatrixBase<Derived>:: CwiseAbs2ReturnType
template<typename Derived>
typedef CwiseUnaryOp<internal::scalar_abs_op<Scalar>, const Derived> Eigen:: MatrixBase<Derived>:: CwiseAbsReturnType
template<typename Derived>
typedef CwiseUnaryOp<internal::scalar_inverse_op<Scalar>, const Derived> Eigen:: MatrixBase<Derived>:: CwiseInverseReturnType
template<typename Derived>
typedef CwiseBinaryOp<internal::scalar_cmp_op<Scalar, Scalar, internal::cmp_EQ>, const Derived, const ConstantReturnType> Eigen:: MatrixBase<Derived>:: CwiseScalarEqualReturnType
template<typename Derived>
typedef CwiseUnaryOp<internal::scalar_sign_op<Scalar>, const Derived> Eigen:: MatrixBase<Derived>:: CwiseSignReturnType
template<typename Derived>
typedef CwiseUnaryOp<internal::scalar_sqrt_op<Scalar>, const Derived> Eigen:: MatrixBase<Derived>:: CwiseSqrtReturnType
Function documentation
template<typename Derived>
static const IdentityReturnType Eigen:: MatrixBase<Derived>:: Identity()
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
template<typename Derived>
static const IdentityReturnType Eigen:: MatrixBase<Derived>:: Identity(Index rows,
Index cols)
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;
Output:
1 0 0 0 1 0 0 0 1 0 0 0
template<typename Derived>
static const BasisReturnType Eigen:: MatrixBase<Derived>:: Unit(Index size,
Index i)
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.
template<typename Derived>
static const BasisReturnType Eigen:: MatrixBase<Derived>:: Unit(Index i)
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.
template<typename Derived>
static const BasisReturnType Eigen:: MatrixBase<Derived>:: UnitW()
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.
template<typename Derived>
static const BasisReturnType Eigen:: MatrixBase<Derived>:: UnitX()
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.
template<typename Derived>
static const BasisReturnType Eigen:: MatrixBase<Derived>:: UnitY()
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.
template<typename Derived>
static const BasisReturnType Eigen:: MatrixBase<Derived>:: UnitZ()
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.
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::
Returns | an expression of the matrix inverse hyperbolic cosine of *this . |
---|
template<typename Derived>
const AdjointReturnType Eigen:: MatrixBase<Derived>:: adjoint() const
Returns | an expression of the adjoint (i.e. conjugate transpose) of *this. |
---|
Example:
Matrix2cf m = Matrix2cf::Random(); cout << "Here is the 2x2 complex matrix m:" << endl << m << endl; cout << "Here is the adjoint of m:" << endl << m.adjoint() << endl;
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)
template<typename Derived>
void Eigen:: MatrixBase<Derived>:: adjointInPlace()
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().
template<typename Derived>
template<typename EssentialPart>
void Eigen:: MatrixBase<Derived>:: applyHouseholderOnTheLeft(const EssentialPart& essential,
const Scalar& tau,
Scalar* workspace)
Parameters | |
---|---|
essential | the essential part of the vector v |
tau | the scaling factor of the Householder transformation |
workspace | a pointer to working space with at least this->cols() entries |
Apply the elementary reflector H given by with from the left to a vector or matrix.
On input:
template<typename Derived>
template<typename EssentialPart>
void Eigen:: MatrixBase<Derived>:: applyHouseholderOnTheRight(const EssentialPart& essential,
const Scalar& tau,
Scalar* workspace)
Parameters | |
---|---|
essential | the essential part of the vector v |
tau | the scaling factor of the Householder transformation |
workspace | a pointer to working space with at least this->rows() entries |
Apply the elementary reflector H given by with from the right to a vector or matrix.
On input:
template<typename Derived>
template<typename OtherDerived>
void Eigen:: MatrixBase<Derived>:: applyOnTheLeft(const EigenBase<OtherDerived>& other)
replaces *this
by other * *this
.
Example:
Matrix3f A = Matrix3f::Random(3,3), B; 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;
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
template<typename Derived>
template<typename OtherScalar>
void Eigen:: MatrixBase<Derived>:: applyOnTheLeft(Index p,
Index q,
const JacobiRotation<OtherScalar>& j)
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 .
template<typename Derived>
template<typename OtherDerived>
void Eigen:: MatrixBase<Derived>:: applyOnTheRight(const EigenBase<OtherDerived>& other)
replaces *this
by *this
* other. It is equivalent to MatrixBase::
Example:
Matrix3f A = Matrix3f::Random(3,3), B; 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
template<typename Derived>
ArrayWrapper<Derived> Eigen:: MatrixBase<Derived>:: array()
Returns | an Array expression of this matrix |
---|
template<typename Derived>
const ArrayWrapper<const Derived> Eigen:: MatrixBase<Derived>:: array() const
Returns | a const Array expression of this matrix |
---|
template<typename Derived>
const DiagonalWrapper<const Derived> Eigen:: MatrixBase<Derived>:: asDiagonal() const
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;
Output:
2 0 0 0 5 0 0 0 6
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::
Returns | an expression of the matrix inverse hyperbolic sine of *this . |
---|
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::
Returns | an expression of the matrix inverse hyperbolic cosine of *this . |
---|
template<typename Derived>
BDCSVD<PlainObject> Eigen:: MatrixBase<Derived>:: bdcSvd(unsigned int computationOptions = 0) const
Returns | the singular value decomposition of *this computed by Divide & Conquer algorithm |
---|
This is defined in the SVD module. #include <Eigen/SVD>
template<typename Derived>
template<typename CustomBinaryOp, typename OtherDerived>
const CwiseBinaryOp<CustomBinaryOp, const Derived, const OtherDerived> Eigen:: MatrixBase<Derived>:: binaryExpr(const Eigen:: MatrixBase<OtherDerived>& other,
const CustomBinaryOp& func = CustomBinaryOp()) const
Returns | an expression of a custom coefficient-wise operator func of *this and other |
---|
The template parameter CustomBinaryOp is the type of the functor of the custom operator (see class CwiseBinaryOp for an example)
Here is an example illustrating the use of custom functors:
#include <Eigen/Core> #include <iostream> using namespace Eigen; using namespace std; // define a custom template binary functor template<typename Scalar> struct MakeComplexOp { EIGEN_EMPTY_STRUCT_CTOR(MakeComplexOp) typedef complex<Scalar> result_type; complex<Scalar> operator()(const Scalar& a, const Scalar& b) const { return complex<Scalar>(a,b); } }; int main(int, char**) { Matrix4d m1 = Matrix4d::Random(), m2 = Matrix4d::Random(); cout << m1.binaryExpr(m2, MakeComplexOp<double>()) << endl; return 0; }
Output:
(0.68,0.271) (0.823,-0.967) (-0.444,-0.687) (-0.27,0.998) (-0.211,0.435) (-0.605,-0.514) (0.108,-0.198) (0.0268,-0.563) (0.566,-0.717) (-0.33,-0.726) (-0.0452,-0.74) (0.904,0.0259) (0.597,0.214) (0.536,0.608) (0.258,-0.782) (0.832,0.678)
template<typename Derived>
RealScalar Eigen:: MatrixBase<Derived>:: blueNorm() const
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.
template<typename Derived>
const ColPivHouseholderQR<PlainObject> Eigen:: MatrixBase<Derived>:: colPivHouseholderQr() const
Returns | the column-pivoting Householder QR decomposition of *this . |
---|
template<typename Derived>
const CompleteOrthogonalDecomposition<PlainObject> Eigen:: MatrixBase<Derived>:: completeOrthogonalDecomposition() const
Returns | the complete orthogonal decomposition of *this . |
---|
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
Parameters | |
---|---|
inverse | Reference to the matrix in which to store the inverse. |
determinant | Reference to the variable in which to store the determinant. |
invertible | Reference to the bool variable in which to store whether the matrix is invertible. |
absDeterminantThreshold | Optional parameter controlling the invertibility check. The matrix will be declared invertible if the absolute value of its determinant is greater than this threshold. |
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.
Example:
Matrix3d m = Matrix3d::Random(); cout << "Here is the matrix m:" << endl << m << endl; Matrix3d inverse; 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; }
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
template<typename Derived>
template<typename ResultType>
void Eigen:: MatrixBase<Derived>:: computeInverseWithCheck(ResultType& inverse,
bool& invertible,
const RealScalar& absDeterminantThreshold = NumTraits<Scalar>::dummy_precision()) const
Parameters | |
---|---|
inverse | Reference to the matrix in which to store the inverse. |
invertible | Reference to the bool variable in which to store whether the matrix is invertible. |
absDeterminantThreshold | Optional parameter controlling the invertibility check. The matrix will be declared invertible if the absolute value of its determinant is greater than this threshold. |
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.
Example:
Matrix3d m = Matrix3d::Random(); cout << "Here is the matrix m:" << endl << m << endl; Matrix3d inverse; 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
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::
Returns | an expression of the matrix cosine of *this . |
---|
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::
Returns | an expression of the matrix hyperbolic cosine of *this . |
---|
template<typename Derived>
const CwiseAbsReturnType Eigen:: MatrixBase<Derived>:: cwiseAbs() const
Returns | an expression of the coefficient-wise absolute value of *this |
---|
Example:
MatrixXd m(2,3); m << 2, -4, 6, -5, 1, 0; cout << m.cwiseAbs() << endl;
Output:
template<typename Derived>
const CwiseAbs2ReturnType Eigen:: MatrixBase<Derived>:: cwiseAbs2() const
Returns | an expression of the coefficient-wise squared absolute value of *this |
---|
Example:
MatrixXd m(2,3); m << 2, -4, 6, -5, 1, 0; cout << m.cwiseAbs2() << endl;
Output:
template<typename Derived>
template<typename OtherDerived>
const CwiseBinaryOp<std::equal_to<Scalar>, const Derived, const OtherDerived> Eigen:: MatrixBase<Derived>:: cwiseEqual(const Eigen:: MatrixBase<OtherDerived>& other) const
Returns | an expression of the coefficient-wise == operator of *this and other |
---|
Example:
MatrixXi m(2,2); m << 1, 0, 1, 1; cout << "Comparing m with identity matrix:" << endl; cout << m.cwiseEqual(MatrixXi::Identity(2,2)) << endl; Index count = m.cwiseEqual(MatrixXi::Identity(2,2)).count(); cout << "Number of coefficients that are equal: " << count << endl;
Output:
Comparing m with identity matrix: 1 1 0 1 Number of coefficients that are equal: 3
template<typename Derived>
const CwiseScalarEqualReturnType Eigen:: MatrixBase<Derived>:: cwiseEqual(const Scalar& s) const
Returns | an expression of the coefficient-wise == operator of *this and a scalar s |
---|
template<typename Derived>
const CwiseInverseReturnType Eigen:: MatrixBase<Derived>:: cwiseInverse() const
Returns | an expression of the coefficient-wise inverse of *this. |
---|
Example:
MatrixXd m(2,3); m << 2, 0.5, 1, 3, 0.25, 1; cout << m.cwiseInverse() << endl;
Output:
template<typename Derived>
template<typename OtherDerived>
const CwiseBinaryOp<internal::scalar_max_op<Scalar, Scalar>, const Derived, const OtherDerived> Eigen:: MatrixBase<Derived>:: cwiseMax(const Eigen:: MatrixBase<OtherDerived>& other) const
Returns | an expression of the coefficient-wise max of *this and other |
---|
Example:
Vector3d v(2,3,4), w(4,2,3); cout << v.cwiseMax(w) << endl;
Output:
4 3 4
template<typename Derived>
const CwiseBinaryOp<internal::scalar_max_op<Scalar, Scalar>, const Derived, const ConstantReturnType> Eigen:: MatrixBase<Derived>:: cwiseMax(const Scalar& other) const
Returns | an expression of the coefficient-wise max of *this and scalar other |
---|
template<typename Derived>
template<typename OtherDerived>
const CwiseBinaryOp<internal::scalar_min_op<Scalar, Scalar>, const Derived, const OtherDerived> Eigen:: MatrixBase<Derived>:: cwiseMin(const Eigen:: MatrixBase<OtherDerived>& other) const
Returns | an expression of the coefficient-wise min of *this and other |
---|
Example:
Vector3d v(2,3,4), w(4,2,3); cout << v.cwiseMin(w) << endl;
Output:
2 2 3
template<typename Derived>
const CwiseBinaryOp<internal::scalar_min_op<Scalar, Scalar>, const Derived, const ConstantReturnType> Eigen:: MatrixBase<Derived>:: cwiseMin(const Scalar& other) const
Returns | an expression of the coefficient-wise min of *this and scalar other |
---|
template<typename Derived>
template<typename OtherDerived>
const CwiseBinaryOp<std::not_equal_to<Scalar>, const Derived, const OtherDerived> Eigen:: MatrixBase<Derived>:: cwiseNotEqual(const Eigen:: MatrixBase<OtherDerived>& other) const
Returns | an expression of the coefficient-wise != operator of *this and other |
---|
Example:
MatrixXi m(2,2); m << 1, 0, 1, 1; cout << "Comparing m with identity matrix:" << endl; cout << m.cwiseNotEqual(MatrixXi::Identity(2,2)) << endl; Index count = m.cwiseNotEqual(MatrixXi::Identity(2,2)).count(); cout << "Number of coefficients that are not equal: " << count << endl;
Output:
Comparing m with identity matrix: 0 0 1 0 Number of coefficients that are not equal: 1
template<typename Derived>
template<typename OtherDerived>
const CwiseBinaryOp<internal::scalar_product_op<Derived ::Scalar, OtherDerived ::Scalar>, const Derived, const OtherDerived> Eigen:: MatrixBase<Derived>:: cwiseProduct(const Eigen:: MatrixBase<OtherDerived>& other) const
Returns | an expression of the Schur product (coefficient wise product) of *this and other |
---|
Example:
Matrix3i a = Matrix3i::Random(), b = Matrix3i::Random(); Matrix3i c = a.cwiseProduct(b); cout << "a:\n" << a << "\nb:\n" << b << "\nc:\n" << c << endl;
Output:
a: 7 6 -3 -2 9 6 6 -6 -5 b: 1 -3 9 0 0 3 3 9 5 c: 7 -18 -27 0 0 18 18 -54 -25
template<typename Derived>
template<typename OtherDerived>
const CwiseBinaryOp<internal::scalar_quotient_op<Scalar>, const Derived, const OtherDerived> Eigen:: MatrixBase<Derived>:: cwiseQuotient(const Eigen:: MatrixBase<OtherDerived>& other) const
Returns | an expression of the coefficient-wise quotient of *this and other |
---|
Example:
Vector3d v(2,3,4), w(4,2,3); cout << v.cwiseQuotient(w) << endl;
Output:
0.5 1.5 1.33
template<typename Derived>
const CwiseSignReturnType Eigen:: MatrixBase<Derived>:: cwiseSign() const
Returns | an expression of the coefficient-wise signum of *this. |
---|
Example:
MatrixXd m(2,3); m << 2, -4, 6, -5, 1, 0; cout << m.cwiseSign() << endl;
Output:
1 -1 1 -1 1 0
template<typename Derived>
const CwiseSqrtReturnType Eigen:: MatrixBase<Derived>:: cwiseSqrt() const
Returns | an expression of the coefficient-wise square root of *this. |
---|
Example:
Vector3d v(1,2,4); cout << v.cwiseSqrt() << endl;
Output:
cwisePow(), cwiseSquare()
template<typename Derived>
Scalar Eigen:: MatrixBase<Derived>:: determinant() const
Returns | the determinant of this matrix |
---|
This is defined in the LU module. #include <Eigen/LU>
template<typename Derived>
DiagonalReturnType Eigen:: MatrixBase<Derived>:: diagonal()
Returns | an expression of the main diagonal of the matrix *this |
---|
*this
is not required to be square.
Example:
Matrix3i m = Matrix3i::Random(); 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
*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:
Matrix4i m = Matrix4i::Random(); 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
template<typename Derived>
ConstDiagonalReturnType Eigen:: MatrixBase<Derived>:: diagonal() const
This is the const version of diagonal().
This is the const version of diagonal<int>().
template<typename Derived>
DiagonalDynamicIndexReturnType Eigen:: MatrixBase<Derived>:: diagonal(Index index)
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:
Matrix4i m = Matrix4i::Random(); 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
template<typename Derived>
ConstDiagonalDynamicIndexReturnType Eigen:: MatrixBase<Derived>:: diagonal(Index index) const
This is the const version of diagonal(Index).
template<typename Derived>
Index Eigen:: MatrixBase<Derived>:: diagonalSize() const
Returns | the size of the main diagonal, which is min(rows(),cols()). |
---|
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
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.
template<typename Derived>
EigenvaluesReturnType Eigen:: MatrixBase<Derived>:: eigenvalues() const
Computes the eigenvalues of a matrix.
Returns | Column vector containing the eigenvalues. |
---|
This is defined in the Eigenvalues module. #include <Eigen/Eigenvalues>
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:
MatrixXd ones = MatrixXd::Ones(3,3); VectorXcd eivals = ones.eigenvalues(); cout << "The eigenvalues of the 3x3 matrix of ones are:" << endl << eivals << endl;
Output:
The eigenvalues of the 3x3 matrix of ones are: (-5.31e-17,0) (3,0) (0,0)
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::
Returns | an expression of the matrix exponential of *this . |
---|
template<typename Derived>
const Derived& Eigen:: MatrixBase<Derived>:: forceAlignedAccess() const
Returns | an expression of *this with forced aligned access |
---|
template<typename Derived>
Derived& Eigen:: MatrixBase<Derived>:: forceAlignedAccess()
Returns | an expression of *this with forced aligned access |
---|
template<typename Derived>
template<bool Enable>
internal::add_const_on_value_type<typename internal::conditional<Enable, ForceAlignedAccess<Derived>, Derived&>::type>::type Eigen:: MatrixBase<Derived>:: forceAlignedAccessIf() const
Returns | an expression of *this with forced aligned access if Enable is true. |
---|
template<typename Derived>
template<bool Enable>
internal::conditional<Enable, ForceAlignedAccess<Derived>, Derived&>::type Eigen:: MatrixBase<Derived>:: forceAlignedAccessIf()
Returns | an expression of *this with forced aligned access if Enable is true. |
---|
template<typename Derived>
const FullPivHouseholderQR<PlainObject> Eigen:: MatrixBase<Derived>:: fullPivHouseholderQr() const
Returns | the full-pivoting Householder QR decomposition of *this . |
---|
template<typename Derived>
const FullPivLU<PlainObject> Eigen:: MatrixBase<Derived>:: fullPivLu() const
Returns | the full-pivoting LU decomposition of *this . |
---|
This is defined in the LU module. #include <Eigen/LU>
template<typename Derived>
const HouseholderQR<PlainObject> Eigen:: MatrixBase<Derived>:: householderQr() const
Returns | the Householder QR decomposition of *this . |
---|
template<typename Derived>
RealScalar Eigen:: MatrixBase<Derived>:: hypotNorm() const
Returns | the l2 norm of *this avoiding undeflow and overflow. This version use a concatenation of hypot() calls, and it is very slow. |
---|
template<typename Derived>
const Inverse<Derived> Eigen:: MatrixBase<Derived>:: inverse() const
Returns | the matrix inverse of this matrix. |
---|
This is defined in the LU module. #include <Eigen/LU>
For small fixed sizes up to 4x4, this method uses cofactors. In the general case, this method uses class PartialPivLU.
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:
Matrix3d m = 10000 * Matrix3d::Identity(); 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;
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
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:
Matrix3d m = Matrix3d::Identity(); 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
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. |
---|
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;
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
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. |
---|
Example:
Matrix3d m = Matrix3d::Identity(); 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
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. |
---|
template<typename Derived>
JacobiSVD<PlainObject> Eigen:: MatrixBase<Derived>:: jacobiSvd(unsigned int computationOptions = 0) const
Returns | the singular value decomposition of *this computed by two-sided Jacobi transformations. |
---|
This is defined in the SVD module. #include <Eigen/SVD>
template<typename Derived>
template<typename OtherDerived>
const Product<Derived, OtherDerived, LazyProduct> Eigen:: MatrixBase<Derived>:: lazyProduct(const MatrixBase<OtherDerived>& other) const
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.
template<typename Derived>
const LDLT<PlainObject> Eigen:: MatrixBase<Derived>:: ldlt() const
Returns | the Cholesky decomposition with full pivoting without square root of *this |
---|
This is defined in the Cholesky module. #include <Eigen/Cholesky>
template<typename Derived>
const LLT<PlainObject> Eigen:: MatrixBase<Derived>:: llt() const
Returns | the LLT decomposition of *this |
---|
This is defined in the Cholesky module. #include <Eigen/Cholesky>
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::
Returns | an expression of the matrix logarithm of *this . |
---|
template<typename Derived>
template<int p>
RealScalar Eigen:: MatrixBase<Derived>:: lpNorm() const
Returns | the coefficient-wise 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::*this . |
---|
In all cases, if *this
is empty, then the value 0 is returned.
template<typename Derived>
const PartialPivLU<PlainObject> Eigen:: MatrixBase<Derived>:: lu() const
Returns | the partial-pivoting LU decomposition of *this . |
---|
This is defined in the LU module. #include <Eigen/LU>
Synonym of partialPivLu().
template<typename Derived>
template<typename EssentialPart>
void Eigen:: MatrixBase<Derived>:: makeHouseholder(EssentialPart& essential,
Scalar& tau,
RealScalar& beta) const
Parameters | |
---|---|
essential | the essential part of the vector v |
tau | the scaling factor of the Householder transformation |
beta | the result of H * *this |
Computes the elementary reflector H such that: where the transformation H is: and the vector v is:
On output:
template<typename Derived>
void Eigen:: MatrixBase<Derived>:: makeHouseholderInPlace(Scalar& tau,
RealScalar& beta)
Parameters | |
---|---|
tau | the scaling factor of the Householder transformation |
beta | the result of H * *this |
Computes the elementary reflector H such that: where the transformation H is: and the vector v is:
The essential part of the vector v
is stored in *this.
On output:
template<typename Derived>
NoAlias<Derived, Eigen:: 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;
template<typename Derived>
RealScalar Eigen:: MatrixBase<Derived>:: norm() const
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. |
---|
template<typename Derived>
void Eigen:: MatrixBase<Derived>:: normalize()
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.
template<typename Derived>
const PlainObject Eigen:: MatrixBase<Derived>:: normalized() const
Returns | an expression of the quotient of *this by its own norm. |
---|
This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.
template<typename Derived>
template<typename OtherDerived>
const CwiseBinaryOp<internal::scalar_boolean_and_op, const Derived, const OtherDerived> Eigen:: MatrixBase<Derived>:: operator&&(const Eigen:: MatrixBase<OtherDerived>& other) const
Returns | an expression of the coefficient-wise boolean and operator of *this and other |
---|
Example:
Array3d v(-1,2,1), w(-3,2,3); cout << ((v<w) && (v<0)) << endl;
Output:
0 0 0
template<typename Derived>
template<typename T>
const CwiseBinaryOp<internal::scalar_product_op<Scalar, T>, Derived, Constant<T>> Eigen:: MatrixBase<Derived>:: operator*(const T& scalar) const
Template parameters | |
---|---|
T | is the scalar type of scalar. It must be compatible with the scalar type of the given expression. |
Returns | an expression of *this scaled by the scalar factor scalar |
template<typename Derived>
template<typename OtherDerived>
const Product<Derived, OtherDerived> Eigen:: MatrixBase<Derived>:: operator*(const MatrixBase<OtherDerived>& other) const
Returns | the matrix product of *this and other. |
---|
template<typename Derived>
template<typename OtherDerived>
Derived& Eigen:: MatrixBase<Derived>:: operator*=(const EigenBase<OtherDerived>& other)
Returns | a reference to *this |
---|
replaces *this
by *this
* other.
Example:
Matrix3f A = Matrix3f::Random(3,3), B; 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
template<typename Derived>
template<typename OtherDerived>
bool Eigen:: MatrixBase<Derived>:: operator!=(const MatrixBase<OtherDerived>& other) const
Returns | true if at least one pair of coefficients of *this and other are not exactly equal to each other. |
---|
template<typename Derived>
template<typename OtherDerived>
const CwiseBinaryOp<sum<Scalar>, const Derived, const OtherDerived> Eigen:: MatrixBase<Derived>:: operator+(const Eigen:: MatrixBase<OtherDerived>& other) const
Returns | an expression of the sum of *this and other |
---|
template<typename Derived>
template<typename OtherDerived>
Derived& Eigen:: MatrixBase<Derived>:: operator+=(const MatrixBase<OtherDerived>& other)
Returns | a reference to *this |
---|
replaces *this
by *this
+ other.
template<typename Derived>
template<typename OtherDerived>
const CwiseBinaryOp<difference<Scalar>, const Derived, const OtherDerived> Eigen:: MatrixBase<Derived>:: operator-(const Eigen:: MatrixBase<OtherDerived>& other) const
Returns | an expression of the difference of *this and other |
---|
template<typename Derived>
template<typename OtherDerived>
Derived& Eigen:: MatrixBase<Derived>:: operator-=(const MatrixBase<OtherDerived>& other)
Returns | a reference to *this |
---|
replaces *this
by *this
- other.
template<typename Derived>
template<typename T>
const CwiseBinaryOp<internal::scalar_quotient_op<Scalar, T>, Derived, Constant<T>> Eigen:: MatrixBase<Derived>:: operator/(const T& scalar) const
Template parameters | |
---|---|
T | is the scalar type of scalar. It must be compatible with the scalar type of the given expression. |
Returns | an expression of *this divided by the scalar value scalar |
template<typename Derived>
Derived& Eigen:: MatrixBase<Derived>:: operator=(const MatrixBase& other)
Special case of the template operator=, in order to prevent the compiler from generating a default operator= (issue hit with g++ 4.1)
template<typename Derived>
template<typename OtherDerived>
bool Eigen:: MatrixBase<Derived>:: operator==(const MatrixBase<OtherDerived>& other) const
Returns | true if each coefficients of *this and other are all exactly equal. |
---|
template<typename Derived>
RealScalar Eigen:: MatrixBase<Derived>:: operatorNorm() const
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 is defined to be
where the maximum is over all vectors and the norm on the right is the Euclidean vector norm. The norm equals the largest singular value, which is the square root of the largest eigenvalue of the positive semi-definite matrix .
The current implementation uses the eigenvalues of , as computed by SelfAdjointView::
Example:
MatrixXd ones = MatrixXd::Ones(3,3); 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
template<typename Derived>
template<typename OtherDerived>
const CwiseBinaryOp<internal::scalar_boolean_or_op, const Derived, const OtherDerived> Eigen:: MatrixBase<Derived>:: operator||(const Eigen:: MatrixBase<OtherDerived>& other) const
Returns | an expression of the coefficient-wise boolean or operator of *this and other |
---|
Example:
Array3d v(-1,2,1), w(-3,2,3); cout << ((v<w) || (v<0)) << endl;
Output:
1 0 1
template<typename Derived>
const PartialPivLU<PlainObject> Eigen:: MatrixBase<Derived>:: partialPivLu() const
Returns | the partial-pivoting LU decomposition of *this . |
---|
This is defined in the LU module. #include <Eigen/LU>
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::
Returns | an expression of the matrix power to p of *this . |
---|
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::
Returns | an expression of the matrix power to p of *this . |
---|
template<typename Derived>
template<unsigned int UpLo>
MatrixBase<Derived>::template ConstSelfAdjointViewReturnType<UpLo>::Type Eigen:: MatrixBase<Derived>:: selfadjointView() const
This is the const version of MatrixBase::selfadjointView()
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:
Matrix3i m = Matrix3i::Random(); 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;
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
template<typename Derived>
Derived& Eigen:: MatrixBase<Derived>:: setIdentity()
Writes the identity expression (not necessarily square) into *this.
Example:
Matrix4i m = Matrix4i::Zero(); m.block<3,3>(1,0).setIdentity(); cout << m << endl;
Output:
0 0 0 0 1 0 0 0 0 1 0 0 0 0 1 0
template<typename Derived>
Derived& Eigen:: MatrixBase<Derived>:: setIdentity(Index rows,
Index cols)
Resizes to the given size, and writes the identity expression (not necessarily square) into *this.
Parameters | |
---|---|
rows | the new number of rows |
cols | the new number of columns |
Example:
MatrixXf m; m.setIdentity(3, 3); cout << m << endl;
Output:
1 0 0 0 1 0 0 0 1
template<typename Derived>
Derived& Eigen:: MatrixBase<Derived>:: setUnit(Index i)
Set the coefficients of *this
to the i-th unit (basis) vector.
Parameters | |
---|---|
i | index of the unique coefficient to be set to 1 |
This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.
template<typename Derived>
Derived& Eigen:: MatrixBase<Derived>:: setUnit(Index newSize,
Index i)
Resizes to the given newSize, and writes the i-th unit (basis) vector into *this.
Parameters | |
---|---|
newSize | the new size of the vector |
i | index of the unique coefficient to be set to 1 |
This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.
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::
Returns | an expression of the matrix sine of *this . |
---|
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::
Returns | an expression of the matrix hyperbolic sine of *this . |
---|
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::
Returns | an expression of the matrix square root of *this . |
---|
template<typename Derived>
RealScalar Eigen:: MatrixBase<Derived>:: squaredNorm() const
Returns | , for vectors, the squared l2 norm of *this , and for matrices the 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. |
---|
template<typename Derived>
RealScalar Eigen:: MatrixBase<Derived>:: stableNorm() const
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 in a standard way |
---|
For architecture/scalar types supporting vectorization, this version is faster than blueNorm(). Otherwise the blueNorm() is much faster.
template<typename Derived>
void Eigen:: MatrixBase<Derived>:: stableNormalize()
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.
template<typename Derived>
const PlainObject Eigen:: MatrixBase<Derived>:: stableNormalized() const
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.
template<typename Derived>
template<unsigned int Mode>
MatrixBase<Derived>::template TriangularViewReturnType<Mode>::Type Eigen:: MatrixBase<Derived>:: triangularView()
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:
Matrix3i m = Matrix3i::Random(); 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)
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
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()
template<typename Derived>
template<typename T>
const CwiseBinaryOp<internal::scalar_product_op<T, Scalar>, Constant<T>, Derived> operator*(const T& scalar,
const StorageBaseType& expr)
Template parameters | |
---|---|
T | is the scalar type of scalar. It must be compatible with the scalar type of the given expression. |
Returns | an expression of expr scaled by the scalar factor scalar |