Eigen::SelfAdjointView< MatrixType_, UpLo > Class Template Reference

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

+ Inheritance diagram for Eigen::SelfAdjointView< MatrixType_, UpLo >:

Public Types

enum  {
  Mode ,
  Flags ,
  TransposeMode
}
 
typedef SelfAdjointView< const typename MatrixType::AdjointReturnType, TransposeModeAdjointReturnType
 
typedef TriangularBase< SelfAdjointViewBase
 
typedef SelfAdjointView< const MatrixConjugateReturnType, UpLo > ConjugateReturnType
 
typedef SelfAdjointView< std::add_const_t< MatrixType >, UpLo > ConstSelfAdjointView
 
typedef SelfAdjointView< const typename MatrixType::ConstTransposeReturnType, TransposeModeConstTransposeReturnType
 
typedef Matrix< RealScalar, internal::traits< MatrixType >::ColsAtCompileTime, 1 > EigenvaluesReturnType
 
typedef internal::remove_all_t< typename MatrixType::ConjugateReturnType > MatrixConjugateReturnType
 
typedef MatrixType_ MatrixType
 
typedef internal::traits< SelfAdjointView >::MatrixTypeNested MatrixTypeNested
 
typedef internal::traits< SelfAdjointView >::MatrixTypeNestedCleaned MatrixTypeNestedCleaned
 
typedef MatrixTypeNestedCleaned NestedExpression
 
typedef MatrixType::PlainObject PlainObject
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef internal::traits< SelfAdjointView >::Scalar Scalar
 The type of coefficients in this matrix. More...
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef SelfAdjointView< typename MatrixType::TransposeReturnType, TransposeModeTransposeReturnType
 
- Public Types inherited from Eigen::TriangularBase< SelfAdjointView< MatrixType_, UpLo > >
enum  
 
typedef internal::traits< SelfAdjointView< MatrixType_, UpLo > >::FullMatrixType DenseMatrixType
 
typedef DenseMatrixType DenseType
 
typedef SelfAdjointView< MatrixType_, UpLo > const & Nested
 
typedef internal::traits< SelfAdjointView< MatrixType_, UpLo > >::Scalar Scalar
 
typedef internal::traits< SelfAdjointView< MatrixType_, UpLo > >::StorageIndex StorageIndex
 
typedef internal::traits< SelfAdjointView< MatrixType_, UpLo > >::StorageKind StorageKind
 
- Public Types inherited from Eigen::EigenBase< Derived >
typedef Eigen::Index Index
 The interface type of indices. More...
 
typedef internal::traits< Derived >::StorageKind StorageKind
 

Public Member Functions

const MatrixTypeNestedCleaned_expression () const
 
const AdjointReturnType adjoint () const
 
Scalar coeff (Index row, Index col) const
 
ScalarcoeffRef (Index row, Index col)
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
const ConjugateReturnType conjugate () const
 
template<bool Cond>
std::conditional_t< Cond, ConjugateReturnType, ConstSelfAdjointViewconjugateIf () const
 
MatrixType::ConstDiagonalReturnType diagonal () const
 
EigenvaluesReturnType eigenvalues () const
 Computes the eigenvalues of a matrix. More...
 
EIGEN_CONSTEXPR Index innerStride () const EIGEN_NOEXCEPT
 
const LDLT< PlainObject, UpLo > ldlt () const
 
const LLT< PlainObject, UpLo > llt () const
 
MatrixTypeNestedCleanednestedExpression ()
 
const MatrixTypeNestedCleanednestedExpression () const
 
template<typename OtherDerived >
const Product< SelfAdjointView, OtherDerived > operator* (const MatrixBase< OtherDerived > &rhs) const
 
RealScalar operatorNorm () const
 Computes the L2 operator norm. More...
 
EIGEN_CONSTEXPR Index outerStride () const EIGEN_NOEXCEPT
 
template<typename DerivedU , typename DerivedV >
SelfAdjointView< MatrixType, UpLo > & rankUpdate (const MatrixBase< DerivedU > &u, const MatrixBase< DerivedV > &v, const Scalar &alpha)
 
template<typename DerivedU , typename DerivedV >
SelfAdjointViewrankUpdate (const MatrixBase< DerivedU > &u, const MatrixBase< DerivedV > &v, const Scalar &alpha=Scalar(1))
 
