diff --git a/src/cholesky.rs b/src/cholesky.rs index 26bb2f83..696b5254 100644 --- a/src/cholesky.rs +++ b/src/cholesky.rs @@ -1,6 +1,36 @@ -//! Cholesky decomposition +//! Cholesky decomposition of Hermitian (or real symmetric) positive definite matrices //! -//! https://en.wikipedia.org/wiki/Cholesky_decomposition +//! See the [Wikipedia page about Cholesky +//! decomposition](https://en.wikipedia.org/wiki/Cholesky_decomposition) for +//! more information. +//! +//! # Example +//! +//! Calculate `L` in the Cholesky decomposition `A = L * L^H`, where `A` is a +//! Hermitian (or real symmetric) positive definite matrix: +//! +//! ``` +//! #[macro_use] +//! extern crate ndarray; +//! extern crate ndarray_linalg; +//! +//! use ndarray::prelude::*; +//! use ndarray_linalg::{CholeskyInto, UPLO}; +//! # fn main() { +//! +//! let a: Array2 = array![ +//! [ 4., 12., -16.], +//! [ 12., 37., -43.], +//! [-16., -43., 98.] +//! ]; +//! let lower = a.cholesky_into(UPLO::Lower).unwrap(); +//! assert!(lower.all_close(&array![ +//! [ 2., 0., 0.], +//! [ 6., 1., 0.], +//! [-8., 5., 3.] +//! ], 1e-9)); +//! # } +//! ``` use ndarray::*; @@ -12,19 +42,44 @@ use super::types::*; pub use lapack_traits::UPLO; -/// Cholesky decomposition of matrix reference +/// Cholesky decomposition of Hermitian (or real symmetric) positive definite matrix reference pub trait Cholesky { type Output; + /// Computes the Cholesky decomposition of the Hermitian (or real + /// symmetric) positive definite matrix. + /// + /// If the argument is `UPLO::Upper`, then computes the decomposition `A = + /// U^H * U` using the upper triangular portion of `A` and returns `U`. + /// Otherwise, if the argument is `UPLO::Lower`, computes the decomposition + /// `A = L * L^H` using the lower triangular portion of `A` and returns + /// `L`. fn cholesky(&self, UPLO) -> Result; } -/// Cholesky decomposition +/// Cholesky decomposition of Hermitian (or real symmetric) positive definite matrix pub trait CholeskyInto: Sized { + /// Computes the Cholesky decomposition of the Hermitian (or real + /// symmetric) positive definite matrix. + /// + /// If the argument is `UPLO::Upper`, then computes the decomposition `A = + /// U^H * U` using the upper triangular portion of `A` and returns `U`. + /// Otherwise, if the argument is `UPLO::Lower`, computes the decomposition + /// `A = L * L^H` using the lower triangular portion of `A` and returns + /// `L`. fn cholesky_into(self, UPLO) -> Result; } -/// Cholesky decomposition of mutable reference of matrix +/// Cholesky decomposition of Hermitian (or real symmetric) positive definite mutable reference of matrix pub trait CholeskyMut { + /// Computes the Cholesky decomposition of the Hermitian (or real + /// symmetric) positive definite matrix, storing the result in `self` and + /// returning it. + /// + /// If the argument is `UPLO::Upper`, then computes the decomposition `A = + /// U^H * U` using the upper triangular portion of `A` and returns `U`. + /// Otherwise, if the argument is `UPLO::Lower`, computes the decomposition + /// `A = L * L^H` using the lower triangular portion of `A` and returns + /// `L`. fn cholesky_mut(&mut self, UPLO) -> Result<&mut Self>; } diff --git a/src/solveh.rs b/src/solveh.rs index c09afd5e..291b251c 100644 --- a/src/solveh.rs +++ b/src/solveh.rs @@ -88,7 +88,7 @@ pub trait SolveH { } /// Represents the Bunch–Kaufman factorization of a Hermitian (or real -/// symmetric) matrix as `A = P * U * D * U^T * P^T`. +/// symmetric) matrix as `A = P * U * D * U^H * P^T`. pub struct FactorizedH { pub a: ArrayBase, pub ipiv: Pivot,