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

Inherits TriangularBase< Derived >.

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 Member Functions

Scalar coeff (Index row, Index col) const
 
ScalarcoeffRef (Index row, Index col)
 
template<typename Other >
void copyCoeff (Index row, Index col, Other &other)
 
Derived & derived ()
 
const Derived & derived () const
 
EigenvaluesReturnType eigenvalues () const
 Computes the eigenvalues of a matrix. More...
 
template<typename DenseDerived >
void evalTo (MatrixBase< DenseDerived > &other) const
 
template<typename DenseDerived >
void evalToLazy (MatrixBase< DenseDerived > &other) const
 
const LDLT< PlainObject, UpLo > ldlt () const
 
const LLT< PlainObject, UpLo > llt () const
 
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. 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))
 
Index size () const
 

Friends

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

Member Typedef Documentation

typedef Matrix<RealScalar, internal::traits<MatrixType>::ColsAtCompileTime, 1> EigenvaluesReturnType

Return type of eigenvalues()

typedef NumTraits<Scalar>::Real RealScalar

Real part of Scalar

Member Function Documentation

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
See Also
MatrixBase::coeffRef()
Warning
the coordinates must fit into the referenced triangular part
void copyCoeff ( Index  row,
Index  col,
Other &  other 
)
inlineinherited
See Also
MatrixBase::copyCoeff(row,col)
Derived& derived ( )
inlineinherited
Returns
a reference to the derived object

Referenced by IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::analyzePattern(), MatrixBase< Derived >::applyOnTheLeft(), MatrixBase< Derived >::applyOnTheRight(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::applyTranspositionOnTheLeft(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::applyTranspositionOnTheRight(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::bottomRows(), EigenBase< SparseSymmetricPermutationProduct< MatrixType, UpLo > >::cols(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::cols(), IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::compute(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::eval(), IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::factorize(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::indices(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::inverse(), IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::IterativeSolverBase(), RotationBase< Derived, 3 >::operator*(), SparseSelfAdjointView< MatrixType, UpLo >::operator*(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::operator*(), SparseMatrixBase< Derived >::operator*(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::operator*(), Transform< Scalar, Dim, Mode, _Options >::operator*(), Eigen::operator*(), MatrixBase< Derived >::operator*=(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::operator=(), DenseBase< Derived >::operator=(), Transform< Scalar, Dim, Mode, _Options >::operator=(), PlainObjectBase< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > >::operator=(), PlainObjectBase< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > >::PlainObjectBase(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::real(), PlainObjectBase< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > >::resizeLike(), EigenBase< SparseSymmetricPermutationProduct< MatrixType, UpLo > >::rows(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::rows(), SimplicialCholeskyBase< SimplicialLDLT< _MatrixType, _UpLo, _Ordering > >::solve(), IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::solve(), SparseLU< _MatrixType, _OrderingType >::solve(), UmfPackLU< _MatrixType >::solve(), CholmodBase< _MatrixType, _UpLo, CholmodSimplicialLLT< _MatrixType, _UpLo > >::solve(), SuperLUBase< _MatrixType, SuperILU< _MatrixType > >::solve(), SparseMatrix< Scalar, RowMajor >::SparseMatrix(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::toDenseMatrix(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::topLeftCorner(), Transform< Scalar, Dim, Mode, _Options >::Transform(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::transpose(), and SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::twistedBy().

const Derived& derived ( ) const
inlineinherited
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:
-3.09e-16
        0
        3
See Also
SelfAdjointEigenSolver::eigenvalues(), MatrixBase::eigenvalues()
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.

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
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()
SelfAdjointView& 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)
SelfAdjointView& 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)
Index size ( ) const
inlineinherited
Returns
the number of coefficients, which is rows()*cols().
See Also
rows(), cols(), SizeAtCompileTime.

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


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