Type Definition nalgebra::base::SquareMatrix

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

A square matrix.

Implementations§

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

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

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

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

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

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

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

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

Computes the transformation equal to self followed by a translation.

Computes the transformation equal to a translation followed by self.

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

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

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

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

Computes the transformation equal to self followed by a translation.

Computes the transformation equal to a translation followed by self.

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

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

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

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

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

Returns true if this matrix is invertible.

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.

Computes the matrix determinant.

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

Computes the Hessenberg decomposition of this matrix using householder reflections.

Attempts to invert this matrix.

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

Computes the Schur decomposition of a square matrix.

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.

Computes the eigenvalues of this matrix.

Computes the eigenvalues of this matrix.

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.

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.

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.

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.

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.

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.

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.

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.

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.

Computes the eigendecomposition of this symmetric matrix.

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

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.

Computes the eigenvalues of this symmetric matrix.

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

Computes the tridiagonalization of this symmetric matrix.

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