use crate::storage::Storage;
use crate::{
Allocator, Bidiagonal, Cholesky, ColPivQR, ComplexField, DefaultAllocator, Dim, DimDiff,
DimMin, DimMinimum, DimSub, FullPivLU, Hessenberg, LU, Matrix, OMatrix, QR, RealField, SVD,
Schur, SymmetricEigen, SymmetricTridiagonal, U1, UDU,
};
impl<T: ComplexField, R: Dim, C: Dim, S: Storage<T, R, C>> Matrix<T, R, C, S> {
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())
}
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())
}
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())
}
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())
}
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())
}
pub fn svd(self, compute_u: bool, compute_v: bool) -> SVD<T, R, C>
where
R: DimMin<C>,
DimMinimum<R, C>: DimSub<U1>, 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)
}
pub fn svd_unordered(self, compute_u: bool, compute_v: bool) -> SVD<T, R, C>
where
R: DimMin<C>,
DimMinimum<R, C>: DimSub<U1>, 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)
}
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>, 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)
}
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>, 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)
}
pub fn polar(self) -> (OMatrix<T, R, R>, OMatrix<T, R, C>)
where
R: DimMin<C>,
DimMinimum<R, C>: DimSub<U1>, 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()
}
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>, 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())
}
}
impl<T: ComplexField, D: Dim, S: Storage<T, D, D>> Matrix<T, D, D, S> {
pub fn cholesky(self) -> Option<Cholesky<T, D>>
where
DefaultAllocator: Allocator<D, D>,
{
Cholesky::new(self.into_owned())
}
pub fn udu(self) -> Option<UDU<T, D>>
where
T: RealField,
DefaultAllocator: Allocator<D> + Allocator<D, D>,
{
UDU::new(self.into_owned())
}
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())
}
pub fn schur(self) -> Schur<T, D>
where
D: DimSub<U1>, DefaultAllocator: Allocator<D, DimDiff<D, U1>>
+ Allocator<DimDiff<D, U1>>
+ Allocator<D, D>
+ Allocator<D>,
{
Schur::new(self.into_owned())
}
pub fn try_schur(self, eps: T::RealField, max_niter: usize) -> Option<Schur<T, D>>
where
D: DimSub<U1>, DefaultAllocator: Allocator<D, DimDiff<D, U1>>
+ Allocator<DimDiff<D, U1>>
+ Allocator<D, D>
+ Allocator<D>,
{
Schur::try_new(self.into_owned(), eps, max_niter)
}
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())
}
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)
}
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())
}
}