Eigne Matrix Class Library

Dependents:   Eigen_test Odometry_test AttitudeEstimation_usingTicker MPU9250_Quaternion_Binary_Serial ... more

Embed: (wiki syntax)

« Back to documentation index

TriangularView< _MatrixType, _Mode > Class Template Reference

TriangularView< _MatrixType, _Mode > Class Template Reference
[Core module]

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

#include <TriangularMatrix.h>

Inherits Eigen::TriangularBase< TriangularView< _MatrixType, _Mode > >.

Public Member Functions

Index rows () const
Index cols () const
template<typename Other >
TriangularViewoperator+= (const DenseBase< Other > &other)
template<typename Other >
TriangularViewoperator-= (const DenseBase< Other > &other)
TriangularViewoperator*= (const typename internal::traits< MatrixType >::Scalar &other)
TriangularViewoperator/= (const typename internal::traits< MatrixType >::Scalar &other)
void fill (const Scalar &value)
TriangularViewsetConstant (const Scalar &value)
TriangularViewsetZero ()
TriangularViewsetOnes ()
Scalar coeff (Index row, Index col) const
Scalar & coeffRef (Index row, Index col)
template<typename OtherDerived >
TriangularViewoperator= (const TriangularBase< OtherDerived > &other)
 Assigns a triangular matrix to a triangular part of a dense matrix.
TriangularView
< MatrixConjugateReturnType,
Mode > 
conjugate ()
const TriangularView
< MatrixConjugateReturnType,
Mode > 
conjugate () const
const TriangularView< const
typename
MatrixType::AdjointReturnType,
TransposeMode > 
adjoint () const
TriangularView< Transpose
< MatrixType >, TransposeMode > 
transpose ()
const TriangularView
< Transpose< MatrixType >
, TransposeMode > 
transpose () const
template<typename OtherDerived >
TriangularProduct< Mode, true,
MatrixType, false,
OtherDerived,
OtherDerived::ColsAtCompileTime==1 > 
operator* (const MatrixBase< OtherDerived > &rhs) const
 Efficient triangular matrix times vector/matrix product.
template<int Side, typename Other >
const
internal::triangular_solve_retval
< Side, TriangularView, Other > 
solve (const MatrixBase< Other > &other) const
template<int Side, typename OtherDerived >
void solveInPlace (const MatrixBase< OtherDerived > &other) const
 "in-place" version of TriangularView::solve() where the result is written in other
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 >
TriangularProduct< Mode, false,
OtherDerived,
OtherDerived::RowsAtCompileTime==1,
MatrixType, false > 
operator* (const MatrixBase< OtherDerived > &lhs, const TriangularView &rhs)
 Efficient vector/matrix times triangular matrix product.

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

Definition at line 158 of file TriangularMatrix.h.


Member Function Documentation

const TriangularView<const typename MatrixType::AdjointReturnType,TransposeMode> adjoint (  ) const
See also:
MatrixBase::adjoint() const

Definition at line 264 of file TriangularMatrix.h.

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

Definition at line 222 of file TriangularMatrix.h.

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

Definition at line 231 of file TriangularMatrix.h.

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

Reimplemented from EigenBase< Derived >.

Definition at line 196 of file TriangularMatrix.h.

TriangularView<MatrixConjugateReturnType,Mode> conjugate (  )
See also:
MatrixBase::conjugate()

Definition at line 257 of file TriangularMatrix.h.

const TriangularView<MatrixConjugateReturnType,Mode> conjugate (  ) const
See also:
MatrixBase::conjugate() const

Definition at line 260 of file TriangularMatrix.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.

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.

void fill ( const Scalar &  value )
See also:
MatrixBase::fill()

Definition at line 210 of file TriangularMatrix.h.

TriangularProduct<Mode, true, MatrixType, false, OtherDerived, OtherDerived::ColsAtCompileTime==1> operator* ( const MatrixBase< OtherDerived > &  rhs ) const

Efficient triangular matrix times vector/matrix product.

Definition at line 282 of file TriangularMatrix.h.

TriangularView& operator*= ( const typename internal::traits< MatrixType >::Scalar &  other )
See also:
MatrixBase::operator*=()

Definition at line 205 of file TriangularMatrix.h.

TriangularView& operator+= ( const DenseBase< Other > &  other )
See also:
MatrixBase::operator+=()

Definition at line 201 of file TriangularMatrix.h.

TriangularView& operator-= ( const DenseBase< Other > &  other )
See also:
MatrixBase::operator-=()

Definition at line 203 of file TriangularMatrix.h.

TriangularView& operator/= ( const typename internal::traits< MatrixType >::Scalar &  other )
See also:
MatrixBase::operator/=()

Definition at line 207 of file TriangularMatrix.h.

TriangularView< MatrixType, Mode > & operator= ( const TriangularBase< OtherDerived > &  other )

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

Definition at line 640 of file TriangularMatrix.h.

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

Reimplemented from EigenBase< Derived >.

Definition at line 195 of file TriangularMatrix.h.

TriangularView& setConstant ( const Scalar &  value )
See also:
MatrixBase::setConstant()

Definition at line 212 of file TriangularMatrix.h.

TriangularView& setOnes (  )
See also:
MatrixBase::setOnes()

Definition at line 217 of file TriangularMatrix.h.

TriangularView& setZero (  )
See also:
MatrixBase::setZero()

Definition at line 215 of file TriangularMatrix.h.

const internal::triangular_solve_retval< Side, TriangularView< Derived, Mode >, Other > solve ( const MatrixBase< Other > &  other ) const
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:

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

Definition at line 216 of file SolveTriangular.h.

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.

Definition at line 174 of file SolveTriangular.h.

const TriangularView<Transpose<MatrixType>,TransposeMode> transpose (  ) const
See also:
MatrixBase::transpose() const

Definition at line 274 of file TriangularMatrix.h.

TriangularView<Transpose<MatrixType>,TransposeMode> transpose (  )
See also:
MatrixBase::transpose()

Definition at line 268 of file TriangularMatrix.h.


Friends And Related Function Documentation

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

Efficient vector/matrix times triangular matrix product.

Definition at line 292 of file TriangularMatrix.h.