Eigen  3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
Eigen::SelfAdjointView< MatrixType_, UpLo > Class Template Reference

Detailed Description

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

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

Template Parameters
MatrixTypethe type of the dense matrix storing the coefficients
TriangularPartcan be either Lower or Upper

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.

See also
class TriangularBase, MatrixBase::selfadjointView()
+ Inheritance diagram for Eigen::SelfAdjointView< MatrixType_, UpLo >:

Public Types

typedef Matrix< RealScalar, internal::traits< MatrixType >::ColsAtCompileTime, 1 > EigenvaluesReturnType
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef internal::traits< SelfAdjointView >::Scalar Scalar
 The type of coefficients in this matrix.
 
- Public Types inherited from Eigen::EigenBase< Derived >
typedef Eigen::Index Index
 The interface type of indices. More...
 

Public Member Functions

const AdjointReturnType adjoint () const
 
Scalar coeff (Index row, Index col) const
 
ScalarcoeffRef (Index row, Index col)
 
const ConjugateReturnType conjugate () const
 
template<bool Cond>
std::conditional_t< Cond, ConjugateReturnType, ConstSelfAdjointViewconjugateIf () const
 
MatrixType::ConstDiagonalReturnType diagonal () const
 
EigenvaluesReturnType eigenvalues () const
 Computes the eigenvalues of a matrix. More...
 
const LDLT< PlainObject, UpLo > ldlt () const
 
const LLT< PlainObject, UpLo > llt () const
 
template<typename OtherDerived >
const Product< SelfAdjointView, OtherDerived > operator* (const MatrixBase< OtherDerived > &rhs) const
 
RealScalar operatorNorm () const
 Computes the L2 operator norm. More...
 
template<typename DerivedU , typename DerivedV >
SelfAdjointViewrankUpdate (const MatrixBase< DerivedU > &u, const MatrixBase< DerivedV > &v, const Scalar &alpha=Scalar(1))
 
template<typename DerivedU >
SelfAdjointViewrankUpdate (const MatrixBase< DerivedU > &u, const Scalar &alpha=Scalar(1))
 
const ConstTransposeReturnType transpose () const
 
template<class Dummy = int>
TransposeReturnType transpose (std::enable_if_t< Eigen::internal::is_lvalue< MatrixType >::value, Dummy * >=nullptr)
 
template<unsigned int TriMode>
std::conditional_t<(TriMode &(Upper|Lower))==(UpLo &(Upper|Lower)), TriangularView< MatrixType, TriMode >, TriangularView< typename MatrixType::AdjointReturnType, TriMode > > triangularView () const
 
- Public Member Functions inherited from Eigen::TriangularBase< SelfAdjointView< MatrixType_, UpLo > >
void copyCoeff (Index row, Index col, Other &other)
 
void evalTo (MatrixBase< DenseDerived > &other) const
 
void evalToLazy (MatrixBase< DenseDerived > &other) const
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & derived ()
 
const Derived & derived () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 

Member Typedef Documentation

◆ EigenvaluesReturnType

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

Return type of eigenvalues()

◆ RealScalar

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

Real part of Scalar

Member Function Documentation

◆ adjoint()

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

◆ coeff()

template<typename MatrixType_ , unsigned int UpLo>
Scalar Eigen::SelfAdjointView< MatrixType_, UpLo >::coeff ( Index  row,
Index  col 
) const
inline
See also
MatrixBase::coeff()
Warning
the coordinates must fit into the referenced triangular part

◆ coeffRef()

template<typename MatrixType_ , unsigned int UpLo>
Scalar& Eigen::SelfAdjointView< MatrixType_, UpLo >::coeffRef ( Index  row,
Index  col 
)
inline
See also
MatrixBase::coeffRef()
Warning
the coordinates must fit into the referenced triangular part

◆ conjugate()

template<typename MatrixType_ , unsigned int UpLo>
const ConjugateReturnType Eigen::SelfAdjointView< MatrixType_, UpLo >::conjugate ( void  ) const
inline
See also
MatrixBase::conjugate() const

◆ conjugateIf()

template<typename MatrixType_ , unsigned int UpLo>
template<bool Cond>
std::conditional_t<Cond,ConjugateReturnType,ConstSelfAdjointView> Eigen::SelfAdjointView< MatrixType_, UpLo >::conjugateIf ( ) const
inline
Returns
an expression of the complex conjugate of *this if Cond==true, returns *this otherwise.

