Eigen  3.3.1
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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.

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

Scalar coeff (Index row, Index col) const
 
ScalarcoeffRef (Index row, Index col)
 
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))
 
template<unsigned int TriMode>
internal::conditional<(TriMode
&(Upper|Lower))==(UpLo &(Upper|Lower)),
TriangularView< MatrixType,
TriMode >, TriangularView
< typename
MatrixType::AdjointReturnType,
TriMode > >::type 
triangularView () const
 
- Public Member Functions inherited from Eigen::TriangularBase< SelfAdjointView< _MatrixType, UpLo > >
void copyCoeff (Index row, Index col, Other &other)
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
Index cols () const
 
Derived & derived ()
 
const Derived & derived () const
 
Index rows () const
 
Index size () const
 

Friends

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

Member 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

Member Function Documentation

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
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
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
template<typename MatrixType , unsigned int UpLo>
SelfAdjointView< MatrixType, UpLo >::EigenvaluesReturnType Eigen::SelfAdjointView< MatrixType, UpLo >::eigenvalues ( ) const
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:

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
See Also
SelfAdjointEigenSolver::eigenvalues(), MatrixBase::eigenvalues()
template<typename MatrixType , unsigned int UpLo>
const LDLT< typename SelfAdjointView< MatrixType, UpLo >::PlainObject, UpLo > Eigen::SelfAdjointView< MatrixType, UpLo >::ldlt ( ) const
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()
template<typename MatrixType , unsigned int UpLo>
const LLT< typename SelfAdjointView< MatrixType, UpLo >::PlainObject, UpLo > Eigen::SelfAdjointView< MatrixType, UpLo >::llt ( ) const
inline

This is defined in the Cholesky module.

#include <Eigen/Cholesky>
Returns
the LLT decomposition of *this
See Also
SelfAdjointView::llt()
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

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

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
See Also
eigenvalues(), MatrixBase::operatorNorm()
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)
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)
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
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

Friends And Related Function Documentation

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

Efficient vector/matrix times triangular matrix product


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