Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
TriangularView< _MatrixType, _Mode > Class Template Reference

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

Inherits TriangularBase< Derived >.

Public Member Functions

const TriangularView< const
typename
MatrixType::AdjointReturnType,
TransposeMode > 
adjoint () const
 
Scalar coeff (Index row, Index col) const
 
Scalar & coeffRef (Index row, Index col)
 
TriangularView
< MatrixConjugateReturnType,
Mode > 
conjugate ()
 
const TriangularView
< MatrixConjugateReturnType,
Mode > 
conjugate () const
 
template<typename Other >
void copyCoeff (Index row, Index col, Other &other)
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename DenseDerived >
void evalTo (MatrixBase< DenseDerived > &other) const
 
template<typename DenseDerived >
void evalToLazy (MatrixBase< DenseDerived > &other) const
 
void fill (const Scalar &value)
 
template<typename OtherDerived >
TriangularProduct< Mode, true,
MatrixType, false,
OtherDerived,
OtherDerived::ColsAtCompileTime==1 > 
operator* (const MatrixBase< OtherDerived > &rhs) const
 
TriangularViewoperator*= (const typename internal::traits< MatrixType >::Scalar &other)
 
template<typename Other >
TriangularViewoperator+= (const DenseBase< Other > &other)
 
template<typename Other >
TriangularViewoperator-= (const DenseBase< Other > &other)
 
TriangularViewoperator/= (const typename internal::traits< MatrixType >::Scalar &other)
 
template<typename OtherDerived >
TriangularViewoperator= (const TriangularBase< OtherDerived > &other)
 
TriangularViewsetConstant (const Scalar &value)
 
TriangularViewsetOnes ()
 
TriangularViewsetZero ()
 
Index size () const
 
template<int Side, typename Other >
const
internal::triangular_solve_retval
< Side, TriangularView
< Derived, Mode >, Other > 
solve (const MatrixBase< Other > &other) const
 
template<int Side, typename OtherDerived >
void solveInPlace (const MatrixBase< OtherDerived > &other) const
 
TriangularView< Transpose
< MatrixType >, TransposeMode > 
transpose ()
 
const TriangularView
< Transpose< MatrixType >
, TransposeMode > 
transpose () const
 

Friends

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

Member Function Documentation

const TriangularView<const typename MatrixType::AdjointReturnType,TransposeMode> adjoint ( ) const
inline
Scalar coeff ( Index  row,
Index  col 
) const
inline
See Also
MatrixBase::coeff()
Warning
the coordinates must fit into the referenced triangular part
Scalar& coeffRef ( Index  row,
Index  col 
)
inline
See Also
MatrixBase::coeffRef()
Warning
the coordinates must fit into the referenced triangular part
TriangularView<MatrixConjugateReturnType,Mode> conjugate ( )
inline
const TriangularView<MatrixConjugateReturnType,Mode> conjugate ( ) const
inline
void copyCoeff ( Index  row,
Index  col,
Other &  other 
)
inlineinherited
See Also
MatrixBase::copyCoeff(row,col)
Derived& derived ( )
inlineinherited
Returns
a reference to the derived object

Referenced by IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::analyzePattern(), MatrixBase< Derived >::applyOnTheLeft(), MatrixBase< Derived >::applyOnTheRight(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::applyTranspositionOnTheLeft(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::applyTranspositionOnTheRight(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::bottomRows(), EigenBase< SparseSymmetricPermutationProduct< MatrixType, UpLo > >::cols(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::cols(), IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::compute(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::eval(), IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::factorize(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::indices(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::inverse(), IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::IterativeSolverBase(), RotationBase< Derived, 3 >::operator*(), SparseSelfAdjointView< MatrixType, UpLo >::operator*(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::operator*(), SparseMatrixBase< Derived >::operator*(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::operator*(), Transform< Scalar, Dim, Mode, _Options >::operator*(), Eigen::operator*(), MatrixBase< Derived >::operator*=(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::operator=(), DenseBase< Derived >::operator=(), Transform< Scalar, Dim, Mode, _Options >::operator=(), PlainObjectBase< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > >::operator=(), PlainObjectBase< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > >::PlainObjectBase(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::real(), PlainObjectBase< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > >::resizeLike(), EigenBase< SparseSymmetricPermutationProduct< MatrixType, UpLo > >::rows(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::rows(), SimplicialCholeskyBase< SimplicialLDLT< _MatrixType, _UpLo, _Ordering > >::solve(), IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::solve(), SparseLU< _MatrixType, _OrderingType >::solve(), UmfPackLU< _MatrixType >::solve(), CholmodBase< _MatrixType, _UpLo, CholmodSimplicialLLT< _MatrixType, _UpLo > >::solve(), SuperLUBase< _MatrixType, SuperILU< _MatrixType > >::solve(), SparseMatrix< Scalar, RowMajor >::SparseMatrix(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::toDenseMatrix(), SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::topLeftCorner(), Transform< Scalar, Dim, Mode, _Options >::Transform(), PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::transpose(), and SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::twistedBy().

const Derived& derived ( ) const
inlineinherited
Returns
a const reference to the derived object
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
TriangularProduct<Mode, true, MatrixType, false, OtherDerived, OtherDerived::ColsAtCompileTime==1> operator* ( const MatrixBase< OtherDerived > &  rhs) const
inline

Efficient triangular matrix times vector/matrix product

TriangularView& operator*= ( const typename internal::traits< MatrixType >::Scalar &  other)
inline
TriangularView& operator+= ( const DenseBase< Other > &  other)
inline
TriangularView& operator-= ( const DenseBase< Other > &  other)
inline
TriangularView& operator/= ( const typename internal::traits< MatrixType >::Scalar &  other)
inline
See Also
MatrixBase::operator/=()
TriangularView& operator= ( const TriangularBase< OtherDerived > &  other)

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

Index size ( ) const
inlineinherited
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
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()
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.

TriangularView<Transpose<MatrixType>,TransposeMode> transpose ( )
inline
const TriangularView<Transpose<MatrixType>,TransposeMode> transpose ( ) const
inline

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


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