Trait ndarray_linalg::types::Lapack

source ·
pub trait Lapack: Scalar {
Show 25 methods // Required methods fn eig( calc_v: bool, l: MatrixLayout, a: &mut [Self], ) -> Result<(Vec<Self::Complex>, Vec<Self::Complex>), Error>; fn eigh( calc_eigenvec: bool, layout: MatrixLayout, uplo: UPLO, a: &mut [Self], ) -> Result<Vec<Self::Real>, Error>; fn eigh_generalized( calc_eigenvec: bool, layout: MatrixLayout, uplo: UPLO, a: &mut [Self], b: &mut [Self], ) -> Result<Vec<Self::Real>, Error>; fn householder(l: MatrixLayout, a: &mut [Self]) -> Result<Vec<Self>, Error>; fn q(l: MatrixLayout, a: &mut [Self], tau: &[Self]) -> Result<(), Error>; fn qr(l: MatrixLayout, a: &mut [Self]) -> Result<Vec<Self>, Error>; fn svd( l: MatrixLayout, calc_u: bool, calc_vt: bool, a: &mut [Self], ) -> Result<SvdOwned<Self>, Error>; fn svddc( layout: MatrixLayout, jobz: JobSvd, a: &mut [Self], ) -> Result<SvdOwned<Self>, Error>; fn least_squares( a_layout: MatrixLayout, a: &mut [Self], b: &mut [Self], ) -> Result<LeastSquaresOwned<Self>, Error>; fn least_squares_nrhs( a_layout: MatrixLayout, a: &mut [Self], b_layout: MatrixLayout, b: &mut [Self], ) -> Result<LeastSquaresOwned<Self>, Error>; fn lu(l: MatrixLayout, a: &mut [Self]) -> Result<Vec<i32>, Error>; fn inv(l: MatrixLayout, a: &mut [Self], p: &Vec<i32>) -> Result<(), Error>; fn solve( l: MatrixLayout, t: Transpose, a: &[Self], p: &Vec<i32>, b: &mut [Self], ) -> Result<(), Error>; fn bk( l: MatrixLayout, uplo: UPLO, a: &mut [Self], ) -> Result<Vec<i32>, Error>; fn invh( l: MatrixLayout, uplo: UPLO, a: &mut [Self], ipiv: &Vec<i32>, ) -> Result<(), Error>; fn solveh( l: MatrixLayout, uplo: UPLO, a: &[Self], ipiv: &Vec<i32>, b: &mut [Self], ) -> Result<(), Error>; fn cholesky( l: MatrixLayout, uplo: UPLO, a: &mut [Self], ) -> Result<(), Error>; fn inv_cholesky( l: MatrixLayout, uplo: UPLO, a: &mut [Self], ) -> Result<(), Error>; fn solve_cholesky( l: MatrixLayout, uplo: UPLO, a: &[Self], b: &mut [Self], ) -> Result<(), Error>; fn rcond( l: MatrixLayout, a: &[Self], anorm: Self::Real, ) -> Result<Self::Real, Error>; fn opnorm(t: NormType, l: MatrixLayout, a: &[Self]) -> Self::Real; fn solve_triangular( al: MatrixLayout, bl: MatrixLayout, uplo: UPLO, d: Diag, a: &[Self], b: &mut [Self], ) -> Result<(), Error>; fn lu_tridiagonal( a: Tridiagonal<Self>, ) -> Result<LUFactorizedTridiagonal<Self>, Error>; fn rcond_tridiagonal( lu: &LUFactorizedTridiagonal<Self>, ) -> Result<Self::Real, Error>; fn solve_tridiagonal( lu: &LUFactorizedTridiagonal<Self>, bl: MatrixLayout, t: Transpose, b: &mut [Self], ) -> Result<(), Error>;
}
Expand description

Trait for primitive types which implements LAPACK subroutines

Required Methods§

source

fn eig( calc_v: bool, l: MatrixLayout, a: &mut [Self], ) -> Result<(Vec<Self::Complex>, Vec<Self::Complex>), Error>

Compute right eigenvalue and eigenvectors for a general matrix

source

fn eigh( calc_eigenvec: bool, layout: MatrixLayout, uplo: UPLO, a: &mut [Self], ) -> Result<Vec<Self::Real>, Error>

