template<typename _MatrixType, unsigned int UpLo>
Eigen::SelfAdjointView class

Expression of a selfadjoint matrix from a triangular part of a dense matrix.

This class is an expression of a sefladjoint matrix from a triangular part of a matrix with given dense storage of the coefficients. It is the return type of MatrixBase::selfadjointView() and most of the time this is the only way that it is used.

Base classes

template<typename Derived>
class TriangularBase
Base class for triangular part in a matrix.

Public types

using EigenvaluesReturnType = Matrix<RealScalar, internal::traits<MatrixType>::ColsAtCompileTime, 1>
using RealScalar = NumTraits<Scalar>::Real
using Scalar = internal::traits<SelfAdjointView>::Scalar
The type of coefficients in this matrix.

Public functions

auto adjoint() const -> const AdjointReturnType
auto coeff(Index row, Index col) const -> Scalar
auto coeffRef(Index row, Index col) -> Scalar&
auto conjugate() const -> const ConjugateReturnType
template<bool Cond>
auto conjugateIf() const -> internal::conditional<Cond, ConjugateReturnType, ConstSelfAdjointView>::type
auto diagonal() const -> MatrixType::ConstDiagonalReturnType
auto eigenvalues() const -> EigenvaluesReturnType
Computes the eigenvalues of a matrix.
auto ldlt() const -> const LDLT<PlainObject, UpLo>
auto llt() const -> const LLT<PlainObject, UpLo>
template<typename OtherDerived>
auto operator*(const MatrixBase<OtherDerived>& rhs) const -> const Product<SelfAdjointView, OtherDerived>
auto operatorNorm() const -> RealScalar
Computes the L2 operator norm.
template<typename DerivedU, typename DerivedV>
auto rankUpdate(const MatrixBase<DerivedU>& u, const MatrixBase<DerivedV>& v, const Scalar& alpha = Scalar(1)) -> SelfAdjointView&
template<typename DerivedU>
auto rankUpdate(const MatrixBase<DerivedU>& u, const Scalar& alpha = Scalar(1)) -> SelfAdjointView&
auto transpose() -> TransposeReturnType
auto transpose() const -> const ConstTransposeReturnType
template<unsigned int TriMode>
auto triangularView() const -> internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), TriangularView<MatrixType, TriMode>, TriangularView<typename MatrixType::AdjointReturnType, TriMode>>::type

Friends

template<typename OtherDerived>
auto operator*(const MatrixBase<OtherDerived>& lhs, const SelfAdjointView& rhs) -> const Product<OtherDerived, SelfAdjointView>

Typedef documentation

template<typename _MatrixType, unsigned int UpLo>
typedef Matrix<RealScalar, internal::traits<MatrixType>::ColsAtCompileTime, 1> Eigen::SelfAdjointView<_MatrixType, UpLo>::EigenvaluesReturnType

Return type of eigenvalues()

template<typename _MatrixType, unsigned int UpLo>
typedef NumTraits<Scalar>::Real Eigen::SelfAdjointView<_MatrixType, UpLo>::RealScalar

Real part of Scalar

Function documentation

template<typename _MatrixType, unsigned int UpLo>
const AdjointReturnType Eigen::SelfAdjointView<_MatrixType, UpLo>::adjoint() const

template<typename _MatrixType, unsigned int UpLo>
Scalar Eigen::SelfAdjointView<_MatrixType, UpLo>::coeff(Index row, Index col) const

template<typename _MatrixType, unsigned int UpLo>
Scalar& Eigen::SelfAdjointView<_MatrixType, UpLo>::coeffRef(Index row, Index col)

template<typename _MatrixType, unsigned int UpLo>
const ConjugateReturnType Eigen::SelfAdjointView<_MatrixType, UpLo>::conjugate() const

template<typename _MatrixType, unsigned int UpLo> template<bool Cond>
internal::conditional<Cond, ConjugateReturnType, ConstSelfAdjointView>::type Eigen::SelfAdjointView<_MatrixType, UpLo>::conjugateIf() const

Returns an expression of the complex conjugate of *this if Cond==true, returns *this otherwise.

template<typename _MatrixType, unsigned int UpLo>
MatrixType::ConstDiagonalReturnType Eigen::SelfAdjointView<_MatrixType, UpLo>::diagonal() const

Returns a const expression of the main diagonal of the matrix *this

This method simply returns the diagonal of the nested expression, thus by-passing the SelfAdjointView decorator.

