Struct nalgebra::linalg::LU
[−]
[src]
pub struct LU<N: Real, R: DimMin<C>, C: Dim> where
DefaultAllocator: Allocator<N, R, C> + Allocator<(usize, usize), DimMinimum<R, C>>, { /* fields omitted */ }
LU decomposition with partial (row) pivoting.
Methods
impl<N: Real, R: DimMin<C>, C: Dim> LU<N, R, C> where
DefaultAllocator: Allocator<N, R, C> + Allocator<(usize, usize), DimMinimum<R, C>>,
[src]
DefaultAllocator: Allocator<N, R, C> + Allocator<(usize, usize), DimMinimum<R, C>>,
fn new(matrix: MatrixMN<N, R, C>) -> Self
Computes the LU decomposition with partial (row) pivoting of matrix
.
fn l(&self) -> MatrixMN<N, R, DimMinimum<R, C>> where
DefaultAllocator: Allocator<N, R, DimMinimum<R, C>>,
DefaultAllocator: Allocator<N, R, DimMinimum<R, C>>,
The lower triangular matrix of this decomposition.
fn l_unpack(self) -> MatrixMN<N, R, DimMinimum<R, C>> where
DefaultAllocator: Reallocator<N, R, C, R, DimMinimum<R, C>>,
DefaultAllocator: Reallocator<N, R, C, R, DimMinimum<R, C>>,
The lower triangular matrix of this decomposition.
fn u(&self) -> MatrixMN<N, DimMinimum<R, C>, C> where
DefaultAllocator: Allocator<N, DimMinimum<R, C>, C>,
DefaultAllocator: Allocator<N, DimMinimum<R, C>, C>,
The upper triangular matrix of this decomposition.
fn p(&self) -> &PermutationSequence<DimMinimum<R, C>>
The row permutations of this decomposition.
fn unpack(
self
) -> (PermutationSequence<DimMinimum<R, C>>, MatrixMN<N, R, DimMinimum<R, C>>, MatrixMN<N, DimMinimum<R, C>, C>) where
DefaultAllocator: Allocator<N, R, DimMinimum<R, C>> + Allocator<N, DimMinimum<R, C>, C> + Reallocator<N, R, C, R, DimMinimum<R, C>>,
self
) -> (PermutationSequence<DimMinimum<R, C>>, MatrixMN<N, R, DimMinimum<R, C>>, MatrixMN<N, DimMinimum<R, C>, C>) where
DefaultAllocator: Allocator<N, R, DimMinimum<R, C>> + Allocator<N, DimMinimum<R, C>, C> + Reallocator<N, R, C, R, DimMinimum<R, C>>,
The row permutations and two triangular matrices of this decomposition: (P, L, U)
.
impl<N: Real, D: DimMin<D, Output = D>> LU<N, D, D> where
DefaultAllocator: Allocator<N, D, D> + Allocator<(usize, usize), D>,
[src]
DefaultAllocator: Allocator<N, D, D> + Allocator<(usize, usize), D>,
fn solve<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
DefaultAllocator: Allocator<N, R2, C2>,
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
DefaultAllocator: Allocator<N, R2, C2>,
Solves the linear system self * x = b
, where x
is the unknown to be determined.
Returns None
if self
is not invertible.
fn solve_mut<R2: Dim, C2: Dim, S2>(&self, b: &mut Matrix<N, R2, C2, S2>) -> bool where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Solves the linear system self * x = b
, where x
is the unknown to be determined.
If the decomposed matrix is not invertible, this returns false
and its input b
may
be overwritten with garbage.
fn try_inverse(&self) -> Option<MatrixN<N, D>>
Computes the inverse of the decomposed matrix.
Returns None
if the matrix is not invertible.
fn try_inverse_to<S2: StorageMut<N, D, D>>(
&self,
out: &mut Matrix<N, D, D, S2>
) -> bool
&self,
out: &mut Matrix<N, D, D, S2>
) -> bool
Computes the inverse of the decomposed matrix and outputs the result to out
.
If the decomposed matrix is not invertible, this returns false
and out
may be
overwritten with garbage.
fn determinant(&self) -> N
Computes the determinant of the decomposed matrix.
fn is_invertible(&self) -> bool
Indicates if the decomposed matrix is invertible.
Trait Implementations
impl<N: Clone + Real, R: Clone + DimMin<C>, C: Clone + Dim> Clone for LU<N, R, C> where
DefaultAllocator: Allocator<N, R, C> + Allocator<(usize, usize), DimMinimum<R, C>>,
[src]
DefaultAllocator: Allocator<N, R, C> + Allocator<(usize, usize), DimMinimum<R, C>>,
fn clone(&self) -> LU<N, R, C>
Returns a copy of the value. Read more
fn clone_from(&mut self, source: &Self)
1.0.0
Performs copy-assignment from source
. Read more
impl<N: Debug + Real, R: Debug + DimMin<C>, C: Debug + Dim> Debug for LU<N, R, C> where
DefaultAllocator: Allocator<N, R, C> + Allocator<(usize, usize), DimMinimum<R, C>>,
[src]
DefaultAllocator: Allocator<N, R, C> + Allocator<(usize, usize), DimMinimum<R, C>>,
impl<N: Real, R: DimMin<C>, C: Dim> Copy for LU<N, R, C> where
DefaultAllocator: Allocator<N, R, C> + Allocator<(usize, usize), DimMinimum<R, C>>,
MatrixMN<N, R, C>: Copy,
PermutationSequence<DimMinimum<R, C>>: Copy,
[src]
DefaultAllocator: Allocator<N, R, C> + Allocator<(usize, usize), DimMinimum<R, C>>,
MatrixMN<N, R, C>: Copy,
PermutationSequence<DimMinimum<R, C>>: Copy,