Compute right eigenvalue and eigenvectors for a symmetric or Hermitian matrix

source

fn eigh_generalized( calc_eigenvec: bool, layout: MatrixLayout, uplo: UPLO, a: &mut [Self], b: &mut [Self], ) -> Result<Vec<Self::Real>, Error>

Compute right eigenvalue and eigenvectors for a symmetric or Hermitian matrix

source

fn householder(l: MatrixLayout, a: &mut [Self]) -> Result<Vec<Self>, Error>

Execute Householder reflection as the first step of QR-decomposition

For C-continuous array, this will call LQ-decomposition of the transposed matrix $ A^T = LQ^T $

source

fn q(l: MatrixLayout, a: &mut [Self], tau: &[Self]) -> Result<(), Error>

Reconstruct Q-matrix from Householder-reflectors

source

fn qr(l: MatrixLayout, a: &mut [Self]) -> Result<Vec<Self>, Error>

Execute QR-decomposition at once

source

fn svd( l: MatrixLayout, calc_u: bool, calc_vt: bool, a: &mut [Self], ) -> Result<SvdOwned<Self>, Error>

Compute singular-value decomposition (SVD)

source

fn svddc( layout: MatrixLayout, jobz: JobSvd, a: &mut [Self], ) -> Result<SvdOwned<Self>, Error>

Compute singular value decomposition (SVD) with divide-and-conquer algorithm

source

fn least_squares( a_layout: MatrixLayout, a: &mut [Self], b: &mut [Self], ) -> Result<LeastSquaresOwned<Self>, Error>

Compute a vector $x$ which minimizes Euclidian norm $| Ax - b|$ for a given matrix $A$ and a vector $b$.

source

fn least_squares_nrhs( a_layout: MatrixLayout, a: &mut [Self], b_layout: MatrixLayout, b: &mut [Self], ) -> Result<LeastSquaresOwned<Self>, Error>

Solve least square problems $\argmin_X | AX - B|$

source

fn lu(l: MatrixLayout, a: &mut [Self]) -> Result<Vec<i32>, Error>

Computes the LU decomposition of a general $m \times n$ matrix with partial pivoting with row interchanges.

For a given matrix $A$, LU decomposition is described as $A = PLU$ where:

  • $L$ is lower matrix
  • $U$ is upper matrix
  • $P$ is permutation matrix represented by Pivot

This is designed as two step computation according to LAPACK API:

  1. Factorize input matrix $A$ into $L$, $U$, and $P$.
  2. Solve linear equation $Ax = b$ by Lapack::solve or compute inverse matrix $A^{-1}$ by Lapack::inv using the output of LU decomposition.
§Output
  • $U$ and $L$ are stored in a after LU decomposition has succeeded.
  • $P$ is returned as Pivot
§Error
source

fn inv(l: MatrixLayout, a: &mut [Self], p: &Vec<i32>) -> Result<(), Error>

Compute inverse matrix $A^{-1}$ from the output of LU-decomposition

source

fn solve( l: MatrixLayout, t: Transpose, a: &[Self], p: &Vec<i32>, b: &mut [Self], ) -> Result<(), Error>

Solve linear equations $Ax = b$ using the output of LU-decomposition

source

fn bk(l: MatrixLayout, uplo: UPLO, a: &mut [Self]) -> Result<Vec<i32>, Error>

Factorize symmetric/Hermitian matrix using Bunch-Kaufman diagonal pivoting method

For a given symmetric matrix $A$, this method factorizes $A = U^T D U$ or $A = L D L^T$ where

  • $U$ (or $L$) are is a product of permutation and unit upper (lower) triangular matrices
  • $D$ is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.

This takes two-step approach based in LAPACK:

  1. Factorize given matrix $A$ into upper ($U$) or lower ($L$) form with diagonal matrix $D$
  2. Then solve linear equation $Ax = b$, and/or calculate inverse matrix $A^{-1}$
source

fn invh( l: MatrixLayout, uplo: UPLO, a: &mut [Self], ipiv: &Vec<i32>, ) -> Result<(), Error>

Compute inverse matrix $A^{-1}$ using the result of Lapack::bk

source

