1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364
use crate::storage::Storage;
use crate::{
Allocator, Bidiagonal, Cholesky, ColPivQR, ComplexField, DefaultAllocator, Dim, DimDiff,
DimMin, DimMinimum, DimSub, FullPivLU, Hessenberg, Matrix, OMatrix, 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. |
/// | Polar (Left Polar) | `P' * U` | `U` is semi-unitary/unitary and `P'` is a positive semi-definite Hermitian Matrix
impl<T: ComplexField, R: Dim, C: Dim, S: Storage<T, R, C>> Matrix<T, R, C, S> {
/// Computes the bidiagonalization using householder reflections.
pub fn bidiagonalize(self) -> Bidiagonal<T, R, C>
where
R: DimMin<C>,
DimMinimum<R, C>: DimSub<U1>,
DefaultAllocator: Allocator<R, C>
+ Allocator<C>
+ Allocator<R>
+ Allocator<DimMinimum<R, C>>
+ Allocator<DimDiff<DimMinimum<R, C>, 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<T, R, C>
where
R: DimMin<C>,
DefaultAllocator: Allocator<R, C> + Allocator<DimMinimum<R, C>>,
{
FullPivLU::new(self.into_owned())
}
/// Computes the LU decomposition with partial (row) pivoting of `matrix`.
pub fn lu(self) -> LU<T, R, C>
where
R: DimMin<C>,
DefaultAllocator: Allocator<R, C> + Allocator<DimMinimum<R, C>>,
{
LU::new(self.into_owned())
}
/// Computes the QR decomposition of this matrix.
pub fn qr(self) -> QR<T, R, C>
where
R: DimMin<C>,
DefaultAllocator: Allocator<R, C> + Allocator<R> + Allocator<DimMinimum<R, C>>,
{
QR::new(self.into_owned())
}
/// Computes the QR decomposition (with column pivoting) of this matrix.
pub fn col_piv_qr(self) -> ColPivQR<T, R, C>
where
R: DimMin<C>,
DefaultAllocator: Allocator<R, C>
+ Allocator<R>
+ Allocator<DimMinimum<R, C>>
+ Allocator<DimMinimum<R, C>>,
{
ColPivQR::new(self.into_owned())
}
/// Computes the Singular Value Decomposition using implicit shift.
/// The singular values are guaranteed to be sorted in descending order.
/// If this order is not required consider using `svd_unordered`.
pub fn svd(self, compute_u: bool, compute_v: bool) -> SVD<T, R, C>
where
R: DimMin<C>,
DimMinimum<R, C>: DimSub<U1>, // for Bidiagonal.
DefaultAllocator: Allocator<R, C>
+ Allocator<C>
+ Allocator<R>
+ Allocator<DimDiff<DimMinimum<R, C>, U1>>
+ Allocator<DimMinimum<R, C>, C>
+ Allocator<R, DimMinimum<R, C>>
+ Allocator<DimMinimum<R, C>>,
{
SVD::new(self.into_owned(), compute_u, compute_v)
}
/// Computes the Singular Value Decomposition using implicit shift.
/// The singular values are not guaranteed to be sorted in any particular order.
/// If a descending order is required, consider using `svd` instead.
pub fn svd_unordered(self, compute_u: bool, compute_v: bool) -> SVD<T, R, C>
where
R: DimMin<C>,
DimMinimum<R, C>: DimSub<U1>, // for Bidiagonal.
DefaultAllocator: Allocator<R, C>
+ Allocator<C>
+ Allocator<R>
+ Allocator<DimDiff<DimMinimum<R, C>, U1>>
+ Allocator<DimMinimum<R, C>, C>
+ Allocator<R, DimMinimum<R, C>>
+ Allocator<DimMinimum<R, C>>,
{
SVD::new_unordered(self.into_owned(), compute_u, compute_v)
}
/// Attempts to compute the Singular Value Decomposition of `matrix` using implicit shift.
/// The singular values are guaranteed to be sorted in descending order.
/// If this order is not required consider using `try_svd_unordered`.
///
/// # 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: T::RealField,
max_niter: usize,
) -> Option<SVD<T, R, C>>
where
R: DimMin<C>,
DimMinimum<R, C>: DimSub<U1>, // for Bidiagonal.
DefaultAllocator: Allocator<R, C>
+ Allocator<C>
+ Allocator<R>
+ Allocator<DimDiff<DimMinimum<R, C>, U1>>
+ Allocator<DimMinimum<R, C>, C>
+ Allocator<R, DimMinimum<R, C>>
+ Allocator<DimMinimum<R, C>>,
{
SVD::try_new(self.into_owned(), compute_u, compute_v, eps, max_niter)
}
/// Attempts to compute the Singular Value Decomposition of `matrix` using implicit shift.
/// The singular values are not guaranteed to be sorted in any particular order.
/// If a descending order is required, consider using `try_svd` instead.
///
/// # 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_unordered(
self,
compute_u: bool,
compute_v: bool,
eps: T::RealField,
max_niter: usize,
) -> Option<SVD<T, R, C>>
where
R: DimMin<C>,
DimMinimum<R, C>: DimSub<U1>, // for Bidiagonal.
DefaultAllocator: Allocator<R, C>
+ Allocator<C>
+ Allocator<R>
+ Allocator<DimDiff<DimMinimum<R, C>, U1>>
+ Allocator<DimMinimum<R, C>, C>
+ Allocator<R, DimMinimum<R, C>>
+ Allocator<DimMinimum<R, C>>
+ Allocator<DimMinimum<R, C>>
+ Allocator<DimDiff<DimMinimum<R, C>, U1>>,
{
SVD::try_new_unordered(self.into_owned(), compute_u, compute_v, eps, max_niter)
}
/// Computes the Polar Decomposition of a `matrix` (indirectly uses SVD).
pub fn polar(self) -> (OMatrix<T, R, R>, OMatrix<T, R, C>)
where
R: DimMin<C>,
DimMinimum<R, C>: DimSub<U1>, // for Bidiagonal.
DefaultAllocator: Allocator<R, C>
+ Allocator<DimMinimum<R, C>, R>
+ Allocator<DimMinimum<R, C>>
+ Allocator<R, R>
+ Allocator<DimMinimum<R, C>, DimMinimum<R, C>>
+ Allocator<C>
+ Allocator<R>
+ Allocator<DimDiff<DimMinimum<R, C>, U1>>
+ Allocator<DimMinimum<R, C>, C>
+ Allocator<R, DimMinimum<R, C>>
+ Allocator<DimMinimum<R, C>>
+ Allocator<DimMinimum<R, C>>
+ Allocator<DimDiff<DimMinimum<R, C>, U1>>,
{
SVD::new_unordered(self.into_owned(), true, true)
.to_polar()
.unwrap()
}
/// Attempts to compute the Polar Decomposition of a `matrix` (indirectly uses SVD).
///
/// # Arguments
///
/// * `eps` − tolerance used to determine when a value converged to 0 when computing the SVD.
/// * `max_niter` − maximum total number of iterations performed by the SVD computation algorithm.
pub fn try_polar(
self,
eps: T::RealField,
max_niter: usize,
) -> Option<(OMatrix<T, R, R>, OMatrix<T, R, C>)>
where
R: DimMin<C>,
DimMinimum<R, C>: DimSub<U1>, // for Bidiagonal.
DefaultAllocator: Allocator<R, C>
+ Allocator<DimMinimum<R, C>, R>
+ Allocator<DimMinimum<R, C>>
+ Allocator<R, R>
+ Allocator<DimMinimum<R, C>, DimMinimum<R, C>>
+ Allocator<C>
+ Allocator<R>
+ Allocator<DimDiff<DimMinimum<R, C>, U1>>
+ Allocator<DimMinimum<R, C>, C>
+ Allocator<R, DimMinimum<R, C>>
+ Allocator<DimMinimum<R, C>>
+ Allocator<DimMinimum<R, C>>
+ Allocator<DimDiff<DimMinimum<R, C>, U1>>,
{
SVD::try_new_unordered(self.into_owned(), true, true, eps, max_niter)
.and_then(|svd| svd.to_polar())
}
}
/// # 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<T: ComplexField, D: Dim, S: Storage<T, D, D>> Matrix<T, D, D, S> {
/// 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<T, D>>
where
DefaultAllocator: Allocator<D, D>,
{
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<UDU<T, D>>
where
T: RealField,
DefaultAllocator: Allocator<D> + Allocator<D, D>,
{
UDU::new(self.into_owned())
}
/// Computes the Hessenberg decomposition of this matrix using householder reflections.
pub fn hessenberg(self) -> Hessenberg<T, D>
where
D: DimSub<U1>,
DefaultAllocator: Allocator<D, D> + Allocator<D> + Allocator<DimDiff<D, U1>>,
{
Hessenberg::new(self.into_owned())
}
/// Computes the Schur decomposition of a square matrix.
pub fn schur(self) -> Schur<T, D>
where
D: DimSub<U1>, // For Hessenberg.
DefaultAllocator: Allocator<D, DimDiff<D, U1>>
+ Allocator<DimDiff<D, U1>>
+ Allocator<D, D>
+ Allocator<D>,
{
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: T::RealField, max_niter: usize) -> Option<Schur<T, D>>
where
D: DimSub<U1>, // For Hessenberg.
DefaultAllocator: Allocator<D, DimDiff<D, U1>>
+ Allocator<DimDiff<D, U1>>
+ Allocator<D, D>
+ Allocator<D>,
{
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<T, D>
where
D: DimSub<U1>,
DefaultAllocator:
Allocator<D, D> + Allocator<DimDiff<D, U1>> + Allocator<D> + Allocator<DimDiff<D, U1>>,
{
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: T::RealField,
max_niter: usize,
) -> Option<SymmetricEigen<T, D>>
where
D: DimSub<U1>,
DefaultAllocator:
Allocator<D, D> + Allocator<DimDiff<D, U1>> + Allocator<D> + Allocator<DimDiff<D, U1>>,
{
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<T, D>
where
D: DimSub<U1>,
DefaultAllocator: Allocator<D, D> + Allocator<DimDiff<D, U1>>,
{
SymmetricTridiagonal::new(self.into_owned())
}
}