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 >:
Inheritance graph
[legend]

List of all members.

Public Types

enum  { Mode = internal::traits<SelfAdjointView>::Mode }
typedef TriangularBase
< SelfAdjointView
Base
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.

Public Member Functions

const MatrixTypeNestedCleaned_expression () const
Scalar coeff (Index row, Index col) const
ScalarcoeffRef (Index row, Index col)
Index cols () const
EigenvaluesReturnType eigenvalues () const
 Computes the eigenvalues of a matrix.
Index innerStride () const
const LDLT< PlainObject, UpLo > ldlt () const
const LLT< PlainObject, UpLo > llt () const
const MatrixTypeNestedCleanednestedExpression () const
MatrixTypeNestedCleanednestedExpression ()
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, const Scalar &alpha=Scalar(1))
template<typename DerivedU >
SelfAdjointViewrankUpdate (const MatrixBase< DerivedU > &u, const Scalar &alpha=Scalar(1))
Index rows () const
 SelfAdjointView (MatrixType &matrix)

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

Definition at line 53 of file SelfAdjointView.h.


Member Typedef Documentation

template<typename MatrixType, unsigned int UpLo>
typedef TriangularBase<SelfAdjointView> Eigen::SelfAdjointView< MatrixType, UpLo >::Base

Definition at line 58 of file SelfAdjointView.h.

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 MatrixType::Index Eigen::SelfAdjointView< MatrixType, UpLo >::Index
template<typename MatrixType, unsigned int UpLo>
typedef internal::traits<SelfAdjointView>::MatrixTypeNested Eigen::SelfAdjointView< MatrixType, UpLo >::MatrixTypeNested

Definition at line 59 of file SelfAdjointView.h.

template<typename MatrixType, unsigned int UpLo>
typedef internal::traits<SelfAdjointView>::MatrixTypeNestedCleaned Eigen::SelfAdjointView< MatrixType, UpLo >::MatrixTypeNestedCleaned

Definition at line 60 of file SelfAdjointView.h.

template<typename MatrixType, unsigned int UpLo>
typedef MatrixType::PlainObject Eigen::SelfAdjointView< MatrixType, UpLo >::PlainObject

Definition at line 70 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.

template<typename MatrixType, unsigned int UpLo>
typedef internal::traits<SelfAdjointView>::Scalar Eigen::SelfAdjointView< MatrixType, UpLo >::Scalar

The type of coefficients in this matrix.

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

Definition at line 63 of file SelfAdjointView.h.


Member Enumeration Documentation

template<typename MatrixType, unsigned int UpLo>
anonymous enum
Enumerator:
Mode 

Definition at line 67 of file SelfAdjointView.h.


Constructor & Destructor Documentation

template<typename MatrixType, unsigned int UpLo>
Eigen::SelfAdjointView< MatrixType, UpLo >::SelfAdjointView ( MatrixType &  matrix) [inline]

Definition at line 72 of file SelfAdjointView.h.


Member Function Documentation

template<typename MatrixType, unsigned int UpLo>
const MatrixTypeNestedCleaned& Eigen::SelfAdjointView< MatrixType, UpLo >::_expression ( ) const [inline]

Definition at line 99 of file SelfAdjointView.h.

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.

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.

template<typename MatrixType, unsigned int UpLo>
Index Eigen::SelfAdjointView< MatrixType, UpLo >::cols ( void  ) const [inline]
Returns:
the number of columns.
See also:
rows(), ColsAtCompileTime

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

Definition at line 76 of file SelfAdjointView.h.

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.

template<typename MatrixType, unsigned int UpLo>
Index Eigen::SelfAdjointView< MatrixType, UpLo >::innerStride ( ) const [inline]
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.

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.

template<typename MatrixType, unsigned int UpLo>
const MatrixTypeNestedCleaned& Eigen::SelfAdjointView< MatrixType, UpLo >::nestedExpression ( ) const [inline]

Definition at line 101 of file SelfAdjointView.h.

template<typename MatrixType, unsigned int UpLo>
MatrixTypeNestedCleaned& Eigen::SelfAdjointView< MatrixType, UpLo >::nestedExpression ( ) [inline]

Definition at line 102 of file SelfAdjointView.h.

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.

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.

template<typename MatrixType, unsigned int UpLo>
Index Eigen::SelfAdjointView< MatrixType, UpLo >::outerStride ( ) const [inline]
template<typename MatrixType , unsigned int UpLo>
template<typename DerivedU , typename DerivedV >
SelfAdjointView< MatrixType, UpLo > & 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)

Definition at line 61 of file SelfadjointRank2Update.h.

template<typename MatrixType , unsigned int UpLo>
template<typename DerivedU >
SelfAdjointView< MatrixType, UpLo > & 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)

Definition at line 114 of file SelfadjointProduct.h.

template<typename MatrixType, unsigned int UpLo>
Index Eigen::SelfAdjointView< MatrixType, UpLo >::rows ( void  ) const [inline]
Returns:
the number of rows.
See also:
cols(), RowsAtCompileTime

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

Definition at line 75 of file SelfAdjointView.h.


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.


Member Data Documentation

template<typename MatrixType, unsigned int UpLo>
MatrixTypeNested Eigen::SelfAdjointView< MatrixType, UpLo >::m_matrix [protected]

Definition at line 189 of file SelfAdjointView.h.


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


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Sat Jun 8 2019 19:40:53