template<typename _MatrixType, unsigned int UpLo>
EigenvaluesReturnType Eigen::SelfAdjointView<_MatrixType, UpLo>::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 SelfAdjointEigenSolver class. The eigenvalues are repeated according to their algebraic multiplicity, so there are as many eigenvalues as rows in the matrix.

Example:

MatrixXd ones = MatrixXd::Ones(3,3);
VectorXd eivals = ones.selfadjointView<Lower>().eigenvalues();
cout << "The eigenvalues of the 3x3 matrix of ones are:" << endl << eivals << endl;

Output:

The eigenvalues of the 3x3 matrix of ones are:
-3.09e-16
        0
        3

template<typename _MatrixType, unsigned int UpLo>
const LDLT<PlainObject, UpLo> Eigen::SelfAdjointView<_MatrixType, UpLo>::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 _MatrixType, unsigned int UpLo>
const LLT<PlainObject, UpLo> Eigen::SelfAdjointView<_MatrixType, UpLo>::llt() const

Returns the LLT decomposition of *this

This is defined in the Cholesky module. #include <Eigen/Cholesky>

template<typename _MatrixType, unsigned int UpLo> template<typename OtherDerived>
const Product<SelfAdjointView, OtherDerived> Eigen::SelfAdjointView<_MatrixType, UpLo>::operator*(const MatrixBase<OtherDerived>& rhs) const

Efficient triangular matrix times vector/matrix product

template<typename _MatrixType, unsigned int UpLo>
RealScalar Eigen::SelfAdjointView<_MatrixType, UpLo>::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 self-adjoint matrix. For a self-adjoint matrix, the operator norm is the largest eigenvalue.

The current implementation uses the eigenvalues of the matrix, as computed by eigenvalues(), to compute the operator norm of the matrix.

Example:

MatrixXd ones = MatrixXd::Ones(3,3);
cout << "The operator norm of the 3x3 matrix of ones is "
     << ones.selfadjointView<Lower>().operatorNorm() << endl;

Output:

The operator norm of the 3x3 matrix of ones is 3

template<typename _MatrixType, unsigned int UpLo> template<typename DerivedU, typename DerivedV>
SelfAdjointView& Eigen::SelfAdjointView<_MatrixType, UpLo>::rankUpdate(const MatrixBase<DerivedU>& u, const MatrixBase<DerivedV>& v, const Scalar& alpha = Scalar(1))

Returns a reference to *this

Perform a symmetric rank 2 update of the selfadjoint matrix *this: $ this = this + \alpha u v^* + conj(\alpha) v u^* $ The vectors u and v must be column vectors, however they can be a adjoint expression without any overhead. Only the meaningful triangular part of the matrix is updated, the rest is left unchanged.

template<typename _MatrixType, unsigned int UpLo> template<typename DerivedU>
SelfAdjointView& Eigen::SelfAdjointView<_MatrixType, UpLo>::rankUpdate(const MatrixBase<DerivedU>& u, const Scalar& alpha = Scalar(1))

Returns a reference to *this

Perform a symmetric rank K update of the selfadjoint matrix *this: $ this = this + \alpha ( u u^* ) $ where u is a vector or matrix.

Note that to perform $ this = this + \alpha ( u^* u ) $ you can simply call this function with u.adjoint().

template<typename _MatrixType, unsigned int UpLo>
TransposeReturnType Eigen::SelfAdjointView<_MatrixType, UpLo>::transpose()

template<typename _MatrixType, unsigned int UpLo>
const ConstTransposeReturnType Eigen::SelfAdjointView<_MatrixType, UpLo>::transpose() const

template<typename _MatrixType, unsigned int UpLo> template<unsigned int TriMode>
internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), TriangularView<MatrixType, TriMode>, TriangularView<typename MatrixType::AdjointReturnType, TriMode>>::type Eigen::SelfAdjointView<_MatrixType, UpLo>::triangularView() const

Returns an expression of a triangular view extracted from the current selfadjoint view of a given triangular part

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

If TriMode references the same triangular part than *this, then this method simply return a TriangularView of the nested expression, otherwise, the nested expression is first transposed, thus returning a TriangularView<Transpose<MatrixType>> object.

template<typename _MatrixType, unsigned int UpLo> template<typename OtherDerived>
const Product<OtherDerived, SelfAdjointView> operator*(const MatrixBase<OtherDerived>& lhs, const SelfAdjointView& rhs)

Efficient vector/matrix times triangular matrix product