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

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

A square matrix.

Methods

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

[src]

Creates a square matrix with its diagonal set to diag and all other entries set to 0.

[src]

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

impl<N: Real, D: Dim, S: Storage<N, D, D>> SquareMatrix<N, D, S> where
    DefaultAllocator: Allocator<N, D, D>, 
[src]

[src]

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

[src]

Returns true if this matrix is invertible.

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

[src]

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

[src]

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

[src]

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

[src]

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

[src]

Computes the transformation equal to self followed by a translation.

[src]

Computes the transformation equal to a translation followed by self.

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

[src]

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

[src]

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

[src]

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

[src]

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

[src]

Computes the transformation equal to self followed by a translation.

[src]

Computes the transformation equal to a translation followed by self.

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

[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 concidered not-zero.

[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 concidered not-zero.

[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 concidered not-zero.

[src]

Solves the linear system self . x = b where x is the unknown and only the lower-triangular part of self is concidered 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.

[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 concidered not-zero.

[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 concidered not-zero.

[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 concidered not-zero.

[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 concidered not-zero.

[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 concidered not-zero.

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

[src]

Computes the matrix determinant.

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

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

[src]

Attempts to invert this matrix.

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

[src]

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

impl<N: Real, 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]

[src]

Computes the Hessenberg decomposition of this matrix using householder reflections.

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

[src]

Computes the tridiagonalization of this symmetric matrix.

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

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

[src]

Attempts to compute the Cholesky decomposition of this matrix.

Returns None if the input matrix is not definite-positive. The intput matrix is assumed to be symmetric and only the lower-triangular part is read.

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

[src]

Computes the Schur decomposition of a square matrix.

[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 − tolerence 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.

[src]

Computes the eigenvalues of this matrix.

[src]

Computes the eigenvalues of this matrix.

impl<N: Real, 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]

[src]

Computes the eigendecomposition of this symmetric matrix.

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

[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.

[src]

Computes the eigenvalues of this symmetric matrix.

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