Trait linxal::prelude::LinxalMatrix
[−]
[src]
pub trait LinxalMatrix<F: LinxalScalar> { fn eigenvalues(&self) -> Result<Array<F::Complex, Ix1>, EigenError>; fn eigenvalues_vectors(
&self,
compute_left: bool,
compute_right: bool
) -> Result<Solution<F, F::Complex>, EigenError>; fn symmetric_eigenvalues(
&self,
uplo: Symmetric
) -> Result<Array<F::RealPart, Ix1>, EigenError>; fn symmetric_eigenvalues_vectors(
&self,
uplo: Symmetric
) -> Result<Solution<F, F::RealPart>, EigenError>; fn solve_linear<D1: Data<Elem = F>>(
&self,
b: &ArrayBase<D1, Ix1>
) -> Result<Array<F, Ix1>, SolveError>; fn solve_symmetric_linear<D1: Data<Elem = F>>(
&self,
b: &ArrayBase<D1, Ix1>,
uplo: Symmetric
) -> Result<Array<F, Ix1>, SolveError>; fn solve_multi_linear<D1: Data<Elem = F>>(
&self,
b: &ArrayBase<D1, Ix2>
) -> Result<Array<F, Ix2>, SolveError>; fn solve_symmetric_multi_linear<D1: Data<Elem = F>>(
&self,
b: &ArrayBase<D1, Ix2>,
uplo: Symmetric
) -> Result<Array<F, Ix2>, SolveError>; fn least_squares<D1, PT>(
&self,
b: &ArrayBase<D1, Ix1>,
problem_type: PT
) -> Result<LeastSquaresSolution<F, Ix1>, LeastSquaresError>
where
D1: Data<Elem = F>,
PT: Into<Option<LeastSquaresType>>; fn multi_least_squares<D1, PT>(
&self,
b: &ArrayBase<D1, Ix2>,
problem_type: PT
) -> Result<LeastSquaresSolution<F, Ix2>, LeastSquaresError>
where
D1: Data<Elem = F>,
PT: Into<Option<LeastSquaresType>>; fn qr(&self) -> Result<QRFactors<F>, QRError>; fn lu(&self) -> Result<LUFactors<F>, LUError>; fn cholesky(&self, uplo: Symmetric) -> Result<Array<F, Ix2>, CholeskyError>; fn svd_full(&self) -> Result<SVDSolution<F>, SVDError>; fn svd_econ(&self) -> Result<SVDSolution<F>, SVDError>; fn singular_values(&self) -> Result<Array<F::RealPart, Ix1>, SVDError>; fn inverse(&self) -> Result<Array<F, Ix2>, Error>; fn conj(&self) -> Array<F, Ix2>; fn conj_t(&self) -> Array<F, Ix2>; fn is_square(&self) -> bool; fn is_identity<T: Into<Option<F::RealPart>>>(&self, tolerance: T) -> bool; fn is_diagonal<T: Into<Option<F::RealPart>>>(&self, tolerance: T) -> bool; fn is_symmetric<T: Into<Option<F::RealPart>>>(&self, tolerance: T) -> bool; fn is_unitary<T: Into<Option<F::RealPart>>>(&self, tolerance: T) -> bool; fn is_triangular<T: Into<Option<F::RealPart>>>(
&self,
uplo: Symmetric,
tolerance: T
) -> bool; }
All-encompassing matrix trait, supporting all of the linear
algebra operations defined for any LinxalScalar
.
Required Methods
fn eigenvalues(&self) -> Result<Array<F::Complex, Ix1>, EigenError>
Compute the eigenvalues of a matrix.
fn eigenvalues_vectors(
&self,
compute_left: bool,
compute_right: bool
) -> Result<Solution<F, F::Complex>, EigenError>
&self,
compute_left: bool,
compute_right: bool
) -> Result<Solution<F, F::Complex>, EigenError>
Compute the eigenvalues and eigenvectors of a matrix.
fn symmetric_eigenvalues(
&self,
uplo: Symmetric
) -> Result<Array<F::RealPart, Ix1>, EigenError>
&self,
uplo: Symmetric
) -> Result<Array<F::RealPart, Ix1>, EigenError>
Compute the eigenvalues of a symmetric matrix.
fn symmetric_eigenvalues_vectors(
&self,
uplo: Symmetric
) -> Result<Solution<F, F::RealPart>, EigenError>
&self,
uplo: Symmetric
) -> Result<Solution<F, F::RealPart>, EigenError>
Compute the eigenvalues and eigenvectors of a symmetric matrix.
fn solve_linear<D1: Data<Elem = F>>(
&self,
b: &ArrayBase<D1, Ix1>
) -> Result<Array<F, Ix1>, SolveError>
&self,
b: &ArrayBase<D1, Ix1>
) -> Result<Array<F, Ix1>, SolveError>
Solve a single system of linear equations.
fn solve_symmetric_linear<D1: Data<Elem = F>>(
&self,
b: &ArrayBase<D1, Ix1>,
uplo: Symmetric
) -> Result<Array<F, Ix1>, SolveError>
&self,
b: &ArrayBase<D1, Ix1>,
uplo: Symmetric
) -> Result<Array<F, Ix1>, SolveError>
Solve a single system of linear equations with a symmetrix coefficient matrix.
fn solve_multi_linear<D1: Data<Elem = F>>(
&self,
b: &ArrayBase<D1, Ix2>
) -> Result<Array<F, Ix2>, SolveError>
&self,
b: &ArrayBase<D1, Ix2>
) -> Result<Array<F, Ix2>, SolveError>
Solve a system of linear equations with multiple RHS vectors.
fn solve_symmetric_multi_linear<D1: Data<Elem = F>>(
&self,
b: &ArrayBase<D1, Ix2>,
uplo: Symmetric
) -> Result<Array<F, Ix2>, SolveError>
&self,
b: &ArrayBase<D1, Ix2>,
uplo: Symmetric
) -> Result<Array<F, Ix2>, SolveError>
Solve a system of linear equations with a symmetrix coefficient matrix for multiple RHS vectors.
Each column of b
is a RHS vector to be solved for.
fn least_squares<D1, PT>(
&self,
b: &ArrayBase<D1, Ix1>,
problem_type: PT
) -> Result<LeastSquaresSolution<F, Ix1>, LeastSquaresError> where
D1: Data<Elem = F>,
PT: Into<Option<LeastSquaresType>>,
&self,
b: &ArrayBase<D1, Ix1>,
problem_type: PT
) -> Result<LeastSquaresSolution<F, Ix1>, LeastSquaresError> where
D1: Data<Elem = F>,
PT: Into<Option<LeastSquaresType>>,
Compute the least squares solution for a single RHS.
fn multi_least_squares<D1, PT>(
&self,
b: &ArrayBase<D1, Ix2>,
problem_type: PT
) -> Result<LeastSquaresSolution<F, Ix2>, LeastSquaresError> where
D1: Data<Elem = F>,
PT: Into<Option<LeastSquaresType>>,
&self,
b: &ArrayBase<D1, Ix2>,
problem_type: PT
) -> Result<LeastSquaresSolution<F, Ix2>, LeastSquaresError> where
D1: Data<Elem = F>,
PT: Into<Option<LeastSquaresType>>,
Compute the least squares solution for a multiple RHS.
fn qr(&self) -> Result<QRFactors<F>, QRError>
Return the QR factorization of the matrix.
See QR::compute.
fn lu(&self) -> Result<LUFactors<F>, LUError>
Return the LU factorization of the matrix.
See LU::compute
fn cholesky(&self, uplo: Symmetric) -> Result<Array<F, Ix2>, CholeskyError>
Return the cholesky factorization of the matrix, via the upper- or lower-triangular matrix defining it.
fn svd_full(&self) -> Result<SVDSolution<F>, SVDError>
Return the full singular value decomposition of the matrix.
The SVDSolution
contains full size matrices u
(m x m) and vt
(n x n).
fn svd_econ(&self) -> Result<SVDSolution<F>, SVDError>
Return the economic singular value decomposition of the matrix.
The SVDSolution
contains sufficient matrices u
(m x p) and
vt
(p x n), where p
is the minimum of m
and n
.
fn singular_values(&self) -> Result<Array<F::RealPart, Ix1>, SVDError>
Return the full singular value decomposition of the matrix.
The SVDSolution
contains full size matrices u
(m x m) and vt
(n x n).
fn inverse(&self) -> Result<Array<F, Ix2>, Error>
Return the inverse of the matrix, if it has one.
fn conj(&self) -> Array<F, Ix2>
Return the conjugate of the matrix.
fn conj_t(&self) -> Array<F, Ix2>
Return the cojugate transpose of the matrix.
fn is_square(&self) -> bool
Returns true iff the matrix is square.
fn is_identity<T: Into<Option<F::RealPart>>>(&self, tolerance: T) -> bool
Returns true iff the matrix is the identity matrix, within an optional tolerance.
fn is_diagonal<T: Into<Option<F::RealPart>>>(&self, tolerance: T) -> bool
Returns true iff the matrix is diagnoal, within a certain tolerance.
fn is_symmetric<T: Into<Option<F::RealPart>>>(&self, tolerance: T) -> bool
Returns true iff the matrix is symmetric (or Hermitian, in the complex case), within an optional tolerance.
fn is_unitary<T: Into<Option<F::RealPart>>>(&self, tolerance: T) -> bool
Returns true iff the matrix is unitary, within an optional tolerance.
fn is_triangular<T: Into<Option<F::RealPart>>>(
&self,
uplo: Symmetric,
tolerance: T
) -> bool
&self,
uplo: Symmetric,
tolerance: T
) -> bool
Returns true iff the matrix is triangular or trapezoidal, with
size specified by uplo
.
Implementors
impl<F: LinxalScalar, D: Data<Elem = F>> LinxalMatrix<F> for ArrayBase<D, Ix2>