template<typename DerivedU >
SelfAdjointView< MatrixType, UpLo > & rankUpdate (const MatrixBase< DerivedU > &u, const Scalar &alpha)
 
template<typename DerivedU >
SelfAdjointViewrankUpdate (const MatrixBase< DerivedU > &u, const Scalar &alpha=Scalar(1))
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
 SelfAdjointView (MatrixType &matrix)
 
const ConstTransposeReturnType transpose () const
 
template<class Dummy = int>
TransposeReturnType transpose (std::enable_if_t< Eigen::internal::is_lvalue< MatrixType >::value, Dummy * >=nullptr)
 
template<unsigned int TriMode>
std::conditional_t<(TriMode &(Upper|Lower))==(UpLo &(Upper|Lower)), TriangularView< MatrixType, TriMode >, TriangularView< typename MatrixType::AdjointReturnType, TriMode > > triangularView () const
 
- Public Member Functions inherited from Eigen::TriangularBase< SelfAdjointView< MatrixType_, UpLo > >
Scalar coeff (Index row, Index col) const
 
ScalarcoeffRef (Index row, Index col)
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
void copyCoeff (Index row, Index col, Other &other)
 
void evalTo (MatrixBase< DenseDerived > &other) const
 
void evalToLazy (MatrixBase< DenseDerived > &other) const
 
EIGEN_CONSTEXPR Index innerStride () const EIGEN_NOEXCEPT
 
Scalaroperator() (Index row, Index col)
 
Scalar operator() (Index row, Index col) const
 
EIGEN_CONSTEXPR Index outerStride () const EIGEN_NOEXCEPT
 
