Eigne Matrix Class Library

Dependents:   Eigen_test Odometry_test AttitudeEstimation_usingTicker MPU9250_Quaternion_Binary_Serial ... more

Embed: (wiki syntax)

« Back to documentation index

SelfAdjointView< MatrixType, UpLo > Class Template Reference

SelfAdjointView< MatrixType, UpLo > Class Template Reference
[Core module]

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

#include <SelfAdjointView.h>

Inherits Eigen::TriangularBase< SelfAdjointView< MatrixType, UpLo > >.

Public Types

typedef internal::traits
< SelfAdjointView >::Scalar 
Scalar
 The type of coefficients in this matrix.
typedef NumTraits< Scalar >::Real RealScalar
 Real part of Scalar.
typedef Matrix< RealScalar,
internal::traits< MatrixType >
::ColsAtCompileTime, 1 > 
EigenvaluesReturnType
 Return type of eigenvalues()

Public Member Functions

Index rows () const
Index cols () const
Scalar coeff (Index row, Index col) const
ScalarcoeffRef (Index row, Index col)
template<typename OtherDerived >
SelfadjointProductMatrix
< MatrixType, Mode, false,
OtherDerived,
0, OtherDerived::IsVectorAtCompileTime > 
operator* (const MatrixBase< OtherDerived > &rhs) const
 Efficient self-adjoint matrix times vector/matrix product.
template<typename DerivedU , typename DerivedV >
SelfAdjointViewrankUpdate (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^* $.
template<typename DerivedU >
SelfAdjointViewrankUpdate (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.
const LLT< PlainObject, UpLo > llt () const
 
const LDLT< PlainObject, UpLo > ldlt () const
 
EigenvaluesReturnType eigenvalues () const
 Computes the eigenvalues of a matrix.
RealScalar operatorNorm () const
 Computes the L2 operator norm.
template<typename Other >
EIGEN_STRONG_INLINE void copyCoeff (Index row, Index col, Other &other)
template<typename DenseDerived >
void evalTo (MatrixBase< DenseDerived > &other) const
 Assigns a triangular or selfadjoint matrix to a dense matrix.
template<typename DenseDerived >
void evalToLazy (MatrixBase< DenseDerived > &other) const
 Assigns a triangular or selfadjoint matrix to a dense matrix.
Index size () const

Friends

template<typename OtherDerived >
SelfadjointProductMatrix
< OtherDerived,
0, OtherDerived::IsVectorAtCompileTime,
MatrixType, Mode, false > 
operator* (const MatrixBase< OtherDerived > &lhs, const SelfAdjointView &rhs)
 Efficient vector/matrix times self-adjoint matrix product.

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

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

Return type of eigenvalues()

Definition at line 160 of file SelfAdjointView.h.

typedef NumTraits<Scalar>::Real RealScalar

Real part of Scalar.

Definition at line 158 of file SelfAdjointView.h.

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

The type of coefficients in this matrix.

Definition at line 63 of file SelfAdjointView.h.


Member Function Documentation

Scalar coeff ( Index  row,
Index  col 
) const
See also:
MatrixBase::coeff()
Warning:
the coordinates must fit into the referenced triangular part

Definition at line 83 of file SelfAdjointView.h.

Scalar& coeffRef ( Index  row,
Index  col 
)
See also:
MatrixBase::coeffRef()
Warning:
the coordinates must fit into the referenced triangular part

Definition at line 92 of file SelfAdjointView.h.

Index cols ( void   ) const
Returns:
the number of columns.
See also:
rows(), ColsAtCompileTime

Reimplemented from EigenBase< Derived >.

Definition at line 76 of file SelfAdjointView.h.

EIGEN_STRONG_INLINE void copyCoeff ( Index  row,
Index  col,
Other &  other 
) [inherited]
See also:
MatrixBase::copyCoeff(row,col)

Definition at line 60 of file TriangularMatrix.h.

SelfAdjointView< MatrixType, UpLo >::EigenvaluesReturnType eigenvalues (  ) const

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.

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.

Definition at line 681 of file TriangularMatrix.h.

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.

Definition at line 697 of file TriangularMatrix.h.

const LDLT< typename SelfAdjointView< MatrixType, UpLo >::PlainObject, UpLo > ldlt (  ) const

Returns:
the Cholesky decomposition with full pivoting without square root of *this

Definition at line 594 of file LDLT.h.

const LLT< typename SelfAdjointView< MatrixType, UpLo >::PlainObject, UpLo > llt (  ) const

Returns:
the LLT decomposition of *this

Definition at line 491 of file LLT.h.

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

Efficient self-adjoint matrix times vector/matrix product.

Definition at line 107 of file SelfAdjointView.h.

SelfAdjointView< MatrixType, UpLo >::RealScalar operatorNorm (  ) const

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.

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.

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.

Index rows ( void   ) const
Returns:
the number of rows.
See also:
cols(), RowsAtCompileTime

Reimplemented from EigenBase< Derived >.

Definition at line 75 of file SelfAdjointView.h.


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.

Definition at line 117 of file SelfAdjointView.h.