◆ diagonal()

template<typename MatrixType_ , unsigned int UpLo>
MatrixType::ConstDiagonalReturnType Eigen::SelfAdjointView< MatrixType_, UpLo >::diagonal ( ) const
inline
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.

See also
MatrixBase::diagonal(), class Diagonal

◆ eigenvalues()

template<typename MatrixType , unsigned int UpLo>
SelfAdjointView< MatrixType, UpLo >::EigenvaluesReturnType Eigen::SelfAdjointView< MatrixType, UpLo >::eigenvalues
inline

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:

VectorXd eivals = ones.selfadjointView<Lower>().eigenvalues();
cout << "The eigenvalues of the 3x3 matrix of ones are:" << endl << eivals << endl;
static const ConstantReturnType Ones()
Definition: CwiseNullaryOp.h:672
EigenvaluesReturnType eigenvalues() const
Computes the eigenvalues of a matrix.
Definition: MatrixBaseEigenvalues.h:90
@ Lower
Definition: Constants.h:211
Matrix< double, Dynamic, 1 > VectorXd
DynamicĂ—1 vector of type double.
Definition: Matrix.h:501
Matrix< double, Dynamic, Dynamic > MatrixXd
DynamicĂ—Dynamic matrix of type double.
Definition: Matrix.h:501

Output:

The eigenvalues of the 3x3 matrix of ones are:
-3.09e-16
        0
        3
See also
SelfAdjointEigenSolver::eigenvalues(), MatrixBase::eigenvalues()

◆ ldlt()

template<typename MatrixType , unsigned int UpLo>
const LDLT< typename SelfAdjointView< MatrixType, UpLo >::PlainObject, UpLo > Eigen::SelfAdjointView< MatrixType, UpLo >::ldlt
inline

This is defined in the Cholesky module.

#include <Eigen/Cholesky>
Returns
the Cholesky decomposition with full pivoting without square root of *this
See also
MatrixBase::ldlt()

◆ llt()

template<typename MatrixType , unsigned int UpLo>
const LLT< typename SelfAdjointView< MatrixType, UpLo >::PlainObject, UpLo > Eigen::SelfAdjointView< MatrixType, UpLo >::llt
inline

This is defined in the Cholesky module.

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

◆ operator*()

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

Efficient triangular matrix times vector/matrix product

◆ operatorNorm()

template<typename MatrixType , unsigned int UpLo>
SelfAdjointView< MatrixType, UpLo >::RealScalar Eigen::SelfAdjointView< MatrixType, UpLo >::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 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:

cout << "The operator norm of the 3x3 matrix of ones is "
<< ones.selfadjointView<Lower>().operatorNorm() << endl;
RealScalar operatorNorm() const
Computes the L2 operator norm.
Definition: MatrixBaseEigenvalues.h:153

Output:

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

◆ rankUpdate() [1/2]

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) 
)

Perform a symmetric rank 2 update of the selfadjoint matrix *this: \( this = this + \alpha u v^* + conj(\alpha) v u^* \)

Returns
a reference to *this

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.

See also
rankUpdate(const MatrixBase<DerivedU>&, Scalar)

◆ rankUpdate() [2/2]

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

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

Returns
a reference to *this

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

See also
rankUpdate(const MatrixBase<DerivedU>&, const MatrixBase<DerivedV>&, Scalar)

◆ transpose() [1/2]

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

◆ transpose() [2/2]

template<typename MatrixType_ , unsigned int UpLo>
template<class Dummy = int>
TransposeReturnType Eigen::SelfAdjointView< MatrixType_, UpLo >::transpose ( std::enable_if_t< Eigen::internal::is_lvalue< MatrixType >::value, Dummy * >  = nullptr)
inline

◆ triangularView()

template<typename MatrixType_ , unsigned int UpLo>
template<unsigned int TriMode>
std::conditional_t<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), TriangularView<MatrixType,TriMode>, TriangularView<typename MatrixType::AdjointReturnType,TriMode> > Eigen::SelfAdjointView< MatrixType_, UpLo >::triangularView ( ) const
inline
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.

See also
MatrixBase::triangularView(), class TriangularView

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