void resize (Index rows, Index cols)
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
DenseMatrixType toDenseMatrix () const
 
 TriangularBase ()
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
template<typename Dest >
void addTo (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheLeft (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheRight (Dest &dst) const
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & const_cast_derived () const
 
const Derived & const_derived () const
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Dest >
void evalTo (Dest &dst) const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
template<typename Dest >
void subTo (Dest &dst) const
 

Protected Attributes

MatrixTypeNested m_matrix
 

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.

Template 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 51 of file SelfAdjointView.h.

Member Typedef Documentation

◆ AdjointReturnType

template<typename MatrixType_ , unsigned int UpLo>
typedef SelfAdjointView<const typename MatrixType::AdjointReturnType,TransposeMode> Eigen::SelfAdjointView< MatrixType_, UpLo >::AdjointReturnType

Definition at line 213 of file SelfAdjointView.h.

◆ Base

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

Definition at line 58 of file SelfAdjointView.h.

◆ ConjugateReturnType

template<typename MatrixType_ , unsigned int UpLo>
typedef SelfAdjointView<const MatrixConjugateReturnType,UpLo> Eigen::SelfAdjointView< MatrixType_, UpLo >::ConjugateReturnType

Definition at line 195 of file SelfAdjointView.h.

◆ ConstSelfAdjointView

template<typename MatrixType_ , unsigned int UpLo>
typedef SelfAdjointView<std::add_const_t<MatrixType>, UpLo> Eigen::SelfAdjointView< MatrixType_, UpLo >::ConstSelfAdjointView

Definition at line 67 of file SelfAdjointView.h.

◆ ConstTransposeReturnType

template<typename MatrixType_ , unsigned int UpLo>
typedef SelfAdjointView<const typename MatrixType::ConstTransposeReturnType,TransposeMode> Eigen::SelfAdjointView< MatrixType_, UpLo >::ConstTransposeReturnType

Definition at line 229 of file SelfAdjointView.h.

◆ EigenvaluesReturnType

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 258 of file SelfAdjointView.h.

◆ MatrixConjugateReturnType

template<typename MatrixType_ , unsigned int UpLo>
typedef internal::remove_all_t<typename MatrixType::ConjugateReturnType> Eigen::SelfAdjointView< MatrixType_, UpLo >::MatrixConjugateReturnType

Definition at line 66 of file SelfAdjointView.h.

◆ MatrixType

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

Definition at line 57 of file SelfAdjointView.h.

◆ MatrixTypeNested

template<typename MatrixType_ , unsigned int UpLo>
typedef internal::traits<SelfAdjointView>::MatrixTypeNested Eigen::SelfAdjointView< MatrixType_, UpLo >::MatrixTypeNested

Definition at line 59 of file SelfAdjointView.h.

◆ MatrixTypeNestedCleaned

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

Definition at line 60 of file SelfAdjointView.h.

◆ NestedExpression

template<typename MatrixType_ , unsigned int UpLo>
typedef MatrixTypeNestedCleaned Eigen::SelfAdjointView< MatrixType_, UpLo >::NestedExpression

Definition at line 61 of file SelfAdjointView.h.

◆ PlainObject

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

Definition at line 74 of file SelfAdjointView.h.

◆ RealScalar

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

Real part of Scalar

Definition at line 256 of file SelfAdjointView.h.

◆ Scalar

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

The type of coefficients in this matrix.

Definition at line 64 of file SelfAdjointView.h.

◆ StorageIndex

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

Definition at line 65 of file SelfAdjointView.h.

◆ TransposeReturnType

template<typename MatrixType_ , unsigned int UpLo>
typedef SelfAdjointView<typename MatrixType::TransposeReturnType,TransposeMode> Eigen::SelfAdjointView< MatrixType_, UpLo >::TransposeReturnType

Definition at line 219 of file SelfAdjointView.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , unsigned int UpLo>
anonymous enum
Enumerator
Mode 
Flags 
TransposeMode 

Definition at line 69 of file SelfAdjointView.h.

69  {
70  Mode = internal::traits<SelfAdjointView>::Mode,
71  Flags = internal::traits<SelfAdjointView>::Flags,
72  TransposeMode = ((int(Mode) & int(Upper)) ? Lower : 0) | ((int(Mode) & int(Lower)) ? Upper : 0)
73  };
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213

Constructor & Destructor Documentation

◆ SelfAdjointView()

template<typename MatrixType_ , unsigned int UpLo>
Eigen::SelfAdjointView< MatrixType_, UpLo >::SelfAdjointView ( MatrixType matrix)
inlineexplicit

Definition at line 77 of file SelfAdjointView.h.

77 : m_matrix(matrix) { }
MatrixTypeNested m_matrix

Member Function Documentation

◆ _expression()

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

Definition at line 111 of file SelfAdjointView.h.

111 { return m_matrix; }

◆ adjoint()

template<typename MatrixType_ , unsigned int UpLo>
const AdjointReturnType Eigen::SelfAdjointView< MatrixType_, UpLo >::adjoint ( ) const
inline
See also
MatrixBase::adjoint() const

Definition at line 216 of file SelfAdjointView.h.

217  { return AdjointReturnType(m_matrix.adjoint()); }
SelfAdjointView< const typename MatrixType::AdjointReturnType, TransposeMode > AdjointReturnType

◆ coeff()

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 92 of file SelfAdjointView.h.

93  {
95  return m_matrix.coeff(row, col);
96  }
RowXpr row(Index i)
This is the const version of row(). *‍/.
ColXpr col(Index i)
This is the const version of col().
void check_coordinates_internal(Index, Index) const

◆ coeffRef()

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 102 of file SelfAdjointView.h.

103  {
106  return m_matrix.coeffRef(row, col);
107  }
#define EIGEN_STATIC_ASSERT_LVALUE(Derived)
Definition: StaticAssert.h:96
SelfAdjointView(MatrixType &matrix)

◆ cols()

template<typename MatrixType_ , unsigned int UpLo>
EIGEN_CONSTEXPR Index Eigen::SelfAdjointView< MatrixType_, UpLo >::cols ( void  ) const
inline

Definition at line 82 of file SelfAdjointView.h.

82 { return m_matrix.cols(); }

◆ conjugate()

template<typename MatrixType_ , unsigned int UpLo>
const ConjugateReturnType Eigen::SelfAdjointView< MatrixType_, UpLo >::conjugate ( void  ) const
inline
See also
MatrixBase::conjugate() const

Definition at line 198 of file SelfAdjointView.h.

199  { return ConjugateReturnType(m_matrix.conjugate()); }
SelfAdjointView< const MatrixConjugateReturnType, UpLo > ConjugateReturnType

◆ conjugateIf()

template<typename MatrixType_ , unsigned int UpLo>
template<bool Cond>
std::conditional_t<Cond,ConjugateReturnType,ConstSelfAdjointView> Eigen::SelfAdjointView< MatrixType_, UpLo >::conjugateIf ( ) const
inline
Returns
an expression of the complex conjugate of *this if Cond==true, returns *this otherwise.

Definition at line 207 of file SelfAdjointView.h.

208  {
209  typedef std::conditional_t<Cond,ConjugateReturnType,ConstSelfAdjointView> ReturnType;
210  return ReturnType(m_matrix.template conjugateIf<Cond>());
211  }

◆ diagonal()

template<typename MatrixType_ , unsigned int UpLo>
MatrixType::ConstDiagonalReturnType Eigen::SelfAdjointView< MatrixType_, UpLo >::diagonal ( ) const
inline
Returns
a const expression of the main diagonal of the matrix *this

This method simply returns the diagonal of the nested expression, thus by-passing the SelfAdjointView decorator.

See also
MatrixBase::diagonal(), class Diagonal

Definition at line 243 of file SelfAdjointView.h.

244  {
246  }
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:67

◆ eigenvalues()

template<typename MatrixType , unsigned int UpLo>
SelfAdjointView< MatrixType, UpLo >::EigenvaluesReturnType Eigen::SelfAdjointView< MatrixType, UpLo >::eigenvalues
inline

Computes the eigenvalues of a matrix.

Returns
Column vector containing the eigenvalues.

This is defined in the Eigenvalues module.

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:

VectorXd eivals = ones.selfadjointView<Lower>().eigenvalues();
cout << "The eigenvalues of the 3x3 matrix of ones are:" << endl << eivals << endl;
VectorXcd eivals
static const ConstantReturnType Ones()
EigenvaluesReturnType eigenvalues() const
Computes the eigenvalues of a matrix.
Matrix< double, Dynamic, 1 > VectorXd
DynamicĂ—1 vector of type double.
Definition: Matrix.h:502
Matrix< double, Dynamic, Dynamic > MatrixXd
DynamicĂ—Dynamic matrix of type double.
Definition: Matrix.h:502

Output:

The eigenvalues of the 3x3 matrix of ones are:
-3.09e-16
        0
        3
See also
SelfAdjointEigenSolver::eigenvalues(), MatrixBase::eigenvalues()

Definition at line 90 of file MatrixBaseEigenvalues.h.

91 {
92  PlainObject thisAsMatrix(*this);
93  return SelfAdjointEigenSolver<PlainObject>(thisAsMatrix, false).eigenvalues();
94 }
MatrixType::PlainObject PlainObject

◆ innerStride()

template<typename MatrixType_ , unsigned int UpLo>
EIGEN_CONSTEXPR Index Eigen::SelfAdjointView< MatrixType_, UpLo >::innerStride ( ) const
inline

Definition at line 86 of file SelfAdjointView.h.

86 { return m_matrix.innerStride(); }

◆ ldlt()

template<typename MatrixType , unsigned int UpLo>
const LDLT< typename SelfAdjointView< MatrixType, UpLo >::PlainObject, UpLo > Eigen::SelfAdjointView< MatrixType, UpLo >::ldlt
inline

This is defined in the Cholesky module.

Returns
the Cholesky decomposition with full pivoting without square root of *this
See also
MatrixBase::ldlt()

Definition at line 667 of file LDLT.h.

668 {
669  return LDLT<PlainObject,UpLo>(m_matrix);
670 }

◆ llt()

template<typename MatrixType , unsigned int UpLo>
const LLT< typename SelfAdjointView< MatrixType, UpLo >::PlainObject, UpLo > Eigen::SelfAdjointView< MatrixType, UpLo >::llt
inline

This is defined in the Cholesky module.

#include <Eigen/Cholesky>
Returns
the LLT decomposition of *this
See also
SelfAdjointView::llt()

Definition at line 548 of file LLT.h.

549 {
550  return LLT<PlainObject,UpLo>(m_matrix);
551 }

◆ nestedExpression() [1/2]

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

Definition at line 116 of file SelfAdjointView.h.

116 { return m_matrix; }

◆ nestedExpression() [2/2]

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

Definition at line 114 of file SelfAdjointView.h.

114 { return m_matrix; }

◆ operator*()

template<typename MatrixType_ , unsigned int UpLo>
template<typename OtherDerived >
const Product<SelfAdjointView,OtherDerived> Eigen::SelfAdjointView< MatrixType_, UpLo >::operator* ( const MatrixBase< OtherDerived > &  rhs) const
inline

Efficient triangular matrix times vector/matrix product

Definition at line 122 of file SelfAdjointView.h.

123  {
124  return Product<SelfAdjointView,OtherDerived>(*this, rhs.derived());
125  }

◆ operatorNorm()

template<typename MatrixType , unsigned int UpLo>
SelfAdjointView< MatrixType, UpLo >::RealScalar Eigen::SelfAdjointView< MatrixType, UpLo >::operatorNorm
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:

cout << "The operator norm of the 3x3 matrix of ones is "
<< ones.selfadjointView<Lower>().operatorNorm() << endl;
RealScalar operatorNorm() const
Computes the L2 operator norm.

Output:

The operator norm of the 3x3 matrix of ones is 3
See also
eigenvalues(), MatrixBase::operatorNorm()

Definition at line 153 of file MatrixBaseEigenvalues.h.

154 {
155  return eigenvalues().cwiseAbs().maxCoeff();
156 }
internal::traits< Derived >::Scalar maxCoeff() const
Definition: Redux.h:533

◆ outerStride()

template<typename MatrixType_ , unsigned int UpLo>
EIGEN_CONSTEXPR Index Eigen::SelfAdjointView< MatrixType_, UpLo >::outerStride ( ) const
inline

Definition at line 84 of file SelfAdjointView.h.

84 { return m_matrix.outerStride(); }

◆ rankUpdate() [1/4]

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 
)

Definition at line 63 of file SelfadjointRank2Update.h.

64 {
65  typedef internal::blas_traits<DerivedU> UBlasTraits;
66  typedef typename UBlasTraits::DirectLinearAccessType ActualUType;
67  typedef internal::remove_all_t<ActualUType> ActualUType_;
68  internal::add_const_on_value_type_t<ActualUType> actualU = UBlasTraits::extract(u.derived());
69 
70  typedef internal::blas_traits<DerivedV> VBlasTraits;
71  typedef typename VBlasTraits::DirectLinearAccessType ActualVType;
72  typedef internal::remove_all_t<ActualVType> ActualVType_;
73  internal::add_const_on_value_type_t<ActualVType> actualV = VBlasTraits::extract(v.derived());
74 
75  // If MatrixType is row major, then we use the routine for lower triangular in the upper triangular case and
76  // vice versa, and take the complex conjugate of all coefficients and vector entries.
77 
78  enum { IsRowMajor = (internal::traits<MatrixType>::Flags&RowMajorBit) ? 1 : 0 };
79  Scalar actualAlpha = alpha * UBlasTraits::extractScalarFactor(u.derived())
80  * numext::conj(VBlasTraits::extractScalarFactor(v.derived()));
81  if (IsRowMajor)
82  actualAlpha = numext::conj(actualAlpha);
83 
84  typedef internal::remove_all_t<typename internal::conj_expr_if<int(IsRowMajor) ^ int(UBlasTraits::NeedToConjugate), ActualUType_>::type> UType;
85  typedef internal::remove_all_t<typename internal::conj_expr_if<int(IsRowMajor) ^ int(VBlasTraits::NeedToConjugate), ActualVType_>::type> VType;
86  internal::selfadjoint_rank2_update_selector<Scalar, Index, UType, VType,
87  (IsRowMajor ? int(UpLo==Upper ? Lower : Upper) : UpLo)>
88  ::run(_expression().const_cast_derived().data(),_expression().outerStride(),UType(actualU),VType(actualV),actualAlpha);
89 
90  return *this;
91 }
Array< int, Dynamic, 1 > v
int data[]
const MatrixTypeNestedCleaned & _expression() const
internal::traits< SelfAdjointView >::Scalar Scalar
The type of coefficients in this matrix.
EIGEN_CONSTEXPR Index outerStride() const EIGEN_NOEXCEPT
const unsigned int RowMajorBit
Definition: Constants.h:68
typename remove_all< T >::type remove_all_t
Definition: Meta.h:119
std::conditional<!Cond, const T &, CwiseUnaryOp< scalar_conjugate_op< typename traits< T >::Scalar >, T > > conj_expr_if
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
Derived & const_cast_derived() const
Definition: EigenBase.h:54

◆ rankUpdate() [2/4]

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)

