Shapeworks Studio  2.1
Shape analysis software suite
List of all members | Public Types | Public Member Functions | Protected Attributes | Friends
Eigen::SelfAdjointView< MatrixType, UpLo > Class Template Reference

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

#include <SelfAdjointView.h>

+ Inheritance diagram for Eigen::SelfAdjointView< MatrixType, UpLo >:
+ Collaboration diagram for Eigen::SelfAdjointView< MatrixType, UpLo >:

Public Types

enum  { Mode = internal::traits<SelfAdjointView>::Mode }
 
typedef TriangularBase< SelfAdjointViewBase
 
typedef internal::traits< SelfAdjointView >::MatrixTypeNested MatrixTypeNested
 
typedef internal::traits< SelfAdjointView >::MatrixTypeNestedCleaned MatrixTypeNestedCleaned
 
typedef internal::traits< SelfAdjointView >::Scalar Scalar
 The type of coefficients in this matrix.
 
typedef MatrixType::Index Index
 
typedef MatrixType::PlainObject PlainObject
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef Matrix< RealScalar, internal::traits< MatrixType >::ColsAtCompileTime, 1 > EigenvaluesReturnType
 
- Public Types inherited from Eigen::TriangularBase< SelfAdjointView< MatrixType, UpLo > >
enum  
 
typedef internal::traits< SelfAdjointView< MatrixType, UpLo > >::Scalar Scalar
 
typedef internal::traits< SelfAdjointView< MatrixType, UpLo > >::StorageKind StorageKind
 
typedef internal::traits< SelfAdjointView< MatrixType, UpLo > >::Index Index
 
typedef internal::traits< SelfAdjointView< MatrixType, UpLo > >::DenseMatrixType DenseMatrixType
 
typedef DenseMatrixType DenseType
 
- Public Types inherited from Eigen::EigenBase< Derived >
typedef internal::traits< Derived >::StorageKind StorageKind
 
typedef internal::traits< Derived >::Index Index
 

Public Member Functions

 SelfAdjointView (MatrixType &matrix)
 
Index rows () const
 
Index cols () const
 
Index outerStride () const
 
Index innerStride () const
 
Scalar coeff (Index row, Index col) const
 
ScalarcoeffRef (Index row, Index col)
 
const MatrixTypeNestedCleaned & _expression () const
 
const MatrixTypeNestedCleaned & nestedExpression () const
 
MatrixTypeNestedCleaned & nestedExpression ()
 
template<typename OtherDerived >
SelfadjointProductMatrix< MatrixType, Mode, false, OtherDerived, 0, OtherDerived::IsVectorAtCompileTime > operator* (const MatrixBase< OtherDerived > &rhs) const
 
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 LLT< PlainObject, UpLo > llt () const
 
const LDLT< PlainObject, UpLo > ldlt () const
 
EigenvaluesReturnType eigenvalues () const
 Computes the eigenvalues of a matrix. More...
 
RealScalar operatorNorm () const
 Computes the L2 operator norm. More...
 
template<typename DerivedU >
SelfAdjointView< MatrixType, UpLo > & rankUpdate (const MatrixBase< DerivedU > &u, const Scalar &alpha)
 
template<typename DerivedU , typename DerivedV >
SelfAdjointView< MatrixType, UpLo > & rankUpdate (const MatrixBase< DerivedU > &u, const MatrixBase< DerivedV > &v, const Scalar &alpha)
 
- Public Member Functions inherited from Eigen::TriangularBase< SelfAdjointView< MatrixType, UpLo > >
Index rows () const
 
Index cols () const
 
Index outerStride () const
 
Index innerStride () const
 
Scalar coeff (Index row, Index col) const
 
Scalar & coeffRef (Index row, Index col)
 
EIGEN_STRONG_INLINE void copyCoeff (Index row, Index col, Other &other)
 
Scalar operator() (Index row, Index col) const
 
Scalar & operator() (Index row, Index col)
 
const SelfAdjointView< MatrixType, UpLo > & derived () const
 
SelfAdjointView< MatrixType, UpLo > & derived ()
 
void evalTo (MatrixBase< DenseDerived > &other) const
 
void evalToLazy (MatrixBase< DenseDerived > &other) const
 
DenseMatrixType toDenseMatrix () const
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
Derived & derived ()
 
const Derived & derived () const
 
Derived & const_cast_derived () const
 
const Derived & const_derived () const
 
Index rows () const
 
Index cols () const
 
Index size () const
 
template<typename Dest >
void evalTo (Dest &dst) const
 
