aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Core/TriangularMatrix.h
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2015-10-09 12:10:58 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2015-10-09 12:10:58 +0200
commit2632b3446cf687e9e3feff2850d53f7928837474 (patch)
treea00410e995ad3d68fd42092afa1958ce352ed619 /Eigen/src/Core/TriangularMatrix.h
parent1429daf85087a6c913a8274fd4c6827f9eb57aef (diff)
Improve documentation of TriangularView.
Diffstat (limited to 'Eigen/src/Core/TriangularMatrix.h')
-rw-r--r--Eigen/src/Core/TriangularMatrix.h75
1 files changed, 66 insertions, 9 deletions
diff --git a/Eigen/src/Core/TriangularMatrix.h b/Eigen/src/Core/TriangularMatrix.h
index e9b34ebdf..438dd4dc9 100644
--- a/Eigen/src/Core/TriangularMatrix.h
+++ b/Eigen/src/Core/TriangularMatrix.h
@@ -222,18 +222,23 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularView
TriangularView& operator=(const TriangularView &other)
{ return Base::operator=(other); }
+ /** \copydoc EigenBase::rows() */
EIGEN_DEVICE_FUNC
inline Index rows() const { return m_matrix.rows(); }
+ /** \copydoc EigenBase::cols() */
EIGEN_DEVICE_FUNC
inline Index cols() const { return m_matrix.cols(); }
+ /** \returns a const reference to the nested expression */
EIGEN_DEVICE_FUNC
const NestedExpression& nestedExpression() const { return m_matrix; }
+
+ /** \returns a reference to the nested expression */
EIGEN_DEVICE_FUNC
NestedExpression& nestedExpression() { return *const_cast<NestedExpression*>(&m_matrix); }
- /** \sa MatrixBase::conjugate() const */
typedef TriangularView<const MatrixConjugateReturnType,Mode> ConjugateReturnType;
+ /** \sa MatrixBase::conjugate() const */
EIGEN_DEVICE_FUNC
inline const ConjugateReturnType conjugate() const
{ return ConjugateReturnType(m_matrix.conjugate()); }
@@ -279,19 +284,28 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularView
using Base::solve;
#endif
+ /** \returns a selfadjoint view of the referenced triangular part which must be either \c #Upper or \c #Lower.
+ *
+ * This is a shortcut for \code this->nestedExpression().selfadjointView<(*this)::Mode>() \endcode
+ * \sa MatrixBase::selfadjointView() */
EIGEN_DEVICE_FUNC
- const SelfAdjointView<MatrixTypeNestedNonRef,Mode> selfadjointView() const
+ SelfAdjointView<MatrixTypeNestedNonRef,Mode> selfadjointView()
{
- EIGEN_STATIC_ASSERT((Mode&UnitDiag)==0,PROGRAMMING_ERROR);
+ EIGEN_STATIC_ASSERT((Mode&(UnitDiag|ZeroDiag))==0,PROGRAMMING_ERROR);
return SelfAdjointView<MatrixTypeNestedNonRef,Mode>(m_matrix);
}
+
+ /** This is the const version of selfadjointView() */
EIGEN_DEVICE_FUNC
- SelfAdjointView<MatrixTypeNestedNonRef,Mode> selfadjointView()
+ const SelfAdjointView<MatrixTypeNestedNonRef,Mode> selfadjointView() const
{
- EIGEN_STATIC_ASSERT((Mode&UnitDiag)==0,PROGRAMMING_ERROR);
+ EIGEN_STATIC_ASSERT((Mode&(UnitDiag|ZeroDiag))==0,PROGRAMMING_ERROR);
return SelfAdjointView<MatrixTypeNestedNonRef,Mode>(m_matrix);
}
+
+ /** \returns the determinant of the triangular matrix
+ * \sa MatrixBase::determinant() */
EIGEN_DEVICE_FUNC
Scalar determinant() const
{
@@ -341,12 +355,16 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularViewImpl<_Mat
Flags = internal::traits<TriangularViewType>::Flags
};
+ /** \returns the outer-stride of the underlying dense matrix
+ * \sa DenseCoeffsBase::outerStride() */
EIGEN_DEVICE_FUNC
inline Index outerStride() const { return derived().nestedExpression().outerStride(); }
+ /** \returns the inner-stride of the underlying dense matrix
+ * \sa DenseCoeffsBase::innerStride() */
EIGEN_DEVICE_FUNC
inline Index innerStride() const { return derived().nestedExpression().innerStride(); }
- /** \sa MatrixBase::operator+=() */
+ /** \sa MatrixBase::operator+=() */
template<typename Other>
EIGEN_DEVICE_FUNC
TriangularViewType& operator+=(const DenseBase<Other>& other) {
@@ -364,7 +382,7 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularViewImpl<_Mat
/** \sa MatrixBase::operator*=() */
EIGEN_DEVICE_FUNC
TriangularViewType& operator*=(const typename internal::traits<MatrixType>::Scalar& other) { return *this = derived().nestedExpression() * other; }
- /** \sa MatrixBase::operator/=() */
+ /** \sa DenseBase::operator/=() */
EIGEN_DEVICE_FUNC
TriangularViewType& operator/=(const typename internal::traits<MatrixType>::Scalar& other) { return *this = derived().nestedExpression() / other; }
@@ -408,21 +426,26 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularViewImpl<_Mat
EIGEN_DEVICE_FUNC
TriangularViewType& operator=(const TriangularBase<OtherDerived>& other);
+ /** Shortcut for\code *this = other.other.triangularView<(*this)::Mode>() \endcode */
template<typename OtherDerived>
EIGEN_DEVICE_FUNC
TriangularViewType& operator=(const MatrixBase<OtherDerived>& other);
+#ifndef EIGEN_PARSED_BY_DOXYGEN
EIGEN_DEVICE_FUNC
TriangularViewType& operator=(const TriangularViewImpl& other)
{ return *this = other.derived().nestedExpression(); }
+ /** \deprecated */
template<typename OtherDerived>
EIGEN_DEVICE_FUNC
void lazyAssign(const TriangularBase<OtherDerived>& other);
+ /** \deprecated */
template<typename OtherDerived>
EIGEN_DEVICE_FUNC
- void lazyAssign(const MatrixBase<OtherDerived>& other);
+ void lazyAssign(const MatrixBase<OtherDerived>& other);
+#endif
/** Efficient triangular matrix times vector/matrix product */
template<typename OtherDerived>
@@ -442,11 +465,39 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularViewImpl<_Mat
return Product<OtherDerived,TriangularViewType>(lhs.derived(),rhs.derived());
}
+ /** \returns the product of the inverse of \c *this with \a other, \a *this being triangular.
+ *
+ * This function computes the inverse-matrix matrix product inverse(\c *this) * \a other if
+ * \a Side==OnTheLeft (the default), or the right-inverse-multiply \a other * inverse(\c *this) if
+ * \a Side==OnTheRight.
+ *
+ * The matrix \c *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 \c *this
+ * is an upper (resp. lower) triangular matrix.
+ *
+ * Example: \include Triangular_solve.cpp
+ * Output: \verbinclude Triangular_solve.out
+ *
+ * This function returns an expression of the inverse-multiply and can works in-place if it is assigned
+ * to the same matrix or vector \a other.
+ *
+ * For users coming from BLAS, this function (and more specifically solveInPlace()) offer
+ * all the operations supported by the \c *TRSV and \c *TRSM BLAS routines.
+ *
+ * \sa TriangularView::solveInPlace()
+ */
template<int Side, typename Other>
EIGEN_DEVICE_FUNC
inline const internal::triangular_solve_retval<Side,TriangularViewType, Other>
solve(const MatrixBase<Other>& other) const;
+ /** "in-place" version of TriangularView::solve() where the result is written in \a 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.
+ */
template<int Side, typename OtherDerived>
EIGEN_DEVICE_FUNC
void solveInPlace(const MatrixBase<OtherDerived>& other) const;
@@ -456,15 +507,21 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularViewImpl<_Mat
void solveInPlace(const MatrixBase<OtherDerived>& other) const
{ return solveInPlace<OnTheLeft>(other); }
+ /** Swaps the coefficients of the common triangular parts of two matrices */
template<typename OtherDerived>
EIGEN_DEVICE_FUNC
+#ifdef EIGEN_PARSED_BY_DOXYGEN
+ void swap(TriangularBase<OtherDerived> &other)
+#else
void swap(TriangularBase<OtherDerived> const & other)
+#endif
{
EIGEN_STATIC_ASSERT_LVALUE(OtherDerived);
call_assignment(derived(), other.const_cast_derived(), internal::swap_assign_op<Scalar>());
}
- // TODO: this overload is ambiguous and it should be deprecated (Gael)
+ /** \deprecated
+ * Shortcut for \code (*this).swap(other.triangularView<(*this)::Mode>()) \endcode */
template<typename OtherDerived>
EIGEN_DEVICE_FUNC
void swap(MatrixBase<OtherDerived> const & other)