Type Definition nalgebra::base::SquareMatrix[][src]

type SquareMatrix<N, D, S> = Matrix<N, D, D, S>;

A square matrix.

Implementations

impl<N, D1: Dim, S: StorageMut<N, D1, D1>> SquareMatrix<N, D1, S> where
    N: Scalar + Zero + One + ClosedAdd + ClosedMul
[src]

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>, 
[src]

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>, 
[src]

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>, 
[src]

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>, 
[src]

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);

impl<N: Scalar + Ring, D: DimName, S: Storage<N, D, D>> SquareMatrix<N, D, S>[src]

pub fn append_scaling(&self, scaling: N) -> MatrixN<N, D> where
    D: DimNameSub<U1>,
    DefaultAllocator: Allocator<N, D, D>, 
[src]

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>, 
[src]

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>, 
[src]

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>, 
[src]

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>, 
[src]

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>>, 
[src]

Computes the transformation equal to a translation followed by self.

impl<N: Scalar + Ring, D: DimName, S: StorageMut<N, D, D>> SquareMatrix<N, D, S>[src]

pub fn append_scaling_mut(&mut self, scaling: N) where
    D: DimNameSub<U1>, 
[src]

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

pub fn prepend_scaling_mut(&mut self, scaling: N) where
    D: DimNameSub<U1>, 
[src]

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
    D: DimNameSub<U1>,
    SB: Storage<N, DimNameDiff<D, U1>>, 
[src]

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
    D: DimNameSub<U1>,
    SB: Storage<N, DimNameDiff<D, U1>>, 
[src]

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
    D: DimNameSub<U1>,
    SB: Storage<N, DimNameDiff<D, U1>>, 
[src]

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>,
    SB: Storage<N, DimNameDiff<D, U1>>,
    DefaultAllocator: Allocator<N, DimNameDiff<D, U1>>, 
[src]

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>>, 
[src]

pub fn transform_vector(
    &self,
    v: &VectorN<N, DimNameDiff<D, U1>>
) -> VectorN<N, DimNameDiff<D, U1>>
[src]

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>>
[src]

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

impl<N: Scalar, D: Dim, S: Storage<N, D, D>> SquareMatrix<N, D, S>[src]

pub fn diagonal(&self) -> VectorN<N, D> where
    DefaultAllocator: Allocator<N, D>, 
[src]

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>, 
[src]

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.

pub fn trace(&self) -> N where
    N: Ring
[src]

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

impl<N: ComplexField, D: Dim, S: Storage<N, D, D>> SquareMatrix<N, D, S>[src]

pub fn symmetric_part(&self) -> MatrixMN<N, D, D> where
    DefaultAllocator: Allocator<N, D, D>, 
[src]

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

pub fn hermitian_part(&self) -> MatrixMN<N, D, D> where
    DefaultAllocator: Allocator<N, D, D>, 
[src]

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>, 
[src]

pub fn is_special_orthogonal(&self, eps: N) -> bool where
    D: DimMin<D, Output = D>,
    DefaultAllocator: Allocator<(usize, usize), D>, 
[src]

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

pub fn is_invertible(&self) -> bool[src]

Returns true if this matrix is invertible.

impl<N: ComplexField, D: DimSub<Dynamic>, S: Storage<N, D, D>> SquareMatrix<N, D, S> where
    DefaultAllocator: Allocator<N, D, D>, 
[src]

pub fn cholesky(self) -> Option<Cholesky<N, D>>[src]

Attempts to compute the Cholesky decomposition of this 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.

impl<N: ComplexField, D: DimMin<D, Output = D>, S: Storage<N, D, D>> SquareMatrix<N, D, S>[src]

pub fn determinant(&self) -> N where
    DefaultAllocator: Allocator<N, D, D> + Allocator<(usize, usize), D>, 
[src]

Computes the matrix determinant.

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

impl<N: ComplexField, D: DimSub<U1>, S: Storage<N, D, D>> SquareMatrix<N, D, S> where
    DefaultAllocator: Allocator<N, D, D> + Allocator<N, D> + Allocator<N, DimDiff<D, U1>>, 
[src]

pub fn hessenberg(self) -> Hessenberg<N, D>[src]

Computes the Hessenberg decomposition of this matrix using householder reflections.

impl<N: ComplexField, D: Dim, S: Storage<N, D, D>> SquareMatrix<N, D, S>[src]

pub fn try_inverse(self) -> Option<MatrixN<N, D>> where
    DefaultAllocator: Allocator<N, D, D>, 
[src]

Attempts to invert this matrix.

impl<N: ComplexField, D: Dim, S: StorageMut<N, D, D>> SquareMatrix<N, D, S>[src]

pub fn try_inverse_mut(&mut self) -> bool where
    DefaultAllocator: Allocator<N, D, D>, 
