Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Friends
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 SelfAdjointView< MatrixType, UpLo >:

List of all members.

Public Types

enum  { Mode }
enum  
typedef TriangularBase
< SelfAdjointView
Base
typedef internal::traits
< SelfAdjointView< MatrixType,
UpLo > >::DenseMatrixType 
DenseMatrixType
typedef DenseMatrixType DenseType
typedef Matrix< RealScalar,
internal::traits< MatrixType >
::ColsAtCompileTime, 1 > 
EigenvaluesReturnType
typedef MatrixType::Index Index
typedef internal::traits
< SelfAdjointView >
::MatrixTypeNested 
MatrixTypeNested
typedef internal::traits
< SelfAdjointView >
::MatrixTypeNestedCleaned 
MatrixTypeNestedCleaned
typedef MatrixType::PlainObject PlainObject
typedef NumTraits< Scalar >::Real RealScalar
typedef internal::traits
< SelfAdjointView >::Scalar 
Scalar
 The type of coefficients in this matrix.
typedef internal::traits
< SelfAdjointView< MatrixType,
UpLo > >::StorageKind 
StorageKind

Public Member Functions

const MatrixTypeNestedCleaned_expression () const
void addTo (Dest &dst) const
void applyThisOnTheLeft (Dest &dst) const
void applyThisOnTheRight (Dest &dst) const
Scalar coeff (Index row, Index col) const
Scalar coeff (Index row, Index col) const
ScalarcoeffRef (Index row, Index col)
ScalarcoeffRef (Index row, Index col)
Index cols () const
SelfAdjointView< MatrixType,
UpLo > & 
const_cast_derived () const
const SelfAdjointView
< MatrixType, UpLo > & 
const_derived () const
void copyCoeff (Index row, Index col, Other &other)
SelfAdjointView< MatrixType,
UpLo > & 
derived ()
const SelfAdjointView
< MatrixType, UpLo > & 
derived () const
EigenvaluesReturnType eigenvalues () const
 Computes the eigenvalues of a matrix.
void evalTo (Dest &dst) const
void evalTo (MatrixBase< DenseDerived > &other) const
void evalToLazy (MatrixBase< DenseDerived > &other) const
Index innerStride () const
const LDLT< PlainObject, UpLo > ldlt () const
const LLT< PlainObject, UpLo > llt () const
const MatrixTypeNestedCleanednestedExpression () const
MatrixTypeNestedCleanednestedExpression ()
Scalar operator() (Index row, Index col) const
Scalaroperator() (Index row, Index col)
template<typename OtherDerived >
SelfadjointProductMatrix
< MatrixType, Mode, false,
OtherDerived,
0, OtherDerived::IsVectorAtCompileTime > 
operator* (const MatrixBase< OtherDerived > &rhs) const
RealScalar operatorNorm () const
 Computes the L2 operator norm.
Index outerStride () const
template<typename DerivedU , typename DerivedV >
SelfAdjointViewrankUpdate (const MatrixBase< DerivedU > &u, const MatrixBase< DerivedV > &v, Scalar alpha=Scalar(1))
template<typename DerivedU >
SelfAdjointViewrankUpdate (const MatrixBase< DerivedU > &u, Scalar alpha=Scalar(1))
Index rows () const
 SelfAdjointView (MatrixType &matrix)
Index size () const
void subTo (Dest &dst) const
DenseMatrixType toDenseMatrix () const

Protected Member Functions

void check_coordinates (Index row, Index col) const
void check_coordinates_internal (Index, Index) 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)

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

Member Typedef Documentation

typedef internal::traits<SelfAdjointView< MatrixType, UpLo > >::DenseMatrixType DenseMatrixType [inherited]
typedef DenseMatrixType DenseType [inherited]
typedef Matrix<RealScalar, internal::traits<MatrixType>::ColsAtCompileTime, 1> EigenvaluesReturnType

Return type of eigenvalues()

typedef MatrixType::Index Index
typedef MatrixType::PlainObject PlainObject
typedef NumTraits<Scalar>::Real RealScalar

Real part of Scalar

typedef internal::traits<SelfAdjointView>::Scalar Scalar

The type of coefficients in this matrix.

Reimplemented from TriangularBase< SelfAdjointView< MatrixType, UpLo > >.

typedef internal::traits<SelfAdjointView< MatrixType, UpLo > >::StorageKind StorageKind [inherited]

Member Enumeration Documentation

anonymous enum
Enumerator:
Mode 
anonymous enum [inherited]

Constructor & Destructor Documentation

SelfAdjointView ( MatrixType &  matrix) [inline]

