Struct nalgebra::linalg::Cholesky [−][src]
pub struct Cholesky<T: SimdComplexField, D: Dim> where
DefaultAllocator: Allocator<T, D, D>, { /* fields omitted */ }
Expand description
The Cholesky decomposition of a symmetric-definite-positive matrix.
Implementations
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.)
Retrieves the lower-triangular factor of the Cholesky decomposition with its strictly upper-triangular part filled with zeros.
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.
Retrieves the lower-triangular factor of the Cholesky decomposition with its strictly uppen-triangular part filled with zeros.
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 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
.
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>,
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.
Computes the determinant of the decomposed matrix.
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 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())
.
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.
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
impl<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>,
impl<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<'de, T: SimdComplexField, D: Dim> Deserialize<'de> for Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
DefaultAllocator: Allocator<T, D>,
OMatrix<T, D, D>: Deserialize<'de>,
impl<'de, T: SimdComplexField, D: Dim> Deserialize<'de> for Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
DefaultAllocator: Allocator<T, D>,
OMatrix<T, D, D>: Deserialize<'de>,
fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error> where
__D: Deserializer<'de>,
fn deserialize<__D>(__deserializer: __D) -> Result<Self, __D::Error> where
__D: Deserializer<'de>,
Deserialize this value from the given Serde deserializer. Read more
impl<T: SimdComplexField, D: Dim> Serialize for Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
DefaultAllocator: Allocator<T, D>,
OMatrix<T, D, D>: Serialize,
impl<T: SimdComplexField, D: Dim> Serialize for Cholesky<T, D> where
DefaultAllocator: Allocator<T, D, D>,
DefaultAllocator: Allocator<T, D>,
OMatrix<T, D, D>: Serialize,
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> !UnwindSafe for Cholesky<T, D>
Blanket Implementations
Mutably borrows from an owned value. Read more
The inverse inclusion map: attempts to construct self
from the equivalent element of its
superset. Read more
Checks if self
is actually part of its subset T
(and can be converted to it).
Use with care! Same as self.to_subset
but without any property checks. Always succeeds.
The inclusion map: converts self
to the equivalent element of its superset.