[src]

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

impl<N: ComplexField, D: Dim, S: Storage<N, D, D>> SquareMatrix<N, D, S> where
    D: DimSub<U1>,
    DefaultAllocator: Allocator<N, D, DimDiff<D, U1>> + Allocator<N, DimDiff<D, U1>> + Allocator<N, D, D> + Allocator<N, D>, 
[src]

pub fn schur(self) -> Schur<N, D>[src]

Computes the Schur decomposition of a square matrix.

pub fn try_schur(
    self,
    eps: N::RealField,
    max_niter: usize
) -> Option<Schur<N, D>>
[src]

Attempts to compute the Schur decomposition of a square matrix.

If only eigenvalues are needed, it is more efficient to call the matrix method .eigenvalues() instead.

Arguments

  • eps − tolerance used to determine when a value converged to 0.
  • max_niter − maximum total number of iterations performed by the algorithm. If this number of iteration is exceeded, None is returned. If niter == 0, then the algorithm continues indefinitely until convergence.

pub fn eigenvalues(&self) -> Option<VectorN<N, D>>[src]

Computes the eigenvalues of this matrix.

pub fn complex_eigenvalues(&self) -> VectorN<NumComplex<N>, D> where
    N: RealField,
    DefaultAllocator: Allocator<NumComplex<N>, D>, 
[src]

Computes the eigenvalues of this matrix.

impl<N: ComplexField, D: Dim, S: Storage<N, D, D>> SquareMatrix<N, D, S>[src]

pub fn solve_lower_triangular<R2: Dim, C2: Dim, S2>(
    &self,
    b: &Matrix<N, R2, C2, S2>
) -> Option<MatrixMN<N, R2, C2>> where
    S2: StorageMut<N, R2, C2>,
    DefaultAllocator: Allocator<N, R2, C2>,
    ShapeConstraint: SameNumberOfRows<R2, D>, 
[src]

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: StorageMut<N, R2, C2>,
    DefaultAllocator: Allocator<N, R2, C2>,
    ShapeConstraint: SameNumberOfRows<R2, D>, 
[src]

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>, 
[src]

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>, 
[src]

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>, 
[src]

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: StorageMut<N, R2, C2>,
    DefaultAllocator: Allocator<N, R2, C2>,
    ShapeConstraint: SameNumberOfRows<R2, D>, 
[src]

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: StorageMut<N, R2, C2>,
    DefaultAllocator: Allocator<N, R2, C2>,
    ShapeConstraint: SameNumberOfRows<R2, D>, 
[src]

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>, 
[src]

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>, 
[src]

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: StorageMut<N, R2, C2>,
    DefaultAllocator: Allocator<N, R2, C2>,
    ShapeConstraint: SameNumberOfRows<R2, D>, 
[src]

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: StorageMut<N, R2, C2>,
    DefaultAllocator: Allocator<N, R2, C2>,
    ShapeConstraint: SameNumberOfRows<R2, D>, 
[src]

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>, 
[src]

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>, 
[src]

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.

impl<N: ComplexField, D: DimSub<U1>, S: Storage<N, D, D>> SquareMatrix<N, D, S> where
    DefaultAllocator: Allocator<N, D, D> + Allocator<N, DimDiff<D, U1>> + Allocator<N::RealField, D> + Allocator<N::RealField, DimDiff<D, U1>>, 
[src]

pub fn symmetric_eigen(self) -> SymmetricEigen<N, D>[src]

Computes the eigendecomposition of this symmetric matrix.

Only the lower-triangular part (including the diagonal) of m is read.

pub fn try_symmetric_eigen(
    self,
    eps: N::RealField,
    max_niter: usize
) -> Option<SymmetricEigen<N, D>>
[src]

Computes the eigendecomposition of the given symmetric matrix with user-specified convergence parameters.

Only the lower-triangular part (including the diagonal) of m is read.

Arguments

  • eps − tolerance used to determine when a value converged to 0.
  • max_niter − maximum total number of iterations performed by the algorithm. If this number of iteration is exceeded, None is returned. If niter == 0, then the algorithm continues indefinitely until convergence.

pub fn symmetric_eigenvalues(&self) -> VectorN<N::RealField, D>[src]

Computes the eigenvalues of this symmetric matrix.

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

impl<N: ComplexField, D: DimSub<U1>, S: Storage<N, D, D>> SquareMatrix<N, D, S> where
    DefaultAllocator: Allocator<N, D, D> + Allocator<N, DimDiff<D, U1>>, 
[src]

pub fn symmetric_tridiagonalize(self) -> SymmetricTridiagonal<N, D>[src]

Computes the tridiagonalization of this symmetric matrix.

Only the lower-triangular part (including the diagonal) of m is read.