◆ rankUpdate() [3/4]

template<typename MatrixType_ , unsigned int UpLo>
template<typename DerivedU >
SelfAdjointView<MatrixType,UpLo>& Eigen::SelfAdjointView< MatrixType_, UpLo >::rankUpdate ( const MatrixBase< DerivedU > &  u,
const Scalar alpha 
)

Definition at line 126 of file SelfadjointProduct.h.

127 {
128  selfadjoint_product_selector<MatrixType,DerivedU,UpLo>::run(_expression().const_cast_derived(), u.derived(), alpha);
129 
130  return *this;
131 }

◆ rankUpdate() [4/4]

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)

◆ rows()

template<typename MatrixType_ , unsigned int UpLo>
EIGEN_CONSTEXPR Index Eigen::SelfAdjointView< MatrixType_, UpLo >::rows ( void  ) const
inline

Definition at line 80 of file SelfAdjointView.h.

80 { return m_matrix.rows(); }

◆ transpose() [1/2]

template<typename MatrixType_ , unsigned int UpLo>
const ConstTransposeReturnType Eigen::SelfAdjointView< MatrixType_, UpLo >::transpose ( ) const
inline
See also
MatrixBase::transpose() const

Definition at line 232 of file SelfAdjointView.h.

233  {
234  return ConstTransposeReturnType(m_matrix.transpose());
235  }
SelfAdjointView< const typename MatrixType::ConstTransposeReturnType, TransposeMode > ConstTransposeReturnType