template<typename Dest >
void addTo (Dest &dst) const
 
template<typename Dest >
void subTo (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheRight (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheLeft (Dest &dst) const
 

Protected Attributes

MatrixTypeNested m_matrix
 

Friends

template<typename OtherDerived >
SelfadjointProductMatrix< OtherDerived, 0, OtherDerived::IsVectorAtCompileTime, MatrixType, Mode, false > operator* (const MatrixBase< OtherDerived > &lhs, const SelfAdjointView &rhs)
 

Additional Inherited Members

- Protected Member Functions inherited from Eigen::TriangularBase< SelfAdjointView< MatrixType, UpLo > >
void check_coordinates (Index row, Index col) const
 
void check_coordinates_internal (Index, Index) const
 

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

Definition at line 53 of file SelfAdjointView.h.

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

Definition at line 160 of file SelfAdjointView.h.

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

Real part of Scalar

Definition at line 158 of file SelfAdjointView.h.

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

Definition at line 83 of file SelfAdjointView.h.

84  {
85  Base::check_coordinates_internal(row, col);
86  return m_matrix.coeff(row, col);
87  }
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

Definition at line 92 of file SelfAdjointView.h.

93  {
94  Base::check_coordinates_internal(row, col);
95  return m_matrix.const_cast_derived().coeffRef(row, col);
96  }
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 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:

Output:

See also
SelfAdjointEigenSolver::eigenvalues(), MatrixBase::eigenvalues()

Definition at line 89 of file MatrixBaseEigenvalues.h.

90 {
91  typedef typename SelfAdjointView<MatrixType, UpLo>::PlainObject PlainObject;
92  PlainObject thisAsMatrix(*this);
93  return SelfAdjointEigenSolver<PlainObject>(thisAsMatrix, false).eigenvalues();
94 }
EigenvaluesReturnType eigenvalues() const
Computes the eigenvalues of a matrix.
template<typename MatrixType , unsigned int UpLo>
const LDLT< typename SelfAdjointView< MatrixType, UpLo >::PlainObject, UpLo > Eigen::SelfAdjointView< MatrixType, UpLo >::ldlt ( ) const
inline
Returns
the Cholesky decomposition with full pivoting without square root of *this

Definition at line 583 of file LDLT.h.

584 {
585  return LDLT<PlainObject,UpLo>(m_matrix);
586 }
template<typename MatrixType , unsigned int UpLo>
const LLT< typename SelfAdjointView< MatrixType, UpLo >::PlainObject, UpLo > Eigen::SelfAdjointView< MatrixType, UpLo >::llt ( ) const
inline
Returns
the LLT decomposition of *this

Definition at line 483 of file LLT.h.

484 {
485  return LLT<PlainObject,UpLo>(m_matrix);
486 }
template<typename MatrixType, unsigned int UpLo>
template<typename OtherDerived >
SelfadjointProductMatrix<MatrixType,Mode,false,OtherDerived,0,OtherDerived::IsVectorAtCompileTime> Eigen::SelfAdjointView< MatrixType, UpLo >::operator* ( const MatrixBase< OtherDerived > &  rhs) const
inline

Efficient self-adjoint matrix times vector/matrix product

Definition at line 107 of file SelfAdjointView.h.

108  {
109  return SelfadjointProductMatrix
110  <MatrixType,Mode,false,OtherDerived,0,OtherDerived::IsVectorAtCompileTime>
111  (m_matrix, rhs.derived());
112  }
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 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:

Output:

See also
eigenvalues(), MatrixBase::operatorNorm()

Definition at line 153 of file MatrixBaseEigenvalues.h.

154 {
155  return eigenvalues().cwiseAbs().maxCoeff();
156 }
EigenvaluesReturnType eigenvalues() const
Computes the eigenvalues of a matrix.
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)

Friends And Related Function Documentation

template<typename MatrixType, unsigned int UpLo>
template<typename OtherDerived >
SelfadjointProductMatrix<OtherDerived,0,OtherDerived::IsVectorAtCompileTime,MatrixType,Mode,false> operator* ( const MatrixBase< OtherDerived > &  lhs,
const SelfAdjointView< MatrixType, UpLo > &  rhs 
)
friend

Efficient vector/matrix times self-adjoint matrix product

Definition at line 117 of file SelfAdjointView.h.

118  {
119  return SelfadjointProductMatrix
120  <OtherDerived,0,OtherDerived::IsVectorAtCompileTime,MatrixType,Mode,false>
121  (lhs.derived(),rhs.m_matrix);
122  }

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