Public Types | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Friends
TriangularView< _MatrixType, _Mode > Class Template Reference

Base class for triangular part in a matrix. More...

#include <TriangularMatrix.h>

+ Inheritance diagram for TriangularView< _MatrixType, _Mode >:

List of all members.

Public Types

enum  
enum  {
  Mode,
  TransposeMode
}
typedef TriangularBase
< TriangularView
Base
typedef internal::traits
< TriangularView >
::DenseMatrixType 
DenseMatrixType
typedef DenseMatrixType DenseType
typedef internal::traits
< TriangularView >::Index 
Index
typedef _MatrixType MatrixType
typedef DenseMatrixType PlainObject
typedef internal::traits
< TriangularView >::Scalar 
Scalar
typedef internal::traits
< TriangularView >
::StorageKind 
StorageKind

Public Member Functions

void addTo (Dest &dst) const
const TriangularView< const
typename
MatrixType::AdjointReturnType,
TransposeMode
adjoint () const
void applyThisOnTheLeft (Dest &dst) const
void applyThisOnTheRight (Dest &dst) const
template<typename ProductDerived , typename _Lhs , typename _Rhs >
TriangularView< MatrixType,
UpLo > & 
assignProduct (const ProductBase< ProductDerived, _Lhs, _Rhs > &prod, const Scalar &alpha)
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
TriangularView
< MatrixConjugateReturnType,
Mode
conjugate ()
const TriangularView
< MatrixConjugateReturnType,
Mode
conjugate () const
TriangularView< _MatrixType,
_Mode > & 
const_cast_derived () const
const TriangularView
< _MatrixType, _Mode > & 
const_derived () const
void copyCoeff (Index row, Index col, Other &other)
TriangularView< _MatrixType,
_Mode > & 
derived ()
const TriangularView
< _MatrixType, _Mode > & 
derived () const
Scalar determinant () const
void evalTo (Dest &dst) const
void evalTo (MatrixBase< DenseDerived > &other) const
void evalToLazy (MatrixBase< DenseDerived > &other) const
void fill (const Scalar &value)
Index innerStride () const
template<typename OtherDerived >
void lazyAssign (const TriangularBase< OtherDerived > &other)
template<typename OtherDerived >
void lazyAssign (const MatrixBase< OtherDerived > &other)
const MatrixTypeNestedCleanednestedExpression () const
MatrixTypeNestedCleanednestedExpression ()
Scalar operator() (Index row, Index col) const
Scalaroperator() (Index row, Index col)
template<typename OtherDerived >
TriangularProduct< Mode, true,
MatrixType, false,
OtherDerived,
OtherDerived::IsVectorAtCompileTime > 
operator* (const MatrixBase< OtherDerived > &rhs) const
TriangularViewoperator*= (const typename internal::traits< MatrixType >::Scalar &other)
template<typename Other >
TriangularViewoperator+= (const DenseBase< Other > &other)
template<typename ProductDerived , typename Lhs , typename Rhs >
TriangularViewoperator+= (const ProductBase< ProductDerived, Lhs, Rhs > &other)
template<typename ProductDerived >
TriangularViewoperator+= (const ScaledProduct< ProductDerived > &other)
template<typename Other >
TriangularViewoperator-= (const DenseBase< Other > &other)
template<typename ProductDerived , typename Lhs , typename Rhs >
TriangularViewoperator-= (const ProductBase< ProductDerived, Lhs, Rhs > &other)
template<typename ProductDerived >
TriangularViewoperator-= (const ScaledProduct< ProductDerived > &other)
TriangularViewoperator/= (const typename internal::traits< MatrixType >::Scalar &other)
template<typename OtherDerived >
TriangularViewoperator= (const TriangularBase< OtherDerived > &other)
template<typename OtherDerived >
TriangularViewoperator= (const MatrixBase< OtherDerived > &other)
TriangularViewoperator= (const TriangularView &other)
template<typename ProductDerived , typename Lhs , typename Rhs >
TriangularViewoperator= (const ProductBase< ProductDerived, Lhs, Rhs > &other)
template<typename ProductDerived >
TriangularViewoperator= (const ScaledProduct< ProductDerived > &other)
Index outerStride () const
Index rows () const
const SelfAdjointView
< MatrixTypeNestedNonRef, Mode
selfadjointView () const
SelfAdjointView
< MatrixTypeNestedNonRef, Mode
selfadjointView ()
TriangularViewsetConstant (const Scalar &value)
TriangularViewsetOnes ()
TriangularViewsetZero ()
Index size () const
template<int Side, typename Other >
const
internal::triangular_solve_retval
< Side, TriangularView, Other > 
solve (const MatrixBase< Other > &other) const
template<typename Other >
const
internal::triangular_solve_retval
< OnTheLeft, TriangularView,
Other > 
solve (const MatrixBase< Other > &other) const
template<int Side, typename OtherDerived >
void solveInPlace (const MatrixBase< OtherDerived > &other) const
template<typename OtherDerived >
void solveInPlace (const MatrixBase< OtherDerived > &other) const
void subTo (Dest &dst) const
template<typename OtherDerived >
void swap (TriangularBase< OtherDerived > const &other)
template<typename OtherDerived >
void swap (MatrixBase< OtherDerived > const &other)
DenseMatrixType toDenseMatrix () const
TriangularView< Transpose
< MatrixType >, TransposeMode
transpose ()
const TriangularView
< Transpose< MatrixType >
, TransposeMode
transpose () const
 TriangularView (const MatrixType &matrix)