Member Function Documentation

const MatrixTypeNestedCleaned& _expression ( ) const [inline]
void addTo ( Dest &  dst) const [inline, inherited]
void applyThisOnTheLeft ( Dest &  dst) const [inline, inherited]
void applyThisOnTheRight ( Dest &  dst) const [inline, inherited]
void check_coordinates ( Index  row,
Index  col 
) const [inline, protected, inherited]
void check_coordinates_internal ( Index  ,
Index   
) const [inline, protected, inherited]
Scalar coeff ( Index  row,
Index  col 
) const [inline, inherited]
Scalar coeff ( Index  row,
Index  col 
) const [inline]
See also:
MatrixBase::coeff()
Warning:
the coordinates must fit into the referenced triangular part
Scalar& coeffRef ( Index  row,
Index  col 
) [inline, inherited]
Scalar& coeffRef ( Index  row,
Index  col 
) [inline]
See also:
MatrixBase::coeffRef()
Warning:
the coordinates must fit into the referenced triangular part
Index cols ( void  ) const [inline]
Returns:
the number of columns.
See also:
rows(), ColsAtCompileTime

Reimplemented from TriangularBase< SelfAdjointView< MatrixType, UpLo > >.

SelfAdjointView< MatrixType, UpLo > & const_cast_derived ( ) const [inline, inherited]
const SelfAdjointView< MatrixType, UpLo > & const_derived ( ) const [inline, inherited]
void copyCoeff ( Index  row,
Index  col,
Other &  other 
) [inline, inherited]
See also:
MatrixBase::copyCoeff(row,col)

References EigenBase< Derived >::derived().

SelfAdjointView< MatrixType, UpLo > & derived ( ) [inline, inherited]
Returns:
a reference to the derived object
const SelfAdjointView< MatrixType, UpLo > & derived ( ) const [inline, inherited]
Returns:
a const reference to the derived object
SelfAdjointView< MatrixType, UpLo >::EigenvaluesReturnType 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:
-2.39e-16
8.66e-17
3
See also:
SelfAdjointEigenSolver::eigenvalues(), MatrixBase::eigenvalues()
void evalTo ( Dest &  dst) const [inline, inherited]
void evalTo ( MatrixBase< DenseDerived > &  other) const [inherited]

Assigns a triangular or selfadjoint matrix to a dense matrix. If the matrix is triangular, the opposite part is set to zero.

void evalToLazy ( MatrixBase< DenseDerived > &  other) const [inherited]

Assigns a triangular or selfadjoint matrix to a dense matrix. If the matrix is triangular, the opposite part is set to zero.

Index innerStride ( ) const [inline]
const LDLT< typename SelfAdjointView< MatrixType, UpLo >::PlainObject, 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
const LLT< typename SelfAdjointView< MatrixType, UpLo >::PlainObject, UpLo > llt ( ) const [inline]

This is defined in the Cholesky module.

 #include <Eigen/Cholesky> 
Returns:
the LLT decomposition of *this
const MatrixTypeNestedCleaned& nestedExpression ( ) const [inline]
Scalar operator() ( Index  row,
Index  col 
) const [inline, inherited]
Scalar& operator() ( Index  row,
Index  col 
) [inline, inherited]
SelfadjointProductMatrix<MatrixType,Mode,false,OtherDerived,0,OtherDerived::IsVectorAtCompileTime> operator* ( const MatrixBase< OtherDerived > &  rhs) const [inline]

Efficient self-adjoint matrix times vector/matrix product

SelfAdjointView< MatrixType, UpLo >::RealScalar 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()
Index outerStride ( ) const [inline]
SelfAdjointView< MatrixType, UpLo > & rankUpdate ( const MatrixBase< DerivedU > &  u,
const MatrixBase< DerivedV > &  v,
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)

References conj(), Eigen::Lower, Eigen::RowMajorBit, and Eigen::Upper.

SelfAdjointView< MatrixType, UpLo > & rankUpdate ( const MatrixBase< DerivedU > &  u,
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)
Index rows ( void  ) const [inline]
Returns:
the number of rows.
See also:
cols(), RowsAtCompileTime

Reimplemented from TriangularBase< SelfAdjointView< MatrixType, UpLo > >.

Index size ( ) const [inline, inherited]
Returns:
the number of coefficients, which is rows()*cols().
See also:
rows(), cols(), SizeAtCompileTime.
void subTo ( Dest &  dst) const [inline, inherited]
DenseMatrixType toDenseMatrix ( ) const [inline, inherited]

Friends And Related Function Documentation

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


Member Data Documentation


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