Type Definition nalgebra::base::SquareMatrix

source ·
pub type SquareMatrix<T, D, S> = Matrix<T, D, D, S>;
Expand description

A square matrix.

Implementations§

source§

impl<T, D1: Dim, S: StorageMut<T, D1, D1>> SquareMatrix<T, D1, S>where T: Scalar + Zero + One + ClosedAdd + ClosedMul,

source

pub fn quadform_tr_with_workspace<D2, S2, R3, C3, S3, D4, S4>( &mut self, work: &mut Vector<T, D2, S2>, alpha: T, lhs: &Matrix<T, R3, C3, S3>, mid: &SquareMatrix<T, D4, S4>, beta: T )where D2: Dim, R3: Dim, C3: Dim, D4: Dim, S2: StorageMut<T, D2>, S3: Storage<T, R3, C3>, S4: Storage<T, 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);
source

pub fn quadform_tr<R3, C3, S3, D4, S4>( &mut self, alpha: T, lhs: &Matrix<T, R3, C3, S3>, mid: &SquareMatrix<T, D4, S4>, beta: T )where R3: Dim, C3: Dim, D4: Dim, S3: Storage<T, R3, C3>, S4: Storage<T, D4, D4>, ShapeConstraint: DimEq<D1, D1> + DimEq<D1, R3> + DimEq<C3, D4>, DefaultAllocator: Allocator<T, 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);
source

pub fn quadform_with_workspace<D2, S2, D3, S3, R4, C4, S4>( &mut self, work: &mut Vector<T, D2, S2>, alpha: T, mid: &SquareMatrix<T, D3, S3>, rhs: &Matrix<T, R4, C4, S4>, beta: T )where D2: Dim, D3: Dim, R4: Dim, C4: Dim, S2: StorageMut<T, D2>, S3: Storage<T, D3, D3>, S4: Storage<T, 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);
source

pub fn quadform<D2, S2, R3, C3, S3>( &mut self, alpha: T, mid: &SquareMatrix<T, D2, S2>, rhs: &Matrix<T, R3, C3, S3>, beta: T )where D2: Dim, R3: Dim, C3: Dim, S2: Storage<T, D2, D2>, S3: Storage<T, R3, C3>, ShapeConstraint: DimEq<D2, R3> + DimEq<D1, C3> + AreMultipliable<C3, R3, D2, U1>, DefaultAllocator: Allocator<T, 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);
source§

impl<T: Scalar + Zero + One + ClosedMul + ClosedAdd, D: DimName, S: Storage<T, D, D>> SquareMatrix<T, D, S>

source

pub fn append_scaling(&self, scaling: T) -> OMatrix<T, D, D>where D: DimNameSub<U1>, DefaultAllocator: Allocator<T, D, D>,

Computes the transformation equal to self followed by an uniform scaling factor.

source

pub fn prepend_scaling(&self, scaling: T) -> OMatrix<T, D, D>where D: DimNameSub<U1>, DefaultAllocator: Allocator<T, D, D>,

Computes the transformation equal to an uniform scaling factor followed by self.

source

pub fn append_nonuniform_scaling<SB>( &self, scaling: &Vector<T, DimNameDiff<D, U1>, SB> ) -> OMatrix<T, D, D>where D: DimNameSub<U1>, SB: Storage<T, DimNameDiff<D, U1>>, DefaultAllocator: Allocator<T, D, D>,

Computes the transformation equal to self followed by a non-uniform scaling factor.

source

pub fn prepend_nonuniform_scaling<SB>( &self, scaling: &Vector<T, DimNameDiff<D, U1>, SB> ) -> OMatrix<T, D, D>where D: DimNameSub<U1>, SB: Storage<T, DimNameDiff<D, U1>>, DefaultAllocator: Allocator<T, D, D>,

Computes the transformation equal to a non-uniform scaling factor followed by self.

source

pub fn append_translation<SB>( &self, shift: &Vector<T, DimNameDiff<D, U1>, SB> ) -> OMatrix<T, D, D>where D: DimNameSub<U1>, SB: Storage<T, DimNameDiff<D, U1>>, DefaultAllocator: Allocator<T, D, D>,