Protected Types

typedef internal::remove_all
< typename
MatrixType::ConjugateReturnType >
::type 
MatrixConjugateReturnType
typedef internal::traits
< TriangularView >
::MatrixTypeNested 
MatrixTypeNested
typedef internal::traits
< TriangularView >
::MatrixTypeNestedCleaned 
MatrixTypeNestedCleaned
typedef internal::traits
< TriangularView >
::MatrixTypeNestedNonRef 
MatrixTypeNestedNonRef

Protected Member Functions

template<typename ProductDerived , typename Lhs , typename Rhs >
TriangularViewassignProduct (const ProductBase< ProductDerived, Lhs, Rhs > &prod, const Scalar &alpha)
void check_coordinates (Index row, Index col) const
void check_coordinates_internal (Index, Index) const

Protected Attributes

MatrixTypeNested m_matrix

Friends

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

Detailed Description

template<typename _MatrixType, unsigned int _Mode>
class Eigen::TriangularView< _MatrixType, _Mode >

Base class for triangular part in a matrix.

Parameters:
MatrixTypethe type of the object in which we are taking the triangular part
Modethe kind of triangular matrix expression to construct. Can be Upper, Lower, UnitUpper, UnitLower, StrictlyUpper, or StrictlyLower. This is in fact a bit field; it must have either Upper or Lower, and additionnaly it may have UnitDiag or ZeroDiag or neither.

This class represents a triangular part of a matrix, not necessarily square. Strictly speaking, for rectangular matrices one should speak of "trapezoid" parts. This class is the return type of MatrixBase::triangularView() and most of the time this is the only way it is used.

See also:
MatrixBase::triangularView()

Member Typedef Documentation

typedef internal::traits<TriangularView>::DenseMatrixType DenseMatrixType
typedef DenseMatrixType DenseType [inherited]
typedef internal::traits<TriangularView>::Index Index
typedef internal::remove_all<typename MatrixType::ConjugateReturnType>::type MatrixConjugateReturnType [protected]
typedef _MatrixType MatrixType
typedef internal::traits<TriangularView>::MatrixTypeNested MatrixTypeNested [protected]
typedef internal::traits<TriangularView>::MatrixTypeNestedCleaned MatrixTypeNestedCleaned [protected]
typedef internal::traits<TriangularView>::MatrixTypeNestedNonRef MatrixTypeNestedNonRef [protected]
typedef internal::traits<TriangularView>::Scalar Scalar
typedef internal::traits<TriangularView>::StorageKind StorageKind

