use crate::storage::Storage; use crate::{ Allocator, Bidiagonal, Cholesky, ColPivQR, ComplexField, DefaultAllocator, Dim, DimDiff, DimMin, DimMinimum, DimSub, FullPivLU, Hessenberg, Matrix, RealField, Schur, SymmetricEigen, SymmetricTridiagonal, LU, QR, SVD, U1, UDU, }; /// # Rectangular matrix decomposition /// /// This section contains the methods for computing some common decompositions of rectangular /// matrices with real or complex components. The following are currently supported: /// /// | Decomposition | Factors | Details | /// | -------------------------|---------------------|--------------| /// | QR | `Q * R` | `Q` is an unitary matrix, and `R` is upper-triangular. | /// | QR with column pivoting | `Q * R * P⁻¹` | `Q` is an unitary matrix, and `R` is upper-triangular. `P` is a permutation matrix. | /// | LU with partial pivoting | `P⁻¹ * L * U` | `L` is lower-triangular with a diagonal filled with `1` and `U` is upper-triangular. `P` is a permutation matrix. | /// | LU with full pivoting | `P⁻¹ * L * U * Q⁻¹` | `L` is lower-triangular with a diagonal filled with `1` and `U` is upper-triangular. `P` and `Q` are permutation matrices. | /// | SVD | `U * Σ * Vᵀ` | `U` and `V` are two orthogonal matrices and `Σ` is a diagonal matrix containing the singular values. | impl> Matrix { /// Computes the bidiagonalization using householder reflections. pub fn bidiagonalize(self) -> Bidiagonal where R: DimMin, DimMinimum: DimSub, DefaultAllocator: Allocator + Allocator + Allocator + Allocator> + Allocator, U1>>, { Bidiagonal::new(self.into_owned()) } /// Computes the LU decomposition with full pivoting of `matrix`. /// /// This effectively computes `P, L, U, Q` such that `P * matrix * Q = LU`. pub fn full_piv_lu(self) -> FullPivLU where R: DimMin, DefaultAllocator: Allocator + Allocator<(usize, usize), DimMinimum>, { FullPivLU::new(self.into_owned()) } /// Computes the LU decomposition with partial (row) pivoting of `matrix`. pub fn lu(self) -> LU where R: DimMin, DefaultAllocator: Allocator + Allocator<(usize, usize), DimMinimum>, { LU::new(self.into_owned()) } /// Computes the QR decomposition of this matrix. pub fn qr(self) -> QR where R: DimMin, DefaultAllocator: Allocator + Allocator + Allocator>, { QR::new(self.into_owned()) } /// Computes the QR decomposition (with column pivoting) of this matrix. pub fn col_piv_qr(self) -> ColPivQR where R: DimMin, DefaultAllocator: Allocator + Allocator + Allocator> + Allocator<(usize, usize), DimMinimum>, { ColPivQR::new(self.into_owned()) } /// Computes the Singular Value Decomposition using implicit shift. pub fn svd(self, compute_u: bool, compute_v: bool) -> SVD where R: DimMin, DimMinimum: DimSub, // for Bidiagonal. DefaultAllocator: Allocator + Allocator + Allocator + Allocator, U1>> + Allocator, C> + Allocator> + Allocator> + Allocator> + Allocator, U1>>, { SVD::new(self.into_owned(), compute_u, compute_v) } /// Attempts to compute the Singular Value Decomposition of `matrix` using implicit shift. /// /// # Arguments /// /// * `compute_u` − set this to `true` to enable the computation of left-singular vectors. /// * `compute_v` − set this to `true` to enable the computation of right-singular vectors. /// * `eps` − tolerance used to determine when a value converged to 0. /// * `max_niter` − maximum total number of iterations performed by the algorithm. If this /// number of iteration is exceeded, `None` is returned. If `niter == 0`, then the algorithm /// continues indefinitely until convergence. pub fn try_svd( self, compute_u: bool, compute_v: bool, eps: N::RealField, max_niter: usize, ) -> Option> where R: DimMin, DimMinimum: DimSub, // for Bidiagonal. DefaultAllocator: Allocator + Allocator + Allocator + Allocator, U1>> + Allocator, C> + Allocator> + Allocator> + Allocator> + Allocator, U1>>, { SVD::try_new(self.into_owned(), compute_u, compute_v, eps, max_niter) } } /// # Square matrix decomposition /// /// This section contains the methods for computing some common decompositions of square /// matrices with real or complex components. The following are currently supported: /// /// | Decomposition | Factors | Details | /// | -------------------------|---------------------------|--------------| /// | Hessenberg | `Q * H * Qᵀ` | `Q` is a unitary matrix and `H` an upper-Hessenberg matrix. | /// | Cholesky | `L * Lᵀ` | `L` is a lower-triangular matrix. | /// | UDU | `U * D * Uᵀ` | `U` is a upper-triangular matrix, and `D` a diagonal matrix. | /// | Schur decomposition | `Q * T * Qᵀ` | `Q` is an unitary matrix and `T` a quasi-upper-triangular matrix. | /// | Symmetric eigendecomposition | `Q ~ Λ ~ Qᵀ` | `Q` is an unitary matrix, and `Λ` is a real diagonal matrix. | /// | Symmetric tridiagonalization | `Q ~ T ~ Qᵀ` | `Q` is an unitary matrix, and `T` is a tridiagonal matrix. | impl> Matrix { /// 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> where DefaultAllocator: Allocator, { Cholesky::new(self.into_owned()) } /// Attempts to compute the UDU decomposition of this matrix. /// /// The input matrix `self` is assumed to be symmetric and this decomposition will only read /// the upper-triangular part of `self`. pub fn udu(self) -> Option> where N: RealField, DefaultAllocator: Allocator + Allocator, { UDU::new(self.into_owned()) } /// Computes the Hessenberg decomposition of this matrix using householder reflections. pub fn hessenberg(self) -> Hessenberg where D: DimSub, DefaultAllocator: Allocator + Allocator + Allocator>, { Hessenberg::new(self.into_owned()) } /// Computes the Schur decomposition of a square matrix. pub fn schur(self) -> Schur where D: DimSub, // For Hessenberg. DefaultAllocator: Allocator> + Allocator> + Allocator + Allocator, { Schur::new(self.into_owned()) } /// Attempts to compute the Schur decomposition of a square matrix. /// /// If only eigenvalues are needed, it is more efficient to call the matrix method /// `.eigenvalues()` instead. /// /// # Arguments /// /// * `eps` − tolerance used to determine when a value converged to 0. /// * `max_niter` − maximum total number of iterations performed by the algorithm. If this /// number of iteration is exceeded, `None` is returned. If `niter == 0`, then the algorithm /// continues indefinitely until convergence. pub fn try_schur(self, eps: N::RealField, max_niter: usize) -> Option> where D: DimSub, // For Hessenberg. DefaultAllocator: Allocator> + Allocator> + Allocator + Allocator, { Schur::try_new(self.into_owned(), eps, max_niter) } /// Computes the eigendecomposition of this symmetric matrix. /// /// Only the lower-triangular part (including the diagonal) of `m` is read. pub fn symmetric_eigen(self) -> SymmetricEigen where D: DimSub, DefaultAllocator: Allocator + Allocator> + Allocator + Allocator>, { SymmetricEigen::new(self.into_owned()) } /// Computes the eigendecomposition of the given symmetric matrix with user-specified /// convergence parameters. /// /// Only the lower-triangular part (including the diagonal) of `m` is read. /// /// # Arguments /// /// * `eps` − tolerance used to determine when a value converged to 0. /// * `max_niter` − maximum total number of iterations performed by the algorithm. If this /// number of iteration is exceeded, `None` is returned. If `niter == 0`, then the algorithm /// continues indefinitely until convergence. pub fn try_symmetric_eigen( self, eps: N::RealField, max_niter: usize, ) -> Option> where D: DimSub, DefaultAllocator: Allocator + Allocator> + Allocator + Allocator>, { SymmetricEigen::try_new(self.into_owned(), eps, max_niter) } /// Computes the tridiagonalization of this symmetric matrix. /// /// Only the lower-triangular part (including the diagonal) of `m` is read. pub fn symmetric_tridiagonalize(self) -> SymmetricTridiagonal where D: DimSub, DefaultAllocator: Allocator + Allocator>, { SymmetricTridiagonal::new(self.into_owned()) } }