#[cfg(feature = "serde-serialize")] use serde::{Deserialize, Serialize}; use alga::general::ComplexField; use crate::allocator::Allocator; use crate::base::{DefaultAllocator, Matrix, MatrixMN, MatrixN, SquareMatrix}; use crate::constraint::{SameNumberOfRows, ShapeConstraint}; use crate::dimension::{Dim, DimAdd, DimSum, DimDiff, DimSub, Dynamic, U1}; use crate::storage::{Storage, StorageMut}; use crate::base::allocator::Reallocator; /// The Cholesky decomposition of a symmetric-definite-positive matrix. #[cfg_attr(feature = "serde-serialize", derive(Serialize, Deserialize))] #[cfg_attr( feature = "serde-serialize", serde(bound(serialize = "DefaultAllocator: Allocator, MatrixN: Serialize")) )] #[cfg_attr( feature = "serde-serialize", serde(bound(deserialize = "DefaultAllocator: Allocator, MatrixN: Deserialize<'de>")) )] #[derive(Clone, Debug)] pub struct Cholesky where DefaultAllocator: Allocator, { chol: MatrixN, } impl Copy for Cholesky where DefaultAllocator: Allocator, MatrixN: Copy, { } impl> Cholesky where DefaultAllocator: Allocator, { /// Attempts to compute the Cholesky decomposition of `matrix`. /// /// Returns `None` if the input matrix is not definite-positive. The input matrix is assumed /// to be symmetric and only the lower-triangular part is read. pub fn new(mut matrix: MatrixN) -> Option { assert!(matrix.is_square(), "The input matrix must be square."); let n = matrix.nrows(); for j in 0..n { for k in 0..j { let factor = unsafe { -*matrix.get_unchecked((j, k)) }; let (mut col_j, col_k) = matrix.columns_range_pair_mut(j, k); let mut col_j = col_j.rows_range_mut(j..); let col_k = col_k.rows_range(j..); col_j.axpy(factor.conjugate(), &col_k, N::one()); } let diag = unsafe { *matrix.get_unchecked((j, j)) }; if !diag.is_zero() { if let Some(denom) = diag.try_sqrt() { unsafe { *matrix.get_unchecked_mut((j, j)) = denom; } let mut col = matrix.slice_range_mut(j + 1.., j); col /= denom; continue; } } // The diagonal element is either zero or its square root could not // be taken (e.g. for negative real numbers). return None; } Some(Cholesky { chol: matrix }) } /// Retrieves the lower-triangular factor of the Cholesky decomposition with its strictly /// upper-triangular part filled with zeros. pub fn unpack(mut self) -> MatrixN { self.chol.fill_upper_triangle(N::zero(), 1); self.chol } /// Retrieves the lower-triangular factor of the Cholesky decomposition, without zeroing-out /// its strict upper-triangular part. /// /// The values of the strict upper-triangular part are garbage and should be ignored by further /// computations. pub fn unpack_dirty(self) -> MatrixN { self.chol } /// Retrieves the lower-triangular factor of the Cholesky decomposition with its strictly /// uppen-triangular part filled with zeros. pub fn l(&self) -> MatrixN { self.chol.lower_triangle() } /// Retrieves the lower-triangular factor of the Cholesky decomposition, without zeroing-out /// its strict upper-triangular part. /// /// This is an allocation-less version of `self.l()`. The values of the strict upper-triangular /// part are garbage and should be ignored by further computations. pub fn l_dirty(&self) -> &MatrixN { &self.chol } /// Solves the system `self * x = b` where `self` is the decomposed matrix and `x` the unknown. /// /// The result is stored on `b`. pub fn solve_mut(&self, b: &mut Matrix) where S2: StorageMut, ShapeConstraint: SameNumberOfRows, { let _ = self.chol.solve_lower_triangular_mut(b); let _ = self.chol.ad_solve_lower_triangular_mut(b); } /// Returns the solution of the system `self * x = b` where `self` is the decomposed matrix and /// `x` the unknown. pub fn solve(&self, b: &Matrix) -> MatrixMN where S2: Storage, DefaultAllocator: Allocator, ShapeConstraint: SameNumberOfRows, { let mut res = b.clone_owned(); self.solve_mut(&mut res); res } /// Computes the inverse of the decomposed matrix. pub fn inverse(&self) -> MatrixN { let shape = self.chol.data.shape(); let mut res = MatrixN::identity_generic(shape.0, shape.1); self.solve_mut(&mut res); res } /// Given the Cholesky decomposition of a matrix `M`, a scalar `sigma` and a vector `v`, /// performs a rank one update such that we end up with the decomposition of `M + sigma * v*v.adjoint()`. pub fn rank_one_update(&mut self, x: &Matrix, sigma: N::RealField) where S2: Storage, DefaultAllocator: Allocator, ShapeConstraint: SameNumberOfRows, { // heavily inspired by Eigen's `llt_rank_update_lower` implementation https://eigen.tuxfamily.org/dox/LLT_8h_source.html let n = x.nrows(); assert_eq!( n, self.chol.nrows(), "The input vector must be of the same size as the factorized matrix." ); let mut x = x.clone_owned(); let mut beta = crate::one::(); for j in 0..n { // updates the diagonal let diag = N::real(unsafe { *self.chol.get_unchecked((j, j)) }); let diag2 = diag * diag; let xj = unsafe { *x.get_unchecked(j) }; let sigma_xj2 = sigma * N::modulus_squared(xj); let gamma = diag2 * beta + sigma_xj2; let new_diag = (diag2 + sigma_xj2 / beta).sqrt(); unsafe { *self.chol.get_unchecked_mut((j, j)) = N::from_real(new_diag) }; beta += sigma_xj2 / diag2; // updates the terms of L let mut xjplus = x.rows_range_mut(j + 1..); let mut col_j = self.chol.slice_range_mut(j + 1.., j); // temp_jplus -= (wj / N::from_real(diag)) * col_j; xjplus.axpy(-xj / N::from_real(diag), &col_j, N::one()); if gamma != crate::zero::() { // col_j = N::from_real(nljj / diag) * col_j + (N::from_real(nljj * sigma / gamma) * N::conjugate(wj)) * temp_jplus; col_j.axpy( N::from_real(new_diag * sigma / gamma) * N::conjugate(xj), &xjplus, N::from_real(new_diag / diag), ); } } } /// Updates the decomposition such that we get the decomposition of a matrix with the given column `c` in the `j`th position. /// Since the matrix is square, an identical row will be added in the `j`th row. pub fn insert_column( self, j: usize, c: &Matrix, ) -> Cholesky> where D: DimAdd, DefaultAllocator: Reallocator> + Reallocator, DimSum, DimSum>, S2: Storage, ShapeConstraint: SameNumberOfRows>, { let n = c.nrows(); assert_eq!( n, self.chol.nrows() + 1, "The new column must have the size of the factored matrix plus one." ); assert!(j < n, "j needs to be within the bound of the new matrix."); // TODO what is the fastest way to produce the new matrix ? let chol= self.chol.clone().insert_column(j, N::zero()).insert_row(j, N::zero()); // TODO see https://en.wikipedia.org/wiki/Cholesky_decomposition#Updating_the_decomposition unimplemented!(); Cholesky { chol } } /// Updates the decomposition such that we get the decomposition of the factored matrix with its `j`th column removed. /// Since the matrix is square, the `j`th row will also be removed. pub fn remove_column( self, j: usize, ) -> Cholesky> where D: DimSub, DefaultAllocator: Reallocator> + Reallocator, DimDiff, DimDiff>, { let n = self.chol.nrows(); assert!(n > 0, "The matrix needs at least one column."); assert!(j < n, "j needs to be within the bound of the matrix."); // TODO what is the fastest way to produce the new matrix ? let mut chol= self.chol.clone().remove_column(j).remove_row(j); // updates the corner let mut corner = chol.slice_range_mut(j.., j..); let colj = self.chol.slice_range(j+1.., j); rank_one_update_helper(&mut corner, &colj, N::real(N::one())); Cholesky { chol } } } impl, S: Storage> SquareMatrix where DefaultAllocator: Allocator, { /// Attempts to compute the Cholesky decomposition of this matrix. /// /// Returns `None` if the input matrix is not definite-positive. The input matrix is assumed /// to be symmetric and only the lower-triangular part is read. pub fn cholesky(self) -> Option> { Cholesky::new(self.into_owned()) } } /// Given the Cholesky decomposition of a matrix `M`, a scalar `sigma` and a vector `v`, /// performs a rank one update such that we end up with the decomposition of `M + sigma * v*v.adjoint()`. fn rank_one_update_helper(chol : &mut Matrix, x: &Matrix, sigma: N::RealField) where N: ComplexField, D: DimSub, R2: Dim, S: StorageMut, S2: Storage, DefaultAllocator: Allocator + Allocator, ShapeConstraint: SameNumberOfRows, { // heavily inspired by Eigen's `llt_rank_update_lower` implementation https://eigen.tuxfamily.org/dox/LLT_8h_source.html let n = x.nrows(); assert_eq!( n, chol.nrows(), "The input vector must be of the same size as the factorized matrix." ); let mut x = x.clone_owned(); let mut beta = crate::one::(); for j in 0..n { // updates the diagonal let diag = N::real(unsafe { *chol.get_unchecked((j, j)) }); let diag2 = diag * diag; let xj = unsafe { *x.get_unchecked(j) }; let sigma_xj2 = sigma * N::modulus_squared(xj); let gamma = diag2 * beta + sigma_xj2; let new_diag = (diag2 + sigma_xj2 / beta).sqrt(); unsafe { *chol.get_unchecked_mut((j, j)) = N::from_real(new_diag) }; beta += sigma_xj2 / diag2; // updates the terms of L let mut xjplus = x.rows_range_mut(j + 1..); let mut col_j = chol.slice_range_mut(j + 1.., j); // temp_jplus -= (wj / N::from_real(diag)) * col_j; xjplus.axpy(-xj / N::from_real(diag), &col_j, N::one()); if gamma != crate::zero::() { // col_j = N::from_real(nljj / diag) * col_j + (N::from_real(nljj * sigma / gamma) * N::conjugate(wj)) * temp_jplus; col_j.axpy( N::from_real(new_diag * sigma / gamma) * N::conjugate(xj), &xjplus, N::from_real(new_diag / diag), ); } } }