pub struct Cholesky<T: SimdComplexField, D: Dim> where
DefaultAllocator: Allocator<T, D, D>, { /* private fields */ }
Expand description
The Cholesky decomposition of a symmetric-definite-positive matrix.
Implementations
sourceimpl<T: SimdComplexField, D: Dim> Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
impl<T: SimdComplexField, D: Dim> Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
sourcepub fn new_unchecked(matrix: OMatrix<T, D, D>) -> Self
pub fn new_unchecked(matrix: OMatrix<T, D, D>) -> Self
Computes the Cholesky decomposition of matrix
without checking that the matrix is definite-positive.
If the input matrix is not definite-positive, the decomposition may contain trash values (Inf, NaN, etc.)
sourcepub fn unpack(self) -> OMatrix<T, D, D>
pub fn unpack(self) -> OMatrix<T, D, D>
Retrieves the lower-triangular factor of the Cholesky decomposition with its strictly upper-triangular part filled with zeros.
sourcepub fn unpack_dirty(self) -> OMatrix<T, D, D>
pub fn unpack_dirty(self) -> OMatrix<T, D, D>
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.
sourcepub fn l(&self) -> OMatrix<T, D, D>
pub fn l(&self) -> OMatrix<T, D, D>
Retrieves the lower-triangular factor of the Cholesky decomposition with its strictly uppen-triangular part filled with zeros.
sourcepub fn l_dirty(&self) -> &OMatrix<T, D, D>
pub fn l_dirty(&self) -> &OMatrix<T, D, D>
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.
sourcepub fn solve_mut<R2: Dim, C2: Dim, S2>(&self, b: &mut Matrix<T, R2, C2, S2>) where
S2: StorageMut<T, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn solve_mut<R2: Dim, C2: Dim, S2>(&self, b: &mut Matrix<T, R2, C2, S2>) where
S2: StorageMut<T, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Solves the system self * x = b
where self
is the decomposed matrix and x
the unknown.
The result is stored on b
.
sourcepub fn solve<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<T, R2, C2, S2>
) -> OMatrix<T, R2, C2> where
S2: Storage<T, R2, C2>,
DefaultAllocator: Allocator<T, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn solve<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<T, R2, C2, S2>
) -> OMatrix<T, R2, C2> where
S2: Storage<T, R2, C2>,
DefaultAllocator: Allocator<T, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Returns the solution of the system self * x = b
where self
is the decomposed matrix and
x
the unknown.
sourcepub fn determinant(&self) -> T::SimdRealField
pub fn determinant(&self) -> T::SimdRealField
Computes the determinant of the decomposed matrix.
sourceimpl<T: ComplexField, D: Dim> Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
impl<T: ComplexField, D: Dim> Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
sourcepub fn new(matrix: OMatrix<T, D, D>) -> Option<Self>
pub fn new(matrix: OMatrix<T, D, D>) -> Option<Self>
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.
sourcepub fn rank_one_update<R2: Dim, S2>(
&mut self,
x: &Vector<T, R2, S2>,
sigma: T::RealField
) where
S2: Storage<T, R2, U1>,
DefaultAllocator: Allocator<T, R2, U1>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn rank_one_update<R2: Dim, S2>(
&mut self,
x: &Vector<T, R2, S2>,
sigma: T::RealField
) where
S2: Storage<T, R2, U1>,
DefaultAllocator: Allocator<T, R2, U1>,
ShapeConstraint: SameNumberOfRows<R2, D>,
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())
.
sourcepub fn insert_column<R2, S2>(
&self,
j: usize,
col: Vector<T, R2, S2>
) -> Cholesky<T, DimSum<D, U1>> where
D: DimAdd<U1>,
R2: Dim,
S2: Storage<T, R2, U1>,
DefaultAllocator: Allocator<T, DimSum<D, U1>, DimSum<D, U1>> + Allocator<T, R2>,
ShapeConstraint: SameNumberOfRows<R2, DimSum<D, U1>>,
pub fn insert_column<R2, S2>(
&self,
j: usize,
col: Vector<T, R2, S2>
) -> Cholesky<T, DimSum<D, U1>> where
D: DimAdd<U1>,
R2: Dim,
S2: Storage<T, R2, U1>,
DefaultAllocator: Allocator<T, DimSum<D, U1>, DimSum<D, U1>> + Allocator<T, R2>,
ShapeConstraint: SameNumberOfRows<R2, DimSum<D, U1>>,
Updates the decomposition such that we get the decomposition of a matrix with the given column col
in the j
th position.
Since the matrix is square, an identical row will be added in the j
th row.
sourcepub fn remove_column(&self, j: usize) -> Cholesky<T, DimDiff<D, U1>> where
D: DimSub<U1>,
DefaultAllocator: Allocator<T, DimDiff<D, U1>, DimDiff<D, U1>> + Allocator<T, D>,
pub fn remove_column(&self, j: usize) -> Cholesky<T, DimDiff<D, U1>> where
D: DimSub<U1>,
DefaultAllocator: Allocator<T, DimDiff<D, U1>, DimDiff<D, U1>> + Allocator<T, D>,
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.
Trait Implementations
sourceimpl<T: Clone + SimdComplexField, D: Clone + Dim> Clone for Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
impl<T: Clone + SimdComplexField, D: Clone + Dim> Clone for Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
sourceimpl<T: Debug + SimdComplexField, D: Debug + Dim> Debug for Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
impl<T: Debug + SimdComplexField, D: Debug + Dim> Debug for Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
impl<T: SimdComplexField, D: Dim> Copy for Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
OMatrix<T, D, D>: Copy,
Auto Trait Implementations
impl<T, D> !RefUnwindSafe for Cholesky<T, D>
impl<T, D> !Send for Cholesky<T, D>
impl<T, D> !Sync for Cholesky<T, D>
impl<T, D> !Unpin for Cholesky<T, D>
impl<T, D> !UnwindSafe for Cholesky<T, D>
Blanket Implementations
sourceimpl<T> BorrowMut<T> for T where
T: ?Sized,
impl<T> BorrowMut<T> for T where
T: ?Sized,
const: unstable · sourcepub fn borrow_mut(&mut self) -> &mut T
pub fn borrow_mut(&mut self) -> &mut T
Mutably borrows from an owned value. Read more
sourceimpl<SS, SP> SupersetOf<SS> for SP where
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SP where
SS: SubsetOf<SP>,
sourcepub fn to_subset(&self) -> Option<SS>
pub fn to_subset(&self) -> Option<SS>
The inverse inclusion map: attempts to construct self
from the equivalent element of its
superset. Read more
sourcepub fn is_in_subset(&self) -> bool
pub fn is_in_subset(&self) -> bool
Checks if self
is actually part of its subset T
(and can be converted to it).
sourcepub fn to_subset_unchecked(&self) -> SS
pub fn to_subset_unchecked(&self) -> SS
Use with care! Same as self.to_subset
but without any property checks. Always succeeds.
sourcepub fn from_subset(element: &SS) -> SP
pub fn from_subset(element: &SS) -> SP
The inclusion map: converts self
to the equivalent element of its superset.
sourceimpl<T> ToOwned for T where
T: Clone,
impl<T> ToOwned for T where
T: Clone,
type Owned = T
type Owned = T
The resulting type after obtaining ownership.
sourcepub fn to_owned(&self) -> T
pub fn to_owned(&self) -> T
Creates owned data from borrowed data, usually by cloning. Read more
sourcepub fn clone_into(&self, target: &mut T)
pub fn clone_into(&self, target: &mut T)
toowned_clone_into
)Uses borrowed data to replace owned data, usually by cloning. Read more