◆ transpose() [2/2]

template<typename MatrixType_ , unsigned int UpLo>
template<class Dummy = int>
TransposeReturnType Eigen::SelfAdjointView< MatrixType_, UpLo >::transpose ( std::enable_if_t< Eigen::internal::is_lvalue< MatrixType >::value, Dummy * >  = nullptr)
inline
See also
MatrixBase::transpose()

Definition at line 223 of file SelfAdjointView.h.

224  {
226  return TransposeReturnType(tmp);
227  }
SelfAdjointView< typename MatrixType::TransposeReturnType, TransposeMode > TransposeReturnType
Expression of the transpose of a matrix.
Definition: Transpose.h:56

◆ triangularView()

template<typename MatrixType_ , unsigned int UpLo>
template<unsigned int TriMode>
std::conditional_t<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), TriangularView<MatrixType,TriMode>, TriangularView<typename MatrixType::AdjointReturnType,TriMode> > Eigen::SelfAdjointView< MatrixType_, UpLo >::triangularView ( ) const
inline
Returns
an expression of a triangular view extracted from the current selfadjoint view of a given triangular part

The parameter TriMode can have the following values: Upper, StrictlyUpper, UnitUpper, Lower, StrictlyLower, UnitLower.

If TriMode references the same triangular part than *this, then this method simply return a TriangularView of the nested expression, otherwise, the nested expression is first transposed, thus returning a TriangularView<Transpose<MatrixType>> object.

See also
MatrixBase::triangularView(), class TriangularView

Definition at line 186 of file SelfAdjointView.h.

187  {
188  std::conditional_t<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), MatrixType&, typename MatrixType::ConstTransposeReturnType> tmp1(m_matrix);
189  std::conditional_t<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), MatrixType&, typename MatrixType::AdjointReturnType> tmp2(tmp1);
190  return std::conditional_t<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)),
191  TriangularView<MatrixType,TriMode>,
192  TriangularView<typename MatrixType::AdjointReturnType,TriMode> >(tmp2);
193  }
Matrix< float, 1, Dynamic > MatrixType

Member Data Documentation

◆ m_matrix

template<typename MatrixType_ , unsigned int UpLo>
MatrixTypeNested Eigen::SelfAdjointView< MatrixType_, UpLo >::m_matrix
protected

Definition at line 266 of file SelfAdjointView.h.


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