Member Enumeration Documentation

anonymous enum [inherited]
anonymous enum
Enumerator:
Mode 
TransposeMode 

Constructor & Destructor Documentation

TriangularView ( const MatrixType matrix) [inline]

Member Function Documentation

void addTo ( Dest &  dst) const [inline, inherited]
const TriangularView<const typename MatrixType::AdjointReturnType,TransposeMode> adjoint ( ) const [inline]
void applyThisOnTheLeft ( Dest &  dst) const [inline, inherited]
void applyThisOnTheRight ( Dest &  dst) const [inline, inherited]
TriangularView<MatrixType,UpLo>& assignProduct ( const ProductBase< ProductDerived, _Lhs, _Rhs > &  prod,
const Scalar alpha 
)
TriangularView& assignProduct ( const ProductBase< ProductDerived, Lhs, Rhs > &  prod,
const Scalar alpha 
) [inline, protected]
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

References TriangularBase< Derived >::check_coordinates_internal(), and TriangularView< _MatrixType, _Mode >::m_matrix.

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

References TriangularBase< Derived >::check_coordinates_internal(), and TriangularView< _MatrixType, _Mode >::m_matrix.

Index cols ( void  ) const [inline]
TriangularView< _MatrixType, _Mode > & const_cast_derived ( ) const [inline, inherited]
const TriangularView< _MatrixType, _Mode > & const_derived ( ) const [inline, inherited]
void copyCoeff ( Index  row,
Index  col,
Other &  other 
) [inline, inherited]
See also:
MatrixBase::copyCoeff(row,col)
TriangularView< _MatrixType, _Mode > & derived ( ) [inline, inherited]
Returns:
a reference to the derived object
const TriangularView< _MatrixType, _Mode > & derived ( ) const [inline, inherited]
Returns:
a const reference to the derived object
Scalar determinant ( ) const [inline]
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.

void fill ( const Scalar value) [inline]
Index innerStride ( ) const [inline]
void lazyAssign ( const TriangularBase< OtherDerived > &  other)
void lazyAssign ( const MatrixBase< OtherDerived > &  other)

References Eigen::Dynamic.

const MatrixTypeNestedCleaned& nestedExpression ( ) const [inline]
Scalar operator() ( Index  row,
Index  col 
) const [inline, inherited]
Scalar& operator() ( Index  row,
Index  col 
) [inline, inherited]
TriangularProduct<Mode,true,MatrixType,false,OtherDerived, OtherDerived::IsVectorAtCompileTime> operator* ( const MatrixBase< OtherDerived > &  rhs) const [inline]

Efficient triangular matrix times vector/matrix product

References TriangularView< _MatrixType, _Mode >::m_matrix, and TriangularView< _MatrixType, _Mode >::Mode.

TriangularView& operator*= ( const typename internal::traits< MatrixType >::Scalar other) [inline]
TriangularView& operator+= ( const DenseBase< Other > &  other) [inline]
TriangularView& operator+= ( const ProductBase< ProductDerived, Lhs, Rhs > &  other) [inline]
TriangularView& operator+= ( const ScaledProduct< ProductDerived > &  other) [inline]
TriangularView& operator-= ( const DenseBase< Other > &  other) [inline]
TriangularView& operator-= ( const ProductBase< ProductDerived, Lhs, Rhs > &  other) [inline]
TriangularView& operator-= ( const ScaledProduct< ProductDerived > &  other) [inline]
TriangularView& operator/= ( const typename internal::traits< MatrixType >::Scalar other) [inline]
TriangularView< MatrixType, Mode > & operator= ( const TriangularBase< OtherDerived > &  other) [inline]