Computes the transformation equal to self followed by a translation.

source

pub fn prepend_translation<SB>( &self, shift: &Vector<T, DimNameDiff<D, U1>, SB> ) -> OMatrix<T, D, D>where D: DimNameSub<U1>, SB: Storage<T, DimNameDiff<D, U1>>, DefaultAllocator: Allocator<T, D, D> + Allocator<T, DimNameDiff<D, U1>>,

Computes the transformation equal to a translation followed by self.

source

pub fn append_scaling_mut(&mut self, scaling: T)where S: StorageMut<T, D, D>, D: DimNameSub<U1>,

Computes in-place the transformation equal to self followed by an uniform scaling factor.

source

pub fn prepend_scaling_mut(&mut self, scaling: T)where S: StorageMut<T, D, D>, D: DimNameSub<U1>,

Computes in-place the transformation equal to an uniform scaling factor followed by self.

source

pub fn append_nonuniform_scaling_mut<SB>( &mut self, scaling: &Vector<T, DimNameDiff<D, U1>, SB> )where S: StorageMut<T, D, D>, D: DimNameSub<U1>, SB: Storage<T, DimNameDiff<D, U1>>,

Computes in-place the transformation equal to self followed by a non-uniform scaling factor.

source

pub fn prepend_nonuniform_scaling_mut<SB>( &mut self, scaling: &Vector<T, DimNameDiff<D, U1>, SB> )where S: StorageMut<T, D, D>, D: DimNameSub<U1>, SB: Storage<T, DimNameDiff<D, U1>>,

Computes in-place the transformation equal to a non-uniform scaling factor followed by self.

source

pub fn append_translation_mut<SB>( &mut self, shift: &Vector<T, DimNameDiff<D, U1>, SB> )where S: StorageMut<T, D, D>, D: DimNameSub<U1>, SB: Storage<T, DimNameDiff<D, U1>>,

Computes the transformation equal to self followed by a translation.

source

pub fn prepend_translation_mut<SB>( &mut self, shift: &Vector<T, DimNameDiff<D, U1>, SB> )where D: DimNameSub<U1>, S: StorageMut<T, D, D>, SB: Storage<T, DimNameDiff<D, U1>>, DefaultAllocator: Allocator<T, DimNameDiff<D, U1>>,

Computes the transformation equal to a translation followed by self.

source§

impl<T: RealField, D: DimNameSub<U1>, S: Storage<T, D, D>> SquareMatrix<T, D, S>where DefaultAllocator: Allocator<T, D, D> + Allocator<T, DimNameDiff<D, U1>> + Allocator<T, DimNameDiff<D, U1>, DimNameDiff<D, U1>>,

source

pub fn transform_vector( &self, v: &OVector<T, DimNameDiff<D, U1>> ) -> OVector<T, DimNameDiff<D, U1>>

Transforms the given vector, assuming the matrix self uses homogeneous coordinates.

source§

impl<T: RealField, S: Storage<T, Const<3>, Const<3>>> SquareMatrix<T, Const<3>, S>

source

pub fn transform_point(&self, pt: &Point<T, 2>) -> Point<T, 2>

Transforms the given point, assuming the matrix self uses homogeneous coordinates.

source§

impl<T: RealField, S: Storage<T, Const<4>, Const<4>>> SquareMatrix<T, Const<4>, S>

source

pub fn transform_point(&self, pt: &Point<T, 3>) -> Point<T, 3>

Transforms the given point, assuming the matrix self uses homogeneous coordinates.

source§

impl<T: Scalar, D: Dim, S: Storage<T, D, D>> SquareMatrix<T, D, S>

source

pub fn diagonal(&self) -> OVector<T, D>where DefaultAllocator: Allocator<T, D>,

The diagonal of this matrix.

source

pub fn map_diagonal<T2: Scalar>(&self, f: impl FnMut(T) -> T2) -> OVector<T2, D>where DefaultAllocator: Allocator<T2, 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.

source

