diff options
author | Gael Guennebaud <g.gael@free.fr> | 2008-07-19 22:59:05 +0000 |
---|---|---|
committer | Gael Guennebaud <g.gael@free.fr> | 2008-07-19 22:59:05 +0000 |
commit | 269f683902231020a910fd4e2c1b74554183e2c8 (patch) | |
tree | f95a531cd1792f1df2122a3879217adc8eb9ec18 /Eigen/src/Cholesky | |
parent | 6e2c53e0562785f2f8c72a67dd8c40451cf2a319 (diff) |
Add cholesky's members to MatrixBase
Various documentation improvements including new snippets (AngleAxis and Cholesky)
Diffstat (limited to 'Eigen/src/Cholesky')
-rw-r--r-- | Eigen/src/Cholesky/Cholesky.h | 17 | ||||
-rw-r--r-- | Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h | 13 |
2 files changed, 26 insertions, 4 deletions
diff --git a/Eigen/src/Cholesky/Cholesky.h b/Eigen/src/Cholesky/Cholesky.h index 8d9b6e8d8..c1f05d768 100644 --- a/Eigen/src/Cholesky/Cholesky.h +++ b/Eigen/src/Cholesky/Cholesky.h @@ -66,7 +66,7 @@ template<typename MatrixType> class Cholesky bool isPositiveDefinite(void) const { return m_isPositiveDefinite; } template<typename Derived> - typename Derived::Eval solve(MatrixBase<Derived> &b); + typename Derived::Eval solve(const MatrixBase<Derived> &b) const; void compute(const MatrixType& matrix); @@ -110,10 +110,14 @@ void Cholesky<MatrixType>::compute(const MatrixType& a) /** \returns the solution of A x = \a b using the current decomposition of A. * In other words, it returns \code A^-1 b \endcode computing * \code L^-* L^1 b \endcode from right to left. + * + * Example: \include Cholesky_solve.cpp + * Output: \verbinclude Cholesky_solve.out + * */ template<typename MatrixType> template<typename Derived> -typename Derived::Eval Cholesky<MatrixType>::solve(MatrixBase<Derived> &b) +typename Derived::Eval Cholesky<MatrixType>::solve(const MatrixBase<Derived> &b) const { const int size = m_matrix.rows(); ei_assert(size==b.size()); @@ -121,5 +125,14 @@ typename Derived::Eval Cholesky<MatrixType>::solve(MatrixBase<Derived> &b) return m_matrix.adjoint().template extract<Upper>().inverseProduct(matrixL().inverseProduct(b)); } +/** \cholesky_module + * \returns the Cholesky decomposition of \c *this + */ +template<typename Derived> +inline const Cholesky<typename ei_eval<Derived>::type> +MatrixBase<Derived>::cholesky() const +{ + return Cholesky<typename ei_eval<Derived>::type>(derived()); +} #endif // EIGEN_CHOLESKY_H diff --git a/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h b/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h index 8905385cc..2d85f78db 100644 --- a/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h +++ b/Eigen/src/Cholesky/CholeskyWithoutSquareRoot.h @@ -77,7 +77,7 @@ template<typename MatrixType> class CholeskyWithoutSquareRoot } template<typename Derived> - typename Derived::Eval solve(MatrixBase<Derived> &b); + typename Derived::Eval solve(const MatrixBase<Derived> &b) const; void compute(const MatrixType& matrix); @@ -127,7 +127,7 @@ void CholeskyWithoutSquareRoot<MatrixType>::compute(const MatrixType& a) */ template<typename MatrixType> template<typename Derived> -typename Derived::Eval CholeskyWithoutSquareRoot<MatrixType>::solve(MatrixBase<Derived> &vecB) +typename Derived::Eval CholeskyWithoutSquareRoot<MatrixType>::solve(const MatrixBase<Derived> &vecB) const { const int size = m_matrix.rows(); ei_assert(size==vecB.size()); @@ -140,5 +140,14 @@ typename Derived::Eval CholeskyWithoutSquareRoot<MatrixType>::solve(MatrixBase<D ); } +/** \cholesky_module + * \returns the Cholesky decomposition without square root of \c *this + */ +template<typename Derived> +inline const CholeskyWithoutSquareRoot<typename ei_eval<Derived>::type> +MatrixBase<Derived>::choleskyNoSqrt() const +{ + return derived(); +} #endif // EIGEN_CHOLESKY_WITHOUT_SQUARE_ROOT_H |