Type Definition nalgebra::base::SquareMatrix [−][src]
pub type SquareMatrix<N, D, S> = Matrix<N, D, D, S>;
Expand description
A square matrix.
Implementations
impl<N, D1: Dim, S: StorageMut<N, D1, D1>> SquareMatrix<N, D1, S> where
N: Scalar + Zero + One + ClosedAdd + ClosedMul,
impl<N, D1: Dim, S: StorageMut<N, D1, D1>> SquareMatrix<N, D1, S> where
N: Scalar + Zero + One + ClosedAdd + ClosedMul,
pub fn quadform_tr_with_workspace<D2, S2, R3, C3, S3, D4, S4>(
&mut self,
work: &mut Vector<N, D2, S2>,
alpha: N,
lhs: &Matrix<N, R3, C3, S3>,
mid: &SquareMatrix<N, D4, S4>,
beta: N
) where
D2: Dim,
R3: Dim,
C3: Dim,
D4: Dim,
S2: StorageMut<N, D2>,
S3: Storage<N, R3, C3>,
S4: Storage<N, D4, D4>,
ShapeConstraint: DimEq<D1, D2> + DimEq<D1, R3> + DimEq<D2, R3> + DimEq<C3, D4>,
pub fn quadform_tr_with_workspace<D2, S2, R3, C3, S3, D4, S4>(
&mut self,
work: &mut Vector<N, D2, S2>,
alpha: N,
lhs: &Matrix<N, R3, C3, S3>,
mid: &SquareMatrix<N, D4, S4>,
beta: N
) where
D2: Dim,
R3: Dim,
C3: Dim,
D4: Dim,
S2: StorageMut<N, D2>,
S3: Storage<N, R3, C3>,
S4: Storage<N, D4, D4>,
ShapeConstraint: DimEq<D1, D2> + DimEq<D1, R3> + DimEq<D2, R3> + DimEq<C3, D4>,
Computes the quadratic form self = alpha * lhs * mid * lhs.transpose() + beta * self
.
This uses the provided workspace work
to avoid allocations for intermediate results.
Examples:
// Note that all those would also work with statically-sized matrices. // We use DMatrix/DVector since that's the only case where pre-allocating the // workspace is actually useful (assuming the same workspace is re-used for // several computations) because it avoids repeated dynamic allocations. let mut mat = DMatrix::identity(2, 2); let lhs = DMatrix::from_row_slice(2, 3, &[1.0, 2.0, 3.0, 4.0, 5.0, 6.0]); let mid = DMatrix::from_row_slice(3, 3, &[0.1, 0.2, 0.3, 0.5, 0.6, 0.7, 0.9, 1.0, 1.1]); // The random shows that values on the workspace do not // matter as they will be overwritten. let mut workspace = DVector::new_random(2); let expected = &lhs * &mid * lhs.transpose() * 10.0 + &mat * 5.0; mat.quadform_tr_with_workspace(&mut workspace, 10.0, &lhs, &mid, 5.0); assert_relative_eq!(mat, expected);
pub fn quadform_tr<R3, C3, S3, D4, S4>(
&mut self,
alpha: N,
lhs: &Matrix<N, R3, C3, S3>,
mid: &SquareMatrix<N, D4, S4>,
beta: N
) where
R3: Dim,
C3: Dim,
D4: Dim,
S3: Storage<N, R3, C3>,
S4: Storage<N, D4, D4>,
ShapeConstraint: DimEq<D1, D1> + DimEq<D1, R3> + DimEq<C3, D4>,
DefaultAllocator: Allocator<N, D1>,
pub fn quadform_tr<R3, C3, S3, D4, S4>(
&mut self,
alpha: N,
lhs: &Matrix<N, R3, C3, S3>,
mid: &SquareMatrix<N, D4, S4>,
beta: N
) where
R3: Dim,
C3: Dim,
D4: Dim,
S3: Storage<N, R3, C3>,
S4: Storage<N, D4, D4>,
ShapeConstraint: DimEq<D1, D1> + DimEq<D1, R3> + DimEq<C3, D4>,
DefaultAllocator: Allocator<N, D1>,
Computes the quadratic form self = alpha * lhs * mid * lhs.transpose() + beta * self
.
This allocates a workspace vector of dimension D1 for intermediate results.
If D1
is a type-level integer, then the allocation is performed on the stack.
Use .quadform_tr_with_workspace(...)
instead to avoid allocations.
Examples:
let mut mat = Matrix2::identity(); let lhs = Matrix2x3::new(1.0, 2.0, 3.0, 4.0, 5.0, 6.0); let mid = Matrix3::new(0.1, 0.2, 0.3, 0.5, 0.6, 0.7, 0.9, 1.0, 1.1); let expected = lhs * mid * lhs.transpose() * 10.0 + mat * 5.0; mat.quadform_tr(10.0, &lhs, &mid, 5.0); assert_relative_eq!(mat, expected);
pub fn quadform_with_workspace<D2, S2, D3, S3, R4, C4, S4>(
&mut self,
work: &mut Vector<N, D2, S2>,
alpha: N,
mid: &SquareMatrix<N, D3, S3>,
rhs: &Matrix<N, R4, C4, S4>,
beta: N
) where
D2: Dim,
D3: Dim,
R4: Dim,
C4: Dim,
S2: StorageMut<N, D2>,
S3: Storage<N, D3, D3>,
S4: Storage<N, R4, C4>,
ShapeConstraint: DimEq<D3, R4> + DimEq<D1, C4> + DimEq<D2, D3> + AreMultipliable<C4, R4, D2, U1>,
pub fn quadform_with_workspace<D2, S2, D3, S3, R4, C4, S4>(
&mut self,
work: &mut Vector<N, D2, S2>,
alpha: N,
mid: &SquareMatrix<N, D3, S3>,
rhs: &Matrix<N, R4, C4, S4>,
beta: N
) where
D2: Dim,
D3: Dim,
R4: Dim,
C4: Dim,
S2: StorageMut<N, D2>,
S3: Storage<N, D3, D3>,
S4: Storage<N, R4, C4>,
ShapeConstraint: DimEq<D3, R4> + DimEq<D1, C4> + DimEq<D2, D3> + AreMultipliable<C4, R4, D2, U1>,
Computes the quadratic form self = alpha * rhs.transpose() * mid * rhs + beta * self
.
This uses the provided workspace work
to avoid allocations for intermediate results.
// Note that all those would also work with statically-sized matrices. // We use DMatrix/DVector since that's the only case where pre-allocating the // workspace is actually useful (assuming the same workspace is re-used for // several computations) because it avoids repeated dynamic allocations. let mut mat = DMatrix::identity(2, 2); let rhs = DMatrix::from_row_slice(3, 2, &[1.0, 2.0, 3.0, 4.0, 5.0, 6.0]); let mid = DMatrix::from_row_slice(3, 3, &[0.1, 0.2, 0.3, 0.5, 0.6, 0.7, 0.9, 1.0, 1.1]); // The random shows that values on the workspace do not // matter as they will be overwritten. let mut workspace = DVector::new_random(3); let expected = rhs.transpose() * &mid * &rhs * 10.0 + &mat * 5.0; mat.quadform_with_workspace(&mut workspace, 10.0, &mid, &rhs, 5.0); assert_relative_eq!(mat, expected);
pub fn quadform<D2, S2, R3, C3, S3>(
&mut self,
alpha: N,
mid: &SquareMatrix<N, D2, S2>,
rhs: &Matrix<N, R3, C3, S3>,
beta: N
) where
D2: Dim,
R3: Dim,
C3: Dim,
S2: Storage<N, D2, D2>,
S3: Storage<N, R3, C3>,
ShapeConstraint: DimEq<D2, R3> + DimEq<D1, C3> + AreMultipliable<C3, R3, D2, U1>,
DefaultAllocator: Allocator<N, D2>,
pub fn quadform<D2, S2, R3, C3, S3>(
&mut self,
alpha: N,
mid: &SquareMatrix<N, D2, S2>,
rhs: &Matrix<N, R3, C3, S3>,
beta: N
) where
D2: Dim,
R3: Dim,
C3: Dim,
S2: Storage<N, D2, D2>,
S3: Storage<N, R3, C3>,
ShapeConstraint: DimEq<D2, R3> + DimEq<D1, C3> + AreMultipliable<C3, R3, D2, U1>,
DefaultAllocator: Allocator<N, D2>,
Computes the quadratic form self = alpha * rhs.transpose() * mid * rhs + beta * self
.
This allocates a workspace vector of dimension D2 for intermediate results.
If D2
is a type-level integer, then the allocation is performed on the stack.
Use .quadform_with_workspace(...)
instead to avoid allocations.
let mut mat = Matrix2::identity(); let rhs = Matrix3x2::new(1.0, 2.0, 3.0, 4.0, 5.0, 6.0); let mid = Matrix3::new(0.1, 0.2, 0.3, 0.5, 0.6, 0.7, 0.9, 1.0, 1.1); let expected = rhs.transpose() * mid * rhs * 10.0 + mat * 5.0; mat.quadform(10.0, &mid, &rhs, 5.0); assert_relative_eq!(mat, expected);
pub fn append_scaling(&self, scaling: N) -> MatrixN<N, D> where
D: DimNameSub<U1>,
DefaultAllocator: Allocator<N, D, D>,
pub fn append_scaling(&self, scaling: N) -> MatrixN<N, D> where
D: DimNameSub<U1>,
DefaultAllocator: Allocator<N, D, D>,
Computes the transformation equal to self
followed by an uniform scaling factor.
pub fn prepend_scaling(&self, scaling: N) -> MatrixN<N, D> where
D: DimNameSub<U1>,
DefaultAllocator: Allocator<N, D, D>,
pub fn prepend_scaling(&self, scaling: N) -> MatrixN<N, D> where
D: DimNameSub<U1>,
DefaultAllocator: Allocator<N, D, D>,
Computes the transformation equal to an uniform scaling factor followed by self
.
pub fn append_nonuniform_scaling<SB>(
&self,
scaling: &Vector<N, DimNameDiff<D, U1>, SB>
) -> MatrixN<N, D> where
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
DefaultAllocator: Allocator<N, D, D>,
pub fn append_nonuniform_scaling<SB>(
&self,
scaling: &Vector<N, DimNameDiff<D, U1>, SB>
) -> MatrixN<N, D> where
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
DefaultAllocator: Allocator<N, D, D>,
Computes the transformation equal to self
followed by a non-uniform scaling factor.
pub fn prepend_nonuniform_scaling<SB>(
&self,
scaling: &Vector<N, DimNameDiff<D, U1>, SB>
) -> MatrixN<N, D> where
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
DefaultAllocator: Allocator<N, D, D>,
pub fn prepend_nonuniform_scaling<SB>(
&self,
scaling: &Vector<N, DimNameDiff<D, U1>, SB>
) -> MatrixN<N, D> where
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
DefaultAllocator: Allocator<N, D, D>,
Computes the transformation equal to a non-uniform scaling factor followed by self
.
pub fn append_translation<SB>(
&self,
shift: &Vector<N, DimNameDiff<D, U1>, SB>
) -> MatrixN<N, D> where
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
DefaultAllocator: Allocator<N, D, D>,
pub fn append_translation<SB>(
&self,
shift: &Vector<N, DimNameDiff<D, U1>, SB>
) -> MatrixN<N, D> where
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
DefaultAllocator: Allocator<N, D, D>,
Computes the transformation equal to self
followed by a translation.
pub fn prepend_translation<SB>(
&self,
shift: &Vector<N, DimNameDiff<D, U1>, SB>
) -> MatrixN<N, D> where
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
DefaultAllocator: Allocator<N, D, D> + Allocator<N, DimNameDiff<D, U1>>,
pub fn prepend_translation<SB>(
&self,
shift: &Vector<N, DimNameDiff<D, U1>, SB>
) -> MatrixN<N, D> where
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
DefaultAllocator: Allocator<N, D, D> + Allocator<N, DimNameDiff<D, U1>>,
Computes the transformation equal to a translation followed by self
.
pub fn append_scaling_mut(&mut self, scaling: N) where
S: StorageMut<N, D, D>,
D: DimNameSub<U1>,
pub fn append_scaling_mut(&mut self, scaling: N) where
S: StorageMut<N, D, D>,
D: DimNameSub<U1>,
Computes in-place the transformation equal to self
followed by an uniform scaling factor.
pub fn prepend_scaling_mut(&mut self, scaling: N) where
S: StorageMut<N, D, D>,
D: DimNameSub<U1>,
pub fn prepend_scaling_mut(&mut self, scaling: N) where
S: StorageMut<N, D, D>,
D: DimNameSub<U1>,
Computes in-place the transformation equal to an uniform scaling factor followed by self
.
pub fn append_nonuniform_scaling_mut<SB>(
&mut self,
scaling: &Vector<N, DimNameDiff<D, U1>, SB>
) where
S: StorageMut<N, D, D>,
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
pub fn append_nonuniform_scaling_mut<SB>(
&mut self,
scaling: &Vector<N, DimNameDiff<D, U1>, SB>
) where
S: StorageMut<N, D, D>,
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
Computes in-place the transformation equal to self
followed by a non-uniform scaling factor.
pub fn prepend_nonuniform_scaling_mut<SB>(
&mut self,
scaling: &Vector<N, DimNameDiff<D, U1>, SB>
) where
S: StorageMut<N, D, D>,
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
pub fn prepend_nonuniform_scaling_mut<SB>(
&mut self,
scaling: &Vector<N, DimNameDiff<D, U1>, SB>
) where
S: StorageMut<N, D, D>,
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
Computes in-place the transformation equal to a non-uniform scaling factor followed by self
.
pub fn append_translation_mut<SB>(
&mut self,
shift: &Vector<N, DimNameDiff<D, U1>, SB>
) where
S: StorageMut<N, D, D>,
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
pub fn append_translation_mut<SB>(
&mut self,
shift: &Vector<N, DimNameDiff<D, U1>, SB>
) where
S: StorageMut<N, D, D>,
D: DimNameSub<U1>,
SB: Storage<N, DimNameDiff<D, U1>>,
Computes the transformation equal to self
followed by a translation.
pub fn prepend_translation_mut<SB>(
&mut self,
shift: &Vector<N, DimNameDiff<D, U1>, SB>
) where
D: DimNameSub<U1>,
S: StorageMut<N, D, D>,
SB: Storage<N, DimNameDiff<D, U1>>,
DefaultAllocator: Allocator<N, DimNameDiff<D, U1>>,
pub fn prepend_translation_mut<SB>(
&mut self,
shift: &Vector<N, DimNameDiff<D, U1>, SB>
) where
D: DimNameSub<U1>,
S: StorageMut<N, D, D>,
SB: Storage<N, DimNameDiff<D, U1>>,
DefaultAllocator: Allocator<N, DimNameDiff<D, U1>>,
Computes the transformation equal to a translation followed by self
.
impl<N: RealField, D: DimNameSub<U1>, S: Storage<N, D, D>> SquareMatrix<N, D, S> where
DefaultAllocator: Allocator<N, D, D> + Allocator<N, DimNameDiff<D, U1>> + Allocator<N, DimNameDiff<D, U1>, DimNameDiff<D, U1>>,
impl<N: RealField, D: DimNameSub<U1>, S: Storage<N, D, D>> SquareMatrix<N, D, S> where
DefaultAllocator: Allocator<N, D, D> + Allocator<N, DimNameDiff<D, U1>> + Allocator<N, DimNameDiff<D, U1>, DimNameDiff<D, U1>>,
pub fn transform_vector(
&self,
v: &VectorN<N, DimNameDiff<D, U1>>
) -> VectorN<N, DimNameDiff<D, U1>>
pub fn transform_vector(
&self,
v: &VectorN<N, DimNameDiff<D, U1>>
) -> VectorN<N, DimNameDiff<D, U1>>
Transforms the given vector, assuming the matrix self
uses homogeneous coordinates.
pub fn transform_point(
&self,
pt: &Point<N, DimNameDiff<D, U1>>
) -> Point<N, DimNameDiff<D, U1>>
pub fn transform_point(
&self,
pt: &Point<N, DimNameDiff<D, U1>>
) -> Point<N, DimNameDiff<D, U1>>
Transforms the given point, assuming the matrix self
uses homogeneous coordinates.
The diagonal of this matrix.
pub fn map_diagonal<N2: Scalar>(&self, f: impl FnMut(N) -> N2) -> VectorN<N2, D> where
DefaultAllocator: Allocator<N2, D>,
pub fn map_diagonal<N2: Scalar>(&self, f: impl FnMut(N) -> N2) -> VectorN<N2, D> where
DefaultAllocator: Allocator<N2, D>,
Apply the given function to this matrix’s diagonal and returns it.
This is a more efficient version of self.diagonal().map(f)
since this
allocates only once.
The symmetric part of self
, i.e., 0.5 * (self + self.transpose())
.
The hermitian part of self
, i.e., 0.5 * (self + self.adjoint())
.
impl<N: RealField, D: Dim, S: Storage<N, D, D>> SquareMatrix<N, D, S> where
DefaultAllocator: Allocator<N, D, D>,
impl<N: RealField, D: Dim, S: Storage<N, D, D>> SquareMatrix<N, D, S> where
DefaultAllocator: Allocator<N, D, D>,
pub fn is_special_orthogonal(&self, eps: N) -> bool where
D: DimMin<D, Output = D>,
DefaultAllocator: Allocator<(usize, usize), D>,
pub fn is_special_orthogonal(&self, eps: N) -> bool where
D: DimMin<D, Output = D>,
DefaultAllocator: Allocator<(usize, usize), D>,
Checks that this matrix is orthogonal and has a determinant equal to 1.
Returns true
if this matrix is invertible.
pub fn determinant(&self) -> N where
DefaultAllocator: Allocator<N, D, D> + Allocator<(usize, usize), D>,
pub fn determinant(&self) -> N where
DefaultAllocator: Allocator<N, D, D> + Allocator<(usize, usize), D>,
Computes the matrix determinant.
If the matrix has a dimension larger than 3, an LU decomposition is used.
Attempts to invert this matrix.
Attempts to invert this matrix in-place. Returns false
and leaves self
untouched if
inversion fails.
Computes the eigenvalues of this matrix.
pub fn complex_eigenvalues(&self) -> VectorN<NumComplex<N>, D> where
N: RealField,
DefaultAllocator: Allocator<NumComplex<N>, D>,
pub fn complex_eigenvalues(&self) -> VectorN<NumComplex<N>, D> where
N: RealField,
DefaultAllocator: Allocator<NumComplex<N>, D>,
Computes the eigenvalues of this matrix.
pub fn solve_lower_triangular<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn solve_lower_triangular<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Computes the solution of the linear system self . x = b
where x
is the unknown and only
the lower-triangular part of self
(including the diagonal) is considered not-zero.
pub fn solve_upper_triangular<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn solve_upper_triangular<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Computes the solution of the linear system self . x = b
where x
is the unknown and only
the upper-triangular part of self
(including the diagonal) is considered not-zero.
pub fn solve_lower_triangular_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>,
pub fn solve_lower_triangular_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>,
Solves the linear system self . x = b
where x
is the unknown and only the
lower-triangular part of self
(including the diagonal) is considered not-zero.
pub fn solve_lower_triangular_with_diag_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>,
diag: N
) -> bool where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn solve_lower_triangular_with_diag_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>,
diag: N
) -> bool where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Solves the linear system self . x = b
where x
is the unknown and only the
lower-triangular part of self
is considered not-zero. The diagonal is never read as it is
assumed to be equal to diag
. Returns false
and does not modify its inputs if diag
is zero.
pub fn solve_upper_triangular_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>,
pub fn solve_upper_triangular_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>,
Solves the linear system self . x = b
where x
is the unknown and only the
upper-triangular part of self
(including the diagonal) is considered not-zero.
pub fn tr_solve_lower_triangular<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn tr_solve_lower_triangular<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Computes the solution of the linear system self.transpose() . x = b
where x
is the unknown and only
the lower-triangular part of self
(including the diagonal) is considered not-zero.
pub fn tr_solve_upper_triangular<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn tr_solve_upper_triangular<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Computes the solution of the linear system self.transpose() . x = b
where x
is the unknown and only
the upper-triangular part of self
(including the diagonal) is considered not-zero.
pub fn tr_solve_lower_triangular_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>,
pub fn tr_solve_lower_triangular_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>,
Solves the linear system self.transpose() . x = b
where x
is the unknown and only the
lower-triangular part of self
(including the diagonal) is considered not-zero.
pub fn tr_solve_upper_triangular_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>,
pub fn tr_solve_upper_triangular_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>,
Solves the linear system self.transpose() . x = b
where x
is the unknown and only the
upper-triangular part of self
(including the diagonal) is considered not-zero.
pub fn ad_solve_lower_triangular<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn ad_solve_lower_triangular<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Computes the solution of the linear system self.adjoint() . x = b
where x
is the unknown and only
the lower-triangular part of self
(including the diagonal) is considered not-zero.
pub fn ad_solve_upper_triangular<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn ad_solve_upper_triangular<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Computes the solution of the linear system self.adjoint() . x = b
where x
is the unknown and only
the upper-triangular part of self
(including the diagonal) is considered not-zero.
pub fn ad_solve_lower_triangular_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>,
pub fn ad_solve_lower_triangular_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>,
Solves the linear system self.adjoint() . x = b
where x
is the unknown and only the
lower-triangular part of self
(including the diagonal) is considered not-zero.
pub fn ad_solve_upper_triangular_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>,
pub fn ad_solve_upper_triangular_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>,
Solves the linear system self.adjoint() . x = b
where x
is the unknown and only the
upper-triangular part of self
(including the diagonal) is considered not-zero.
pub fn solve_lower_triangular_unchecked<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> MatrixMN<N, R2, C2> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn solve_lower_triangular_unchecked<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> MatrixMN<N, R2, C2> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Computes the solution of the linear system self . x = b
where x
is the unknown and only
the lower-triangular part of self
(including the diagonal) is considered not-zero.
pub fn solve_upper_triangular_unchecked<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> MatrixMN<N, R2, C2> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn solve_upper_triangular_unchecked<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> MatrixMN<N, R2, C2> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Computes the solution of the linear system self . x = b
where x
is the unknown and only
the upper-triangular part of self
(including the diagonal) is considered not-zero.
pub fn solve_lower_triangular_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn solve_lower_triangular_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Solves the linear system self . x = b
where x
is the unknown and only the
lower-triangular part of self
(including the diagonal) is considered not-zero.
pub fn solve_lower_triangular_with_diag_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>,
diag: N
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn solve_lower_triangular_with_diag_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>,
diag: N
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Solves the linear system self . x = b
where x
is the unknown and only the
lower-triangular part of self
is considered not-zero. The diagonal is never read as it is
assumed to be equal to diag
. Returns false
and does not modify its inputs if diag
is zero.
pub fn solve_upper_triangular_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn solve_upper_triangular_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Solves the linear system self . x = b
where x
is the unknown and only the
upper-triangular part of self
(including the diagonal) is considered not-zero.
pub fn tr_solve_lower_triangular_unchecked<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> MatrixMN<N, R2, C2> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn tr_solve_lower_triangular_unchecked<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> MatrixMN<N, R2, C2> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Computes the solution of the linear system self.transpose() . x = b
where x
is the unknown and only
the lower-triangular part of self
(including the diagonal) is considered not-zero.
pub fn tr_solve_upper_triangular_unchecked<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> MatrixMN<N, R2, C2> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn tr_solve_upper_triangular_unchecked<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> MatrixMN<N, R2, C2> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Computes the solution of the linear system self.transpose() . x = b
where x
is the unknown and only
the upper-triangular part of self
(including the diagonal) is considered not-zero.
pub fn tr_solve_lower_triangular_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn tr_solve_lower_triangular_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Solves the linear system self.transpose() . x = b
where x
is the unknown and only the
lower-triangular part of self
(including the diagonal) is considered not-zero.
pub fn tr_solve_upper_triangular_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn tr_solve_upper_triangular_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Solves the linear system self.transpose() . x = b
where x
is the unknown and only the
upper-triangular part of self
(including the diagonal) is considered not-zero.
pub fn ad_solve_lower_triangular_unchecked<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> MatrixMN<N, R2, C2> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn ad_solve_lower_triangular_unchecked<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> MatrixMN<N, R2, C2> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Computes the solution of the linear system self.adjoint() . x = b
where x
is the unknown and only
the lower-triangular part of self
(including the diagonal) is considered not-zero.
pub fn ad_solve_upper_triangular_unchecked<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> MatrixMN<N, R2, C2> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn ad_solve_upper_triangular_unchecked<R2: Dim, C2: Dim, S2>(
&self,
b: &Matrix<N, R2, C2, S2>
) -> MatrixMN<N, R2, C2> where
S2: Storage<N, R2, C2>,
DefaultAllocator: Allocator<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Computes the solution of the linear system self.adjoint() . x = b
where x
is the unknown and only
the upper-triangular part of self
(including the diagonal) is considered not-zero.
pub fn ad_solve_lower_triangular_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn ad_solve_lower_triangular_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Solves the linear system self.adjoint() . x = b
where x
is the unknown and only the
lower-triangular part of self
(including the diagonal) is considered not-zero.
pub fn ad_solve_upper_triangular_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
pub fn ad_solve_upper_triangular_unchecked_mut<R2: Dim, C2: Dim, S2>(
&self,
b: &mut Matrix<N, R2, C2, S2>
) where
S2: StorageMut<N, R2, C2>,
ShapeConstraint: SameNumberOfRows<R2, D>,
Solves the linear system self.adjoint() . x = b
where x
is the unknown and only the
upper-triangular part of self
(including the diagonal) is considered not-zero.
Computes the eigenvalues of this symmetric matrix.
Only the lower-triangular part of the matrix is read.