pub fn trace(&self) -> Twhere T: Scalar + Zero + ClosedAdd,

Computes a trace of a square matrix, i.e., the sum of its diagonal elements.

source§

impl<T: SimdComplexField, D: Dim, S: Storage<T, D, D>> SquareMatrix<T, D, S>

source

pub fn symmetric_part(&self) -> OMatrix<T, D, D>where DefaultAllocator: Allocator<T, D, D>,

The symmetric part of self, i.e., 0.5 * (self + self.transpose()).

source

pub fn hermitian_part(&self) -> OMatrix<T, D, D>where DefaultAllocator: Allocator<T, D, D>,

The hermitian part of self, i.e., 0.5 * (self + self.adjoint()).

source§

impl<T: RealField, D: Dim, S: Storage<T, D, D>> SquareMatrix<T, D, S>where DefaultAllocator: Allocator<T, D, D>,

source

pub fn is_special_orthogonal(&self, eps: T) -> boolwhere D: DimMin<D, Output = D>, DefaultAllocator: Allocator<(usize, usize), D>,

Checks that this matrix is orthogonal and has a determinant equal to 1.

source

pub fn is_invertible(&self) -> bool

Returns true if this matrix is invertible.

source§

impl<T: ComplexField, D: DimMin<D, Output = D>, S: Storage<T, D, D>> SquareMatrix<T, D, S>

source

pub fn determinant(&self) -> Twhere DefaultAllocator: Allocator<T, D, D> + Allocator<(usize, usize), D>,

Computes the matrix determinant.

If the matrix has a dimension larger than 3, an LU decomposition is used.

source§

impl<T: ComplexField, D: Dim, S: Storage<T, D, D>> SquareMatrix<T, D, S>

source

pub fn try_inverse(self) -> Option<OMatrix<T, D, D>>where DefaultAllocator: Allocator<T, D, D>,

Attempts to invert this matrix.

source§

impl<T: ComplexField, D: Dim, S: StorageMut<T, D, D>> SquareMatrix<T, D, S>

source

pub fn try_inverse_mut(&mut self) -> boolwhere DefaultAllocator: Allocator<T, D, D>,

Attempts to invert this matrix in-place. Returns false and leaves self untouched if inversion fails.

source§

impl<T: ComplexField, D, S: Storage<T, D, D>> SquareMatrix<T, D, S>where D: DimSub<U1> + Dim, DefaultAllocator: Allocator<T, D, DimDiff<D, U1>> + Allocator<T, DimDiff<D, U1>> + Allocator<T, D, D> + Allocator<T, D>,

source

pub fn eigenvalues(&self) -> Option<OVector<T, D>>

Computes the eigenvalues of this matrix.

source

pub fn complex_eigenvalues(&self) -> OVector<NumComplex<T>, D>where T: RealField, DefaultAllocator: Allocator<NumComplex<T>, D>,

Computes the eigenvalues of this matrix.

source§

impl<T: ComplexField, D: Dim, S: Storage<T, D, D>> SquareMatrix<T, D, S>

source

pub fn solve_lower_triangular<R2: Dim, C2: Dim, S2>( &self, b: &Matrix<T, R2, C2, S2> ) -> Option<OMatrix<T, R2, C2>>where S2: Storage<T, R2, C2>, DefaultAllocator: Allocator<T, 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.

source

pub fn solve_upper_triangular<R2: Dim, C2: Dim, S2>( &self, b: &Matrix<T, R2, C2, S2> ) -> Option<OMatrix<T, R2, C2>>where S2: Storage<T, R2, C2>, DefaultAllocator: Allocator<T, 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.

source

pub fn solve_lower_triangular_mut<R2: Dim, C2: Dim, S2>( &self, b: &mut Matrix<T, R2, C2, S2> ) -> boolwhere S2: StorageMut<T, 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.

source

pub fn solve_lower_triangular_with_diag_mut<R2: Dim, C2: Dim, S2>( &self, b: &mut Matrix<T, R2, C2, S2>, diag: T ) -> boolwhere S2: StorageMut<T, 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.