fn solveh( l: MatrixLayout, uplo: UPLO, a: &[Self], ipiv: &Vec<i32>, b: &mut [Self], ) -> Result<(), Error>

Solve symmetric/Hermitian linear equation $Ax = b$ using the result of Lapack::bk

source

fn cholesky(l: MatrixLayout, uplo: UPLO, a: &mut [Self]) -> Result<(), Error>

Solve symmetric/Hermitian positive-definite linear equations using Cholesky decomposition

For a given positive definite matrix $A$, Cholesky decomposition is described as $A = U^T U$ or $A = LL^T$ where

  • $L$ is lower matrix
  • $U$ is upper matrix

This is designed as two step computation according to LAPACK API

  1. Factorize input matrix $A$ into $L$ or $U$
  2. Solve linear equation $Ax = b$ by Lapack::solve_cholesky or compute inverse matrix $A^{-1}$ by Lapack::inv_cholesky
source

fn inv_cholesky( l: MatrixLayout, uplo: UPLO, a: &mut [Self], ) -> Result<(), Error>

Compute inverse matrix $A^{-1}$ using $U$ or $L$ calculated by Lapack::cholesky

source

fn solve_cholesky( l: MatrixLayout, uplo: UPLO, a: &[Self], b: &mut [Self], ) -> Result<(), Error>

Solve linear equation $Ax = b$ using $U$ or $L$ calculated by Lapack::cholesky

source

fn rcond( l: MatrixLayout, a: &[Self], anorm: Self::Real, ) -> Result<Self::Real, Error>

Estimates the the reciprocal of the condition number of the matrix in 1-norm.

anorm should be the 1-norm of the matrix a.

source

fn opnorm(t: NormType, l: MatrixLayout, a: &[Self]) -> Self::Real

Compute norm of matrices

For a $n \times m$ matrix $$ A = \begin{pmatrix} a_{11} & \cdots & a_{1m} \\ \vdots & \ddots & \vdots \\ a_{n1} & \cdots & a_{nm} \end{pmatrix} $$ LAPACK can compute three types of norms:

  • Operator norm based on 1-norm for its domain linear space: $$ \Vert A \Vert_1 = \sup_{\Vert x \Vert_1 = 1} \Vert Ax \Vert_1 = \max_{1 \le j \le m } \sum_{i=1}^n |a_{ij}| $$ where $\Vert x\Vert_1 = \sum_{j=1}^m |x_j|$ is 1-norm for a vector $x$.

  • Operator norm based on $\infty$-norm for its domain linear space: $$ \Vert A \Vert_\infty = \sup_{\Vert x \Vert_\infty = 1} \Vert Ax \Vert_\infty = \max_{1 \le i \le n } \sum_{j=1}^m |a_{ij}| $$ where $\Vert x\Vert_\infty = \max_{j=1}^m |x_j|$ is $\infty$-norm for a vector $x$.

  • Frobenious norm $$ \Vert A \Vert_F = \sqrt{\mathrm{Tr} \left(AA^\dagger\right)} = \sqrt{\sum_{i=1}^n \sum_{j=1}^m |a_{ij}|^2} $$

source

fn solve_triangular( al: MatrixLayout, bl: MatrixLayout, uplo: UPLO, d: Diag, a: &[Self], b: &mut [Self], ) -> Result<(), Error>

source

fn lu_tridiagonal( a: Tridiagonal<Self>, ) -> Result<LUFactorizedTridiagonal<Self>, Error>

Computes the LU factorization of a tridiagonal m x n matrix a using partial pivoting with row interchanges.

source

fn rcond_tridiagonal( lu: &LUFactorizedTridiagonal<Self>, ) -> Result<Self::Real, Error>

source

fn solve_tridiagonal( lu: &LUFactorizedTridiagonal<Self>, bl: MatrixLayout, t: Transpose, b: &mut [Self], ) -> Result<(), Error>

Object Safety§

This trait is not object safe.

Implementations on Foreign Types§

source§

impl Lapack for f32

source§

fn eig( calc_v: bool, l: MatrixLayout, a: &mut [f32], ) -> Result<(Vec<<f32 as Scalar>::Complex>, Vec<<f32 as Scalar>::Complex>), Error>

source§

fn eigh( calc_eigenvec: bool, layout: MatrixLayout, uplo: UPLO, a: &mut [f32], ) -> Result<Vec<<f32 as Scalar>::Real>, Error>

source§

fn eigh_generalized( calc_eigenvec: bool, layout: MatrixLayout, uplo: UPLO, a: &mut [f32], b: &mut [f32], ) -> Result<Vec<<f32 as Scalar>::Real>, Error>

source§

fn householder(l: MatrixLayout, a: &mut [f32]) -> Result<Vec<f32>, Error>

source§

fn q(l: MatrixLayout, a: &mut [f32], tau: &[f32]) -> Result<(), Error>

source§

fn qr(l: MatrixLayout, a: &mut [f32]) -> Result<Vec<f32>, Error>

source§

fn svd( l: MatrixLayout, calc_u: bool, calc_vt: bool, a: &mut [f32], ) -> Result<SvdOwned<f32>, Error>

source§

fn svddc( layout: MatrixLayout, jobz: JobSvd, a: &mut [f32], ) -> Result<SvdOwned<f32>, Error>

source§

fn least_squares( l: MatrixLayout, a: &mut [f32], b: &mut [f32], ) -> Result<LeastSquaresOwned<f32>, Error>

source§

fn least_squares_nrhs( a_layout: MatrixLayout, a: &mut [f32], b_layout: MatrixLayout, b: &mut [f32], ) -> Result<LeastSquaresOwned<f32>, Error>

source§

fn lu(l: MatrixLayout, a: &mut [f32]) -> Result<Vec<i32>, Error>

source§

fn inv(l: MatrixLayout, a: &mut [f32], p: &Vec<i32>) -> Result<(), Error>

source§

fn solve( l: MatrixLayout, t: Transpose, a: &[f32], p: &Vec<i32>, b: &mut [f32], ) -> Result<(), Error>

source§

fn bk(l: MatrixLayout, uplo: UPLO, a: &mut [f32]) -> Result<Vec<i32>, Error>

source§

fn invh( l: MatrixLayout, uplo: UPLO, a: &mut [f32], ipiv: &Vec<i32>, ) -> Result<(), Error>

source§

fn solveh( l: MatrixLayout, uplo: UPLO, a: &[f32], ipiv: &Vec<i32>, b: &mut [f32], ) -> Result<(), Error>

source§

fn cholesky(l: MatrixLayout, uplo: UPLO, a: &mut [f32]) -> Result<(), Error>

source§

fn inv_cholesky(l: MatrixLayout, uplo: UPLO, a: &mut [f32]) -> Result<(), Error>

source§

fn solve_cholesky( l: MatrixLayout, uplo: UPLO, a: &[f32], b: &mut [f32], ) -> Result<(), Error>

source§

fn rcond( l: MatrixLayout, a: &[f32], anorm: <f32 as Scalar>::Real, ) -> Result<<f32 as Scalar>::Real, Error>

source§

fn opnorm(t: NormType, l: MatrixLayout, a: &[f32]) -> <f32 as Scalar>::Real

source§

fn solve_triangular( al: MatrixLayout, bl: MatrixLayout, uplo: UPLO, d: Diag, a: &[f32], b: &mut [f32], ) -> Result<(), Error>

source§

fn lu_tridiagonal( a: Tridiagonal<f32>, ) -> Result<LUFactorizedTridiagonal<f32>, Error>

source§

fn rcond_tridiagonal( lu: &LUFactorizedTridiagonal<f32>, ) -> Result<<f32 as Scalar>::Real, Error>

source§

fn solve_tridiagonal( lu: &LUFactorizedTridiagonal<f32>, bl: MatrixLayout, t: Transpose, b: &mut [f32], ) -> Result<(), Error>

source§

impl Lapack for f64

source§

fn eig( calc_v: bool, l: MatrixLayout, a: &mut [f64], ) -> Result<(Vec<<f64 as Scalar>::Complex>, Vec<<f64 as Scalar>::Complex>), Error>

source§

fn eigh( calc_eigenvec: bool, layout: MatrixLayout, uplo: UPLO, a: &mut [f64], ) -> Result<Vec<<f64 as Scalar>::Real>, Error>

source§

fn eigh_generalized( calc_eigenvec: bool, layout: MatrixLayout, uplo: UPLO, a: &mut [f64], b: &mut [f64], ) -> Result<Vec<<f64 as Scalar>::Real>, Error>

source§

fn householder(l: MatrixLayout, a: &mut [f64]) -> Result<Vec<f64>, Error>

source§

fn q(l: MatrixLayout, a: &mut [f64], tau: &[f64]) -> Result<(), Error>

source§

fn qr(l: MatrixLayout, a: &mut [f64]) -> Result<Vec<f64>, Error>

source§

fn svd( l: MatrixLayout, calc_u: bool, calc_vt: bool, a: &mut [f64], ) -> Result<SvdOwned<f64>, Error>

source§

fn svddc( layout: MatrixLayout, jobz: JobSvd, a: &mut [f64], ) -> Result<SvdOwned<f64>, Error>

source§

fn least_squares( l: MatrixLayout, a: &mut [f64], b: &mut [f64], ) -> Result<LeastSquaresOwned<f64>, Error>

source§

fn least_squares_nrhs( a_layout: MatrixLayout, a: &mut [f64], b_layout: MatrixLayout, b: &mut [f64], ) -> Result<LeastSquaresOwned<f64>, Error>

source§

fn lu(l: MatrixLayout, a: &mut [f64]) -> Result<Vec<i32>, Error>

source§

fn inv(l: MatrixLayout, a: &mut [f64], p: &Vec<i32>) -> Result<(), Error>

source§

fn solve( l: MatrixLayout, t: Transpose, a: &[f64], p: &Vec<i32>, b: &mut [f64], ) -> Result<(), Error>

source§

fn bk(l: MatrixLayout, uplo: UPLO, a: &mut [f64]) -> Result<Vec<i32>, Error>

source§

fn invh( l: MatrixLayout, uplo: UPLO, a: &mut [f64], ipiv: &Vec<i32>, ) -> Result<(), Error>

source§

fn solveh( l: MatrixLayout, uplo: UPLO, a: &[f64], ipiv: &Vec<i32>, b: &mut [f64], ) -> Result<(), Error>

source§

fn cholesky(l: MatrixLayout, uplo: UPLO, a: &mut [f64]) -> Result<(), Error>

source§

fn inv_cholesky(l: MatrixLayout, uplo: UPLO, a: &mut [f64]) -> Result<(), Error>

source§

fn solve_cholesky( l: MatrixLayout, uplo: UPLO, a: &[f64], b: &mut [f64], ) -> Result<(), Error>

source§

fn rcond( l: MatrixLayout, a: &[f64], anorm: <f64 as Scalar>::Real, ) -> Result<<f64 as Scalar>::Real, Error>

source§

fn opnorm(t: NormType, l: MatrixLayout, a: &[f64]) -> <f64 as Scalar>::Real

source§

fn solve_triangular( al: MatrixLayout, bl: MatrixLayout, uplo: UPLO, d: Diag, a: &[f64], b: &mut [f64], ) -> Result<(), Error>

source§

fn lu_tridiagonal( a: Tridiagonal<f64>, ) -> Result<LUFactorizedTridiagonal<f64>, Error>

source§

fn rcond_tridiagonal( lu: &LUFactorizedTridiagonal<f64>, ) -> Result<<f64 as Scalar>::Real, Error>

source§

fn solve_tridiagonal( lu: &LUFactorizedTridiagonal<f64>, bl: MatrixLayout, t: Transpose, b: &mut [f64], ) -> Result<(), Error>

source§

impl Lapack for Complex<f32>

source§

fn eig( calc_v: bool, l: MatrixLayout, a: &mut [Complex<f32>], ) -> Result<(Vec<<Complex<f32> as Scalar>::Complex>, Vec<<Complex<f32> as Scalar>::Complex>), Error>

source§

fn eigh( calc_eigenvec: bool, layout: MatrixLayout, uplo: UPLO, a: &mut [Complex<f32>], ) -> Result<Vec<<Complex<f32> as Scalar>::Real>, Error>

source§

fn eigh_generalized( calc_eigenvec: bool, layout: MatrixLayout, uplo: UPLO, a: &mut [Complex<f32>], b: &mut [Complex<f32>], ) -> Result<Vec<<Complex<f32> as Scalar>::Real>, Error>

source§

fn householder( l: MatrixLayout, a: &mut [Complex<f32>], ) -> Result<Vec<Complex<f32>>, Error>

source§

fn q( l: MatrixLayout, a: &mut [Complex<f32>], tau: &[Complex<f32>], ) -> Result<(), Error>

source§

fn qr( l: MatrixLayout, a: &mut [Complex<f32>], ) -> Result<Vec<Complex<f32>>, Error>

source§

fn svd( l: MatrixLayout, calc_u: bool, calc_vt: bool, a: &mut [Complex<f32>], ) -> Result<SvdOwned<Complex<f32>>, Error>

source§

fn svddc( layout: MatrixLayout, jobz: JobSvd, a: &mut [Complex<f32>], ) -> Result<SvdOwned<Complex<f32>>, Error>

source§

fn least_squares( l: MatrixLayout, a: &mut [Complex<f32>], b: &mut [Complex<f32>], ) -> Result<LeastSquaresOwned<Complex<f32>>, Error>

source§

fn least_squares_nrhs( a_layout: MatrixLayout, a: &mut [Complex<f32>], b_layout: MatrixLayout, b: &mut [Complex<f32>], ) -> Result<LeastSquaresOwned<Complex<f32>>, Error>

source§

fn lu(l: MatrixLayout, a: &mut [Complex<f32>]) -> Result<Vec<i32>, Error>

source§

fn inv( l: MatrixLayout, a: &mut [Complex<f32>], p: &Vec<i32>, ) -> Result<(), Error>

source§

fn solve( l: MatrixLayout, t: Transpose, a: &[Complex<f32>], p: &Vec<i32>, b: &mut [Complex<f32>], ) -> Result<(), Error>

source§

fn bk( l: MatrixLayout, uplo: UPLO, a: &mut [Complex<f32>], ) -> Result<Vec<i32>, Error>

source§

fn invh( l: MatrixLayout, uplo: UPLO, a: &mut [Complex<f32>], ipiv: &Vec<i32>, ) -> Result<(), Error>

source§

fn solveh( l: MatrixLayout, uplo: UPLO, a: &[Complex<f32>], ipiv: &Vec<i32>, b: &mut [Complex<f32>], ) -> Result<(), Error>

source§

fn cholesky( l: MatrixLayout, uplo: UPLO, a: &mut [Complex<f32>], ) -> Result<(), Error>

source§

fn inv_cholesky( l: MatrixLayout, uplo: UPLO, a: &mut [Complex<f32>], ) -> Result<(), Error>

source§

fn solve_cholesky( l: MatrixLayout, uplo: UPLO, a: &[Complex<f32>], b: &mut [Complex<f32>], ) -> Result<(), Error>

source§

fn rcond( l: MatrixLayout, a: &[Complex<f32>], anorm: <Complex<f32> as Scalar>::Real, ) -> Result<<Complex<f32> as Scalar>::Real, Error>

source§

fn opnorm( t: NormType, l: MatrixLayout, a: &[Complex<f32>], ) -> <Complex<f32> as Scalar>::Real

source§

fn solve_triangular( al: MatrixLayout, bl: MatrixLayout, uplo: UPLO, d: Diag, a: &[Complex<f32>], b: &mut [Complex<f32>], ) -> Result<(), Error>

source§

fn lu_tridiagonal( a: Tridiagonal<Complex<f32>>, ) -> Result<LUFactorizedTridiagonal<Complex<f32>>, Error>

source§

fn rcond_tridiagonal( lu: &LUFactorizedTridiagonal<Complex<f32>>, ) -> Result<<Complex<f32> as Scalar>::Real, Error>

source§

fn solve_tridiagonal( lu: &LUFactorizedTridiagonal<Complex<f32>>, bl: MatrixLayout, t: Transpose, b: &mut [Complex<f32>], ) -> Result<(), Error>

source§

impl Lapack for Complex<f64>

source§

fn eig( calc_v: bool, l: MatrixLayout, a: &mut [Complex<f64>], ) -> Result<(Vec<<Complex<f64> as Scalar>::Complex>, Vec<<Complex<f64> as Scalar>::Complex>), Error>

source§

fn eigh( calc_eigenvec: bool, layout: MatrixLayout, uplo: UPLO, a: &mut [Complex<f64>], ) -> Result<Vec<<Complex<f64> as Scalar>::Real>, Error>

source§

fn eigh_generalized( calc_eigenvec: bool, layout: MatrixLayout, uplo: UPLO, a: &mut [Complex<f64>], b: &mut [Complex<f64>], ) -> Result<Vec<<Complex<f64> as Scalar>::Real>, Error>

source§

fn householder( l: MatrixLayout, a: &mut [Complex<f64>], ) -> Result<Vec<Complex<f64>>, Error>

source§

fn q( l: MatrixLayout, a: &mut [Complex<f64>], tau: &[Complex<f64>], ) -> Result<(), Error>

source§

fn qr( l: MatrixLayout, a: &mut [Complex<f64>], ) -> Result<Vec<Complex<f64>>, Error>

source§

fn svd( l: MatrixLayout, calc_u: bool, calc_vt: bool, a: &mut [Complex<f64>], ) -> Result<SvdOwned<Complex<f64>>, Error>

source§

fn svddc( layout: MatrixLayout, jobz: JobSvd, a: &mut [Complex<f64>], ) -> Result<SvdOwned<Complex<f64>>, Error>

source§

fn least_squares( l: MatrixLayout, a: &mut [Complex<f64>], b: &mut [Complex<f64>], ) -> Result<LeastSquaresOwned<Complex<f64>>, Error>

source§

fn least_squares_nrhs( a_layout: MatrixLayout, a: &mut [Complex<f64>], b_layout: MatrixLayout, b: &mut [Complex<f64>], ) -> Result<LeastSquaresOwned<Complex<f64>>, Error>

source§

fn lu(l: MatrixLayout, a: &mut [Complex<f64>]) -> Result<Vec<i32>, Error>

source§

fn inv( l: MatrixLayout, a: &mut [Complex<f64>], p: &Vec<i32>, ) -> Result<(), Error>

source§

fn solve( l: MatrixLayout, t: Transpose, a: &[Complex<f64>], p: &Vec<i32>, b: &mut [Complex<f64>], ) -> Result<(), Error>

source§

fn bk( l: MatrixLayout, uplo: UPLO, a: &mut [Complex<f64>], ) -> Result<Vec<i32>, Error>

source§

fn invh( l: MatrixLayout, uplo: UPLO, a: &mut [Complex<f64>], ipiv: &Vec<i32>, ) -> Result<(), Error>

source§

fn solveh( l: MatrixLayout, uplo: UPLO, a: &[Complex<f64>], ipiv: &Vec<i32>, b: &mut [Complex<f64>], ) -> Result<(), Error>

source§

fn cholesky( l: MatrixLayout, uplo: UPLO, a: &mut [Complex<f64>], ) -> Result<(), Error>

source§

fn inv_cholesky( l: MatrixLayout, uplo: UPLO, a: &mut [Complex<f64>], ) -> Result<(), Error>

source§

fn solve_cholesky( l: MatrixLayout, uplo: UPLO, a: &[Complex<f64>], b: &mut [Complex<f64>], ) -> Result<(), Error>

source§

fn rcond( l: MatrixLayout, a: &[Complex<f64>], anorm: <Complex<f64> as Scalar>::Real, ) -> Result<<Complex<f64> as Scalar>::Real, Error>

source§

fn opnorm( t: NormType, l: MatrixLayout, a: &[Complex<f64>], ) -> <Complex<f64> as Scalar>::Real

source§

fn solve_triangular( al: MatrixLayout, bl: MatrixLayout, uplo: UPLO, d: Diag, a: &[Complex<f64>], b: &mut [Complex<f64>], ) -> Result<(), Error>

source§

fn lu_tridiagonal( a: Tridiagonal<Complex<f64>>, ) -> Result<LUFactorizedTridiagonal<Complex<f64>>, Error>

source§

fn rcond_tridiagonal( lu: &LUFactorizedTridiagonal<Complex<f64>>, ) -> Result<<Complex<f64> as Scalar>::Real, Error>

source§

fn solve_tridiagonal( lu: &LUFactorizedTridiagonal<Complex<f64>>, bl: MatrixLayout, t: Transpose, b: &mut [Complex<f64>], ) -> Result<(), Error>

Implementors§