Assigns a triangular matrix to a triangular part of a dense matrix

References TriangularBase< Derived >::cols(), EigenBase< Derived >::derived(), Eigen::EvalBeforeAssigningBit, and TriangularBase< Derived >::rows().

TriangularView< MatrixType, Mode > & operator= ( const MatrixBase< OtherDerived > &  other) [inline]
TriangularView& operator= ( const TriangularView< _MatrixType, _Mode > &  other) [inline]
TriangularView& operator= ( const ProductBase< ProductDerived, Lhs, Rhs > &  other) [inline]
TriangularView& operator= ( const ScaledProduct< ProductDerived > &  other) [inline]
Index outerStride ( ) const [inline]
Index rows ( void  ) const [inline]
TriangularView& setConstant ( const Scalar value) [inline]
TriangularView& setOnes ( ) [inline]
TriangularView& setZero ( ) [inline]
Index size ( ) const [inline, inherited]
Returns:
the number of coefficients, which is rows()*cols().
See also:
rows(), cols(), SizeAtCompileTime.
const internal::triangular_solve_retval< Side, TriangularView< Derived, Mode >, Other > solve ( const MatrixBase< Other > &  other) const [inline]
Returns:
the product of the inverse of *this with other, *this being triangular.

This function computes the inverse-matrix matrix product inverse(*this) * other if Side==OnTheLeft (the default), or the right-inverse-multiply other * inverse(*this) if Side==OnTheRight.

The matrix *this must be triangular and invertible (i.e., all the coefficients of the diagonal must be non zero). It works as a forward (resp. backward) substitution if *this is an upper (resp. lower) triangular matrix.

Example:

#ifndef _MSC_VER
  #warning deprecated
#endif
/*
Matrix3d m = Matrix3d::Zero();
m.part<Eigen::UpperTriangular>().setOnes();
cout << "Here is the matrix m:" << endl << m << endl;
Matrix3d n = Matrix3d::Ones();
n.part<Eigen::LowerTriangular>() *= 2;
cout << "Here is the matrix n:" << endl << n << endl;
cout << "And now here is m.inverse()*n, taking advantage of the fact that"
        " m is upper-triangular:" << endl
     << m.marked<Eigen::UpperTriangular>().solveTriangular(n);
*/

Output:

This function returns an expression of the inverse-multiply and can works in-place if it is assigned to the same matrix or vector other.

For users coming from BLAS, this function (and more specifically solveInPlace()) offer all the operations supported by the *TRSV and *TRSM BLAS routines.

See also:
TriangularView::solveInPlace()
const internal::triangular_solve_retval<OnTheLeft,TriangularView, Other> solve ( const MatrixBase< Other > &  other) const [inline]
void solveInPlace ( const MatrixBase< OtherDerived > &  _other) const

"in-place" version of TriangularView::solve() where the result is written in other

Warning:
The parameter is only marked 'const' to make the C++ compiler accept a temporary expression here. This function will const_cast it, so constness isn't honored here.

See TriangularView:solve() for the details.

References Eigen::Lower, Eigen::OnTheLeft, Eigen::OnTheRight, Eigen::Upper, and Eigen::ZeroDiag.

void solveInPlace ( const MatrixBase< OtherDerived > &  other) const [inline]
void subTo ( Dest &  dst) const [inline, inherited]
void swap ( TriangularBase< OtherDerived > const &  other) [inline]
void swap ( MatrixBase< OtherDerived > const &  other) [inline]
DenseMatrixType toDenseMatrix ( ) const [inline, inherited]

Friends And Related Function Documentation

TriangularProduct<Mode,false,OtherDerived,OtherDerived::IsVectorAtCompileTime,MatrixType,false> operator* ( const MatrixBase< OtherDerived > &  lhs,
const TriangularView< _MatrixType, _Mode > &  rhs 
) [friend]

Efficient vector/matrix times triangular matrix product


Member Data Documentation


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