source

pub fn solve_upper_triangular_mut<R2: Dim, C2: Dim, S2>( &self, b: &mut Matrix<T, R2, C2, S2> ) -> boolwhere S2: StorageMut<T, 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.

source

pub fn tr_solve_lower_triangular<R2: Dim, C2: Dim, S2>( &self, b: &Matrix<T, R2, C2, S2> ) -> Option<OMatrix<T, R2, C2>>where S2: Storage<T, R2, C2>, DefaultAllocator: Allocator<T, 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.

source

pub fn tr_solve_upper_triangular<R2: Dim, C2: Dim, S2>( &self, b: &Matrix<T, R2, C2, S2> ) -> Option<OMatrix<T, R2, C2>>where S2: Storage<T, R2, C2>, DefaultAllocator: Allocator<T, 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.

source

pub fn tr_solve_lower_triangular_mut<R2: Dim, C2: Dim, S2>( &self, b: &mut Matrix<T, R2, C2, S2> ) -> boolwhere S2: StorageMut<T, 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.

source

pub fn tr_solve_upper_triangular_mut<R2: Dim, C2: Dim, S2>( &self, b: &mut Matrix<T, R2, C2, S2> ) -> boolwhere S2: StorageMut<T, 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.

source

pub fn ad_solve_lower_triangular<R2: Dim, C2: Dim, S2>( &self, b: &Matrix<T, R2, C2, S2> ) -> Option<OMatrix<T, R2, C2>>where S2: Storage<T, R2, C2>, DefaultAllocator: Allocator<T, 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.

source

pub fn ad_solve_upper_triangular<R2: Dim, C2: Dim, S2>( &self, b: &Matrix<T, R2, C2, S2> ) -> Option<OMatrix<T, R2, C2>>where S2: Storage<T, R2, C2>, DefaultAllocator: Allocator<T, 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.

source

pub fn ad_solve_lower_triangular_mut<R2: Dim, C2: Dim, S2>( &self, b: &mut Matrix<T, R2, C2, S2> ) -> boolwhere S2: StorageMut<T, 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.

source

pub fn ad_solve_upper_triangular_mut<R2: Dim, C2: Dim, S2>( &self, b: &mut Matrix<T, R2, C2, S2> ) -> boolwhere S2: StorageMut<T, 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.

source§

impl<T: SimdComplexField, D: Dim, S: Storage<T, D, D>> SquareMatrix<T, D, S>

source

pub fn solve_lower_triangular_unchecked<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>,

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.

source

pub fn solve_upper_triangular_unchecked<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>,

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.

source

pub fn solve_lower_triangular_unchecked_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 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.

source

pub fn solve_lower_triangular_with_diag_unchecked_mut<R2: Dim, C2: Dim, S2>( &self, b: &mut Matrix<T, R2, C2, S2>, diag: T )where S2: StorageMut<T, 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.

source

pub fn solve_upper_triangular_unchecked_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 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.

source

pub fn tr_solve_lower_triangular_unchecked<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>,

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.

source

pub fn tr_solve_upper_triangular_unchecked<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>,

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.

source

pub fn tr_solve_lower_triangular_unchecked_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 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.

source

pub fn tr_solve_upper_triangular_unchecked_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 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.

source

pub fn ad_solve_lower_triangular_unchecked<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>,

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.

source

pub fn ad_solve_upper_triangular_unchecked<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>,

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.

source

pub fn ad_solve_lower_triangular_unchecked_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 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.

source

pub fn ad_solve_upper_triangular_unchecked_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 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.

source§

impl<T: ComplexField, D: DimSub<U1>, S: Storage<T, D, D>> SquareMatrix<T, D, S>where DefaultAllocator: Allocator<T, D, D> + Allocator<T, DimDiff<D, U1>> + Allocator<T::RealField, D> + Allocator<T::RealField, DimDiff<D, U1>>,

source

pub fn symmetric_eigenvalues(&self) -> OVector<T::RealField, D>

Computes the eigenvalues of this symmetric matrix.

Only the lower-triangular part of the matrix is read.