lapack #
The V Linear Algebra Package
This package implements Linear Algebra routines in V.
| Backend | Description | Status | Compilation Flags |
|---|---|---|---|
| BLAS | Pure V implementation | WIP | NONE |
| LAPACKE | LAPACKE is a C interface to LAPACK. It is a standard part of the LAPACK distribution. | ||
| Check the section LAPACKE Backend for more information. | Stable | -d vsl_lapack_lapacke |
Therefore, its routines are a little more lower level than the ones in the package vsl.la.
LAPACKE Backend
We provide a backend for the LAPACKE library. This backend is probably the fastest one for all platforms but it requires the installation of the LAPACKE library.
Use the compilation flag -d vsl_lapack_lapacke to use the LAPACKE backend instead of the pure V implementation and make sure that the LAPACKE library is installed in your system.
Check the section below for more information about installing the LAPACKE library.
Install dependencies
Homebrew (macOS)
brew install lapack
Debian/Ubuntu GNU Linux
sudo apt-get install -y --no-install-recommends
Arch Linux/Manjaro GNU Linux
The best way of installing OpenBLAS is using blas-openblas.
sudo pacman -S blas-openblasfn dgeev #
fn dgeev(calc_vl LeftEigenVectorsJob, calc_vr RightEigenVectorsJob, n int, mut a []f64, lda int, mut wr []f64, mut wi []f64, mut vl []f64, ldvl int, mut vr []f64, ldvr int) int
dgeev computes for an N-by-N real nonsymmetric matrix A, the eigenvalues and, optionally, the left and/or right eigenvectors.
See: http://www.netlib.org/lapack/explore-html/d9/d28/dgeev_8f.html
See: https://software.intel.com/en-us/mkl-developer-reference-c-geev
See: https://www.nag.co.uk/numeric/fl/nagdoc_fl26/html/f08/f08naf.html
The right eigenvector v(j) of A satisfies
A * v(j) = lambda(j) * v(j)
where lambda(j) is its eigenvalue.
The left eigenvector u(j) of A satisfies
u(j)**H * A = lambda(j) * u(j)**H
where u(j)**H denotes the conjugate-transpose of u(j).
The computed eigenvectors are normalized to have Euclidean norm equal to 1 and largest component real.
fn dgeev_standardized #
fn dgeev_standardized(jobvl LeftEigenVectorsJob, jobvr RightEigenVectorsJob, n int, mut a []f64, lda int, mut wr []f64, mut wi []f64, mut vl []f64, ldvl int, mut vr []f64, ldvr int) int
dgeev_standardized - Standardized dgeev wrapper
fn dgeqrf #
fn dgeqrf(m int, n int, mut a []f64, lda int, mut tau []f64) int
dgeqrf - Direct wrapper for dgeqrf with standardized interface
fn dgeqrf_standardized #
fn dgeqrf_standardized(m int, n int, mut a []f64, lda int, mut tau []f64) int
dgeqrf_standardized - Standardized dgeqrf wrapper (placeholder until implemented)
fn dgesv #
fn dgesv(n int, nrhs int, mut a []f64, lda int, mut ipiv []int, mut b []f64, ldb int) int
dgesv computes the solution to a real system of linear equations.
See: http://www.netlib.org/lapack/explore-html/d8/d72/dgesv_8f.html
See: https://software.intel.com/en-us/mkl-developer-reference-c-gesv
The system is:
A * X = B,
where A is an N-by-N matrix and X and B are N-by-NRHS matrices.
The LU decomposition with partial pivoting and row interchanges is used to factor A as
A = P * L * U,
where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. The factored form of A is then used to solve the system of equations A * X = B.
Note: matrix 'a' will be modified
fn dgesvd #
fn dgesvd(jobu SVDJob, jobvt SVDJob, m int, n int, mut a []f64, lda int, s []f64, mut u []f64, ldu int, mut vt []f64, ldvt int, superb []f64) int
dgesvd computes the singular value decomposition (SingularValueDecomposition) of a real M-by-N matrix A, optionally computing the left and/or right singular vectors.
See: http://www.netlib.org/lapack/explore-html/d8/d2d/dgesvd_8f.html
See: https://software.intel.com/en-us/mkl-developer-reference-c-gesvd
The SingularValueDecomposition is written
A = U * SIGMA * transpose(V)
where SIGMA is an M-by-N matrix which is zero except for its min(m,n) diagonal elements, U is an M-by-M orthogonal matrix, and V is an N-by-N orthogonal matrix. The diagonal elements of SIGMA are the singular values of A; they are real and non-negative, and are returned in descending order. The first min(m,n) columns of U and V are the left and right singular vectors of A.
Note that the routine returns V**T, not V.
Note: matrix 'a' will be modified
fn dgetrf #
fn dgetrf(m int, n int, mut a []f64, lda int, mut ipiv []int) int
dgetrf computes an LU factorization of a general M-by-N matrix A using partial pivoting with row interchanges.
See: http://www.netlib.org/lapack/explore-html/d3/d6a/dgetrf_8f.html
See: https://software.intel.com/en-us/mkl-developer-reference-c-getrf
The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).
This is the right-looking Level 3 BLAS version of the algorithm.
Note: (1) matrix 'a' will be modified(2) ipiv indices are 1-based (i.e. Fortran)
fn dgetri #
fn dgetri(n int, mut a []f64, lda int, mut ipiv []int) int
dgetri computes the inverse of a matrix using the LU factorization computed by DGETRF.
See: http://www.netlib.org/lapack/explore-html/df/da4/dgetri_8f.html
See: https://software.intel.com/en-us/mkl-developer-reference-c-getri
This method inverts U and then computes inv(A) by solving the system inv(A)*L = inv(U) for inv(A).
fn dlange #
fn dlange(norm rune, m int, n int, a []f64, lda int, work []f64) f64
fn dorgqr #
fn dorgqr(m int, n int, k int, mut a []f64, lda int, tau []f64) int
dorgqr - Direct wrapper for dorgqr with standardized interface
fn dorgqr_standardized #
fn dorgqr_standardized(m int, n int, k int, mut a []f64, lda int, tau []f64) int
dorgqr_standardized - Standardized dorgqr wrapper
fn dpotrf #
fn dpotrf(uplo blas.Uplo, n int, mut a []f64, lda int) int
dpotrf computes the Cholesky factorization of a real symmetric positive definite matrix A.
See: http://www.netlib.org/lapack/explore-html/d0/d8a/dpotrf_8f.html
See: https://software.intel.com/en-us/mkl-developer-reference-c-potrf
The factorization has the form
A = U**T * U, if UPLO = 'U'
or
A = L * L**T, if UPLO = 'L'
where U is an upper triangular matrix and L is lower triangular.
This is the block version of the algorithm, calling Level 3 BLAS.
fn dpotrf_standardized #
fn dpotrf_standardized(uplo blas.Uplo, n int, mut a []f64, lda int) int
dpotrf_standardized - Standardized dpotrf wrapper
fn dsyev #
fn dsyev(jobz EigenVectorsJob, uplo blas.Uplo, n int, mut a []f64, lda int, mut w []f64) int
dsyev - Direct wrapper for dsyev with standardized interface
fn dsyev_standardized #
fn dsyev_standardized(jobz EigenVectorsJob, uplo blas.Uplo, n int, mut a []f64, lda int, mut w []f64) int
dsyev_standardized - Standardized dsyev wrapper
fn from_lapack64_apply_ortho #
fn from_lapack64_apply_ortho(apply lapack64.ApplyOrtho) ApplyOrtho
Convert lapack64.ApplyOrtho to VSL ApplyOrtho
fn from_lapack64_balance_job #
fn from_lapack64_balance_job(job lapack64.BalanceJob) BalanceJob
Convert lapack64.BalanceJob to VSL BalanceJob
fn from_lapack64_direct #
fn from_lapack64_direct(direct lapack64.Direct) Direct
Convert lapack64.Direct to VSL Direct
fn from_lapack64_eigen_vectors_comp #
fn from_lapack64_eigen_vectors_comp(comp lapack64.EigenVectorsComp) EigenVectorsComp
Convert lapack64.EigenVectorsComp to VSL EigenVectorsComp
fn from_lapack64_eigen_vectors_how_many #
fn from_lapack64_eigen_vectors_how_many(how_many lapack64.EigenVectorsHowMany) EigenVectorsHowMany
Convert lapack64.EigenVectorsHowMany to VSL EigenVectorsHowMany
fn from_lapack64_eigen_vectors_job #
fn from_lapack64_eigen_vectors_job(job lapack64.EigenVectorsJob) EigenVectorsJob
Convert lapack64.EigenVectorsJob to VSL EigenVectorsJob
fn from_lapack64_eigen_vectors_side #
fn from_lapack64_eigen_vectors_side(side lapack64.EigenVectorsSide) EigenVectorsSide
Convert lapack64.EigenVectorsSide to VSL EigenVectorsSide
fn from_lapack64_gen_ortho #
fn from_lapack64_gen_ortho(gen lapack64.GenOrtho) GenOrtho
Convert lapack64.GenOrtho to VSL GenOrtho
fn from_lapack64_gsvd_job #
fn from_lapack64_gsvd_job(job lapack64.GSVDJob) GSVDJob
Convert lapack64.GSVDJob to VSL GSVDJob
fn from_lapack64_left_eigen_vectors_job #
fn from_lapack64_left_eigen_vectors_job(job lapack64.LeftEigenVectorsJob) LeftEigenVectorsJob
Convert lapack64.LeftEigenVectorsJob to VSL LeftEigenVectorsJob
fn from_lapack64_matrix_norm #
fn from_lapack64_matrix_norm(norm lapack64.MatrixNorm) MatrixNorm
Convert lapack64.MatrixNorm to VSL MatrixNorm
fn from_lapack64_matrix_type #
fn from_lapack64_matrix_type(mtype lapack64.MatrixType) MatrixType
Convert lapack64.MatrixType to VSL MatrixType
fn from_lapack64_maximize_norm_x_job #
fn from_lapack64_maximize_norm_x_job(job lapack64.MaximizeNormXJob) MaximizeNormXJob
Convert lapack64.MaximizeNormXJob to VSL MaximizeNormXJob
fn from_lapack64_ortho_comp #
fn from_lapack64_ortho_comp(comp lapack64.OrthoComp) OrthoComp
Convert lapack64.OrthoComp to VSL OrthoComp
fn from_lapack64_pivot #
fn from_lapack64_pivot(pivot lapack64.Pivot) Pivot
Convert lapack64.Pivot to VSL Pivot
fn from_lapack64_right_eigen_vectors_job #
fn from_lapack64_right_eigen_vectors_job(job lapack64.RightEigenVectorsJob) RightEigenVectorsJob
Convert lapack64.RightEigenVectorsJob to VSL RightEigenVectorsJob
fn from_lapack64_schur_comp #
fn from_lapack64_schur_comp(comp lapack64.SchurComp) SchurComp
Convert lapack64.SchurComp to VSL SchurComp
fn from_lapack64_schur_job #
fn from_lapack64_schur_job(job lapack64.SchurJob) SchurJob
Convert lapack64.SchurJob to VSL SchurJob
fn from_lapack64_sort #
fn from_lapack64_sort(sort lapack64.Sort) Sort
Convert lapack64.Sort to VSL Sort
fn from_lapack64_storev #
fn from_lapack64_storev(storev lapack64.StoreV) StoreV
Convert lapack64.StoreV to VSL StoreV
fn from_lapack64_svd_job #
fn from_lapack64_svd_job(job lapack64.SVDJob) SVDJob
Convert lapack64.SVDJob to VSL SVDJob
fn from_lapack64_update_schur_comp #
fn from_lapack64_update_schur_comp(comp lapack64.UpdateSchurComp) UpdateSchurComp
Convert lapack64.UpdateSchurComp to VSL UpdateSchurComp
fn geev #
fn geev(a [][]f64, jobvl LeftEigenVectorsJob, jobvr RightEigenVectorsJob) !([]f64, []f64, [][]f64, [][]f64)
geev computes the eigenvalues and, optionally, the left and/or right eigenvectors for a real nonsymmetric matrix A.
fn geqrf #
fn geqrf(mut a [][]f64) ![]f64
geqrf computes a QR factorization of a real m×n matrix A.
fn gesv #
fn gesv(mut a [][]f64, mut b [][]f64) !
gesv solves the linear system A*X = B using LU factorization with partial pivoting. A is an n×n matrix and B is an n×nrhs matrix. On return, A is overwritten with the LU factorization and B contains the solution X.
fn gesvd #
fn gesvd(a [][]f64, jobu SVDJob, jobvt SVDJob) !([]f64, [][]f64, [][]f64)
gesvd computes the singular value decomposition (SVD) of a real m×n matrix A.
fn getrf #
fn getrf(mut a [][]f64) ![]int
getrf computes the LU factorization of a general m×n matrix A using partial pivoting. On return, A contains the LU factorization and the function returns the pivot indices.
fn getri #
fn getri(mut a [][]f64, mut ipiv []int) !
getri computes the inverse of a matrix using the LU factorization computed by getrf. A should contain the LU factorization from getrf, and ipiv should be the pivot indices.
fn orgqr #
fn orgqr(mut a [][]f64, tau []f64) !
orgqr generates the m×n matrix Q with orthonormal columns defined as the first n columns of a product of k elementary reflectors of order m.
fn potrf #
fn potrf(mut a [][]f64, uplo blas.Uplo) !
potrf computes the Cholesky factorization of a symmetric positive definite matrix. uplo specifies whether the upper or lower triangle of A is stored.
fn syev #
fn syev(mut a [][]f64, jobz EigenVectorsJob, uplo blas.Uplo) ![]f64
syev computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.
fn to_lapack64_apply_ortho #
fn to_lapack64_apply_ortho(apply ApplyOrtho) lapack64.ApplyOrtho
Convert VSL ApplyOrtho to lapack64.ApplyOrtho
fn to_lapack64_balance_job #
fn to_lapack64_balance_job(job BalanceJob) lapack64.BalanceJob
Convert VSL BalanceJob to lapack64.BalanceJob
fn to_lapack64_direct #
fn to_lapack64_direct(direct Direct) lapack64.Direct
Convert VSL Direct to lapack64.Direct
fn to_lapack64_eigen_vectors_comp #
fn to_lapack64_eigen_vectors_comp(comp EigenVectorsComp) lapack64.EigenVectorsComp
Convert VSL EigenVectorsComp to lapack64.EigenVectorsComp
fn to_lapack64_eigen_vectors_how_many #
fn to_lapack64_eigen_vectors_how_many(how_many EigenVectorsHowMany) lapack64.EigenVectorsHowMany
Convert VSL EigenVectorsHowMany to lapack64.EigenVectorsHowMany
fn to_lapack64_eigen_vectors_job #
fn to_lapack64_eigen_vectors_job(job EigenVectorsJob) lapack64.EigenVectorsJob
Convert VSL EigenVectorsJob to lapack64.EigenVectorsJob
fn to_lapack64_eigen_vectors_side #
fn to_lapack64_eigen_vectors_side(side EigenVectorsSide) lapack64.EigenVectorsSide
Convert VSL EigenVectorsSide to lapack64.EigenVectorsSide
fn to_lapack64_gen_ortho #
fn to_lapack64_gen_ortho(gen GenOrtho) lapack64.GenOrtho
Convert VSL GenOrtho to lapack64.GenOrtho
fn to_lapack64_gsvd_job #
fn to_lapack64_gsvd_job(job GSVDJob) lapack64.GSVDJob
Convert VSL GSVDJob to lapack64.GSVDJob
fn to_lapack64_left_eigen_vectors_job #
fn to_lapack64_left_eigen_vectors_job(job LeftEigenVectorsJob) lapack64.LeftEigenVectorsJob
Convert VSL LeftEigenVectorsJob to lapack64.LeftEigenVectorsJob
fn to_lapack64_matrix_norm #
fn to_lapack64_matrix_norm(norm MatrixNorm) lapack64.MatrixNorm
Convert VSL MatrixNorm to lapack64.MatrixNorm
fn to_lapack64_matrix_type #
fn to_lapack64_matrix_type(mtype MatrixType) lapack64.MatrixType
Convert VSL MatrixType to lapack64.MatrixType
fn to_lapack64_maximize_norm_x_job #
fn to_lapack64_maximize_norm_x_job(job MaximizeNormXJob) lapack64.MaximizeNormXJob
Convert VSL MaximizeNormXJob to lapack64.MaximizeNormXJob
fn to_lapack64_ortho_comp #
fn to_lapack64_ortho_comp(comp OrthoComp) lapack64.OrthoComp
Convert VSL OrthoComp to lapack64.OrthoComp
fn to_lapack64_pivot #
fn to_lapack64_pivot(pivot Pivot) lapack64.Pivot
Convert VSL Pivot to lapack64.Pivot
fn to_lapack64_right_eigen_vectors_job #
fn to_lapack64_right_eigen_vectors_job(job RightEigenVectorsJob) lapack64.RightEigenVectorsJob
Convert VSL RightEigenVectorsJob to lapack64.RightEigenVectorsJob
fn to_lapack64_schur_comp #
fn to_lapack64_schur_comp(comp SchurComp) lapack64.SchurComp
Convert VSL SchurComp to lapack64.SchurComp
fn to_lapack64_schur_job #
fn to_lapack64_schur_job(job SchurJob) lapack64.SchurJob
Convert VSL SchurJob to lapack64.SchurJob
fn to_lapack64_sort #
fn to_lapack64_sort(sort Sort) lapack64.Sort
Convert VSL Sort to lapack64.Sort
fn to_lapack64_storev #
fn to_lapack64_storev(storev StoreV) lapack64.StoreV
Convert VSL StoreV to lapack64.StoreV
fn to_lapack64_svd_job #
fn to_lapack64_svd_job(job SVDJob) lapack64.SVDJob
Convert VSL SVDJob to lapack64.SVDJob
fn to_lapack64_update_schur_comp #
fn to_lapack64_update_schur_comp(comp UpdateSchurComp) lapack64.UpdateSchurComp
Convert VSL UpdateSchurComp to lapack64.UpdateSchurComp
enum ApplyOrtho #
enum ApplyOrtho as u8 {
// Apply P or Pᵀ.
apply_p = u8(`P`)
// Apply Q or Qᵀ.
apply_q = u8(`Q`)
}
ApplyOrtho specifies which orthogonal matrix is applied in Dormbr.
enum BalanceJob #
enum BalanceJob as u8 {
permute = u8(`P`)
scale = u8(`S`)
permute_scale = u8(`B`)
balance_none = u8(`N`)
}
BalanceJob specifies matrix balancing operation.
enum Direct #
enum Direct as u8 {
// Reflectors are right-multiplied, H_0 * H_1 * ... * H_{k-1}.
forward = u8(`F`)
// Reflectors are left-multiplied, H_{k-1} * ... * H_1 * H_0.
backward = u8(`B`)
}
Direct specifies the direction of the multiplication for the Householder matrix.
enum EigenVectorsComp #
enum EigenVectorsComp as u8 {
// Compute eigenvectors of the original symmetric matrix.
ev_orig = u8(`V`)
// Compute eigenvectors of the tridiagonal matrix.
ev_tridiag = u8(`I`)
// Do not compute eigenvectors.
ev_comp_none = u8(`N`)
}
EigenVectorsComp specifies how eigenvectors are computed in Dsteqr.
enum EigenVectorsHowMany #
enum EigenVectorsHowMany as u8 {
// Compute all right and/or left eigenvectors.
ev_all = u8(`A`)
// Compute all right and/or left eigenvectors multiplied by an input matrix.
ev_all_mul_q = u8(`B`)
// Compute selected right and/or left eigenvectors.
ev_selected = u8(`S`)
}
EigenVectorsHowMany specifies which eigenvectors are computed in Dtrevc3 and how.
enum EigenVectorsJob #
enum EigenVectorsJob as u8 {
// Compute eigenvectors.
ev_compute = u8(`V`)
// Do not compute eigenvectors.
ev_none = u8(`N`)
}
EigenVectorsJob specifies whether eigenvectors are computed in Dsyev.
enum EigenVectorsSide #
enum EigenVectorsSide as u8 {
// Compute only right eigenvectors.
ev_right = u8(`R`)
// Compute only left eigenvectors.
ev_left = u8(`L`)
// Compute both right and left eigenvectors.
ev_both = u8(`B`)
}
EigenVectorsSide specifies what eigenvectors are computed in Dtrevc3.
enum GSVDJob #
enum GSVDJob as u8 {
// Compute orthogonal matrix U.
gsvd_u = u8(`U`)
// Compute orthogonal matrix V.
gsvd_v = u8(`V`)
// Compute orthogonal matrix Q.
gsvd_q = u8(`Q`)
// Use unit-initialized matrix.
gsvd_unit = u8(`I`)
// Do not compute orthogonal matrix.
gsvd_none = u8(`N`)
}
GSVDJob specifies the singular vector computation type for Generalized SingularValueDecomposition.
enum GenOrtho #
enum GenOrtho as u8 {
// Generate Pᵀ.
generate_pt = u8(`P`)
// Generate Q.
generate_q = u8(`Q`)
}
GenOrtho specifies which orthogonal matrix is generated in Dorgbr.
enum LeftEigenVectorsJob #
enum LeftEigenVectorsJob as u8 {
// Compute left eigenvectors.
left_ev_compute = u8(`V`)
// Do not compute left eigenvectors.
left_ev_none = u8(`N`)
}
LeftEigenVectorsJob specifies whether left eigenvectors are computed in Dgeev.
enum MatrixNorm #
enum MatrixNorm as u8 {
// max(abs(A(i,j)))
max_abs = u8(`M`)
// Maximum absolute column sum (one norm)
max_column_sum = u8(`O`)
// Maximum absolute row sum (infinity norm)
max_row_sum = u8(`I`)
// Frobenius norm (sqrt of sum of squares)
frobenius = u8(`F`)
}
MatrixNorm represents the kind of matrix norm to compute.
enum MatrixType #
enum MatrixType as u8 {
// A general dense matrix.
general = u8(`G`)
// An upper triangular matrix.
upper_tri = u8(`U`)
// A lower triangular matrix.
lower_tri = u8(`L`)
}
MatrixType represents the kind of matrix represented in the data.
enum MaximizeNormXJob #
enum MaximizeNormXJob as u8 {
// Solve Z*x=h-f where h is a vector of ±1.
local_look_ahead = 0
// Compute an approximate null-vector e of Z, normalize e and solve Z*x=±e-f.
normalized_null_vector = 2
}
MaximizeNormXJob specifies the heuristic method for computing a contribution to the reciprocal Dif-estimate in Dlatdf.
enum OrthoComp #
enum OrthoComp as u8 {
// Do not compute the orthogonal matrix.
ortho_none = u8(`N`)
// The orthogonal matrix is formed explicitly and returned in the argument.
ortho_explicit = u8(`I`)
// The orthogonal matrix is post-multiplied into the matrix stored in the argument on entry.
ortho_postmul = u8(`V`)
}
OrthoComp specifies whether and how the orthogonal matrix is computed in Dgghrd.
enum Pivot #
enum Pivot as u8 {
variable = u8(`V`)
top = u8(`T`)
bottom = u8(`B`)
}
Pivot specifies the pivot type for plane rotations.
enum RightEigenVectorsJob #
enum RightEigenVectorsJob as u8 {
// Compute right eigenvectors.
right_ev_compute = u8(`V`)
// Do not compute right eigenvectors.
right_ev_none = u8(`N`)
}
RightEigenVectorsJob specifies whether right eigenvectors are computed in Dgeev.
enum SVDJob #
enum SVDJob as u8 {
// Compute all columns of the orthogonal matrix U or V.
svd_all = u8(`A`)
// Compute the singular vectors and store them in the orthogonal matrix U or V.
svd_store = u8(`S`)
// Compute the singular vectors and overwrite them on the input matrix A.
svd_overwrite = u8(`O`)
// Do not compute singular vectors.
svd_none = u8(`N`)
}
SVDJob specifies the singular vector computation type for SingularValueDecomposition.
enum SchurComp #
enum SchurComp as u8 {
// Compute Schur vectors of the original matrix.
schur_orig = u8(`V`)
// Compute Schur vectors of the upper Hessenberg matrix.
schur_hess = u8(`I`)
// Do not compute Schur vectors.
schur_none = u8(`N`)
}
SchurComp specifies whether and how the Schur vectors are computed in Dhseqr.
enum SchurJob #
enum SchurJob as u8 {
eigenvalues_only = u8(`E`)
eigenvalues_and_schur = u8(`S`)
}
SchurJob specifies whether the Schur form is computed in Dhseqr.
enum Sort #
enum Sort as u8 {
sort_increasing = u8(`I`)
sort_decreasing = u8(`D`)
}
Sort is the sorting order.
enum StoreV #
enum StoreV as u8 {
// Reflector stored in a column of the matrix.
column_wise = u8(`C`)
// Reflector stored in a row of the matrix.
row_wise = u8(`R`)
}
StoreV indicates the storage direction of elementary reflectors.
enum UpdateSchurComp #
enum UpdateSchurComp as u8 {
// Update the matrix of Schur vectors.
update_schur = u8(`V`)
// Do not update the matrix of Schur vectors.
update_schur_none = u8(`N`)
}
UpdateSchurComp specifies whether the matrix of Schur vectors is updated in Dtrexc.
- README
- fn dgeev
- fn dgeev_standardized
- fn dgeqrf
- fn dgeqrf_standardized
- fn dgesv
- fn dgesvd
- fn dgetrf
- fn dgetri
- fn dlange
- fn dorgqr
- fn dorgqr_standardized
- fn dpotrf
- fn dpotrf_standardized
- fn dsyev
- fn dsyev_standardized
- fn from_lapack64_apply_ortho
- fn from_lapack64_balance_job
- fn from_lapack64_direct
- fn from_lapack64_eigen_vectors_comp
- fn from_lapack64_eigen_vectors_how_many
- fn from_lapack64_eigen_vectors_job
- fn from_lapack64_eigen_vectors_side
- fn from_lapack64_gen_ortho
- fn from_lapack64_gsvd_job
- fn from_lapack64_left_eigen_vectors_job
- fn from_lapack64_matrix_norm
- fn from_lapack64_matrix_type
- fn from_lapack64_maximize_norm_x_job
- fn from_lapack64_ortho_comp
- fn from_lapack64_pivot
- fn from_lapack64_right_eigen_vectors_job
- fn from_lapack64_schur_comp
- fn from_lapack64_schur_job
- fn from_lapack64_sort
- fn from_lapack64_storev
- fn from_lapack64_svd_job
- fn from_lapack64_update_schur_comp
- fn geev
- fn geqrf
- fn gesv
- fn gesvd
- fn getrf
- fn getri
- fn orgqr
- fn potrf
- fn syev
- fn to_lapack64_apply_ortho
- fn to_lapack64_balance_job
- fn to_lapack64_direct
- fn to_lapack64_eigen_vectors_comp
- fn to_lapack64_eigen_vectors_how_many
- fn to_lapack64_eigen_vectors_job
- fn to_lapack64_eigen_vectors_side
- fn to_lapack64_gen_ortho
- fn to_lapack64_gsvd_job
- fn to_lapack64_left_eigen_vectors_job
- fn to_lapack64_matrix_norm
- fn to_lapack64_matrix_type
- fn to_lapack64_maximize_norm_x_job
- fn to_lapack64_ortho_comp
- fn to_lapack64_pivot
- fn to_lapack64_right_eigen_vectors_job
- fn to_lapack64_schur_comp
- fn to_lapack64_schur_job
- fn to_lapack64_sort
- fn to_lapack64_storev
- fn to_lapack64_svd_job
- fn to_lapack64_update_schur_comp
- enum ApplyOrtho
- enum BalanceJob
- enum Direct
- enum EigenVectorsComp
- enum EigenVectorsHowMany
- enum EigenVectorsJob
- enum EigenVectorsSide
- enum GSVDJob
- enum GenOrtho
- enum LeftEigenVectorsJob
- enum MatrixNorm
- enum MatrixType
- enum MaximizeNormXJob
- enum OrthoComp
- enum Pivot
- enum RightEigenVectorsJob
- enum SVDJob
- enum SchurComp
- enum SchurJob
- enum Sort
- enum StoreV
- enum UpdateSchurComp