Teuchos - Trilinos Tools Package
Version of the Day
|
The Templated LAPACK Wrapper Class. More...
#include <Teuchos_LAPACK.hpp>
Public Member Functions | |
Constructors/Destructors. | |
LAPACK (void) | |
Default Constructor. More... | |
LAPACK (const LAPACK< OrdinalType, ScalarType > &lapack) | |
Copy Constructor. More... | |
virtual | ~LAPACK (void) |
Destructor. More... | |
Symmetric Positive Definite Linear System Routines. | |
void | PTTRF (const OrdinalType &n, ScalarType *d, ScalarType *e, OrdinalType *info) const |
Computes the L*D*L' factorization of a Hermitian/symmetric positive definite tridiagonal matrix A . More... | |
void | PTTRS (const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *d, const ScalarType *e, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const |
Solves a tridiagonal system A*X=B using the *D*L' factorization of A computed by PTTRF. More... | |
void | POTRF (const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *info) const |
Computes Cholesky factorization of a real symmetric positive definite matrix A . More... | |
void | POTRS (const char &UPLO, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const |
Solves a system of linear equations A*X=B , where A is a symmetric positive definite matrix factored by POTRF and the nrhs solutions are returned in B . More... | |
void | POTRI (const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *info) const |
Computes the inverse of a real symmetric positive definite matrix A using the Cholesky factorization A from POTRF. More... | |
void | POCON (const char &UPLO, const OrdinalType &n, const ScalarType *A, const OrdinalType &lda, const ScalarType &anorm, ScalarType *rcond, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const |
Estimates the reciprocal of the condition number (1-norm) of a real symmetric positive definite matrix A using the Cholesky factorization from POTRF. More... | |
void | POSV (const char &UPLO, const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const |
Computes the solution to a real system of linear equations A*X=B , where A is a symmetric positive definite matrix and the nrhs solutions are returned in B . More... | |
void | POEQU (const OrdinalType &n, const ScalarType *A, const OrdinalType &lda, MagnitudeType *S, MagnitudeType *scond, MagnitudeType *amax, OrdinalType *info) const |
Computes row and column scalings intended to equilibrate a symmetric positive definite matrix A and reduce its condition number (w.r.t. 2-norm). More... | |
void | PORFS (const char &UPLO, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const ScalarType *AF, const OrdinalType &ldaf, const ScalarType *B, const OrdinalType &ldb, ScalarType *X, const OrdinalType &ldx, ScalarType *FERR, ScalarType *BERR, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const |
Improves the computed solution to a system of linear equations when the coefficient matrix is symmetric positive definite, and provides error bounds and backward error estimates for the solution. More... | |
void | POSVX (const char &FACT, const char &UPLO, const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, ScalarType *AF, const OrdinalType &ldaf, char *EQUED, ScalarType *S, ScalarType *B, const OrdinalType &ldb, ScalarType *X, const OrdinalType &ldx, ScalarType *rcond, ScalarType *FERR, ScalarType *BERR, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const |
Uses the Cholesky factorization to compute the solution to a real system of linear equations A*X=B , where A is symmetric positive definite. System can be equilibrated by POEQU and iteratively refined by PORFS, if requested. More... | |
General Linear System Routines. | |
void | GELS (const char &TRANS, const OrdinalType &m, const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Solves an over/underdetermined real m by n linear system A using QR or LQ factorization of A. More... | |
void | GELSS (const OrdinalType &m, const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, MagnitudeType *S, const MagnitudeType rcond, OrdinalType *rank, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const |
Use the SVD to solve a possibly rank-deficient linear least-squares problem. More... | |
void | GELSS (const OrdinalType &m, const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, ScalarType *S, const ScalarType &rcond, OrdinalType *rank, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Legacy GELSS interface for real-valued ScalarType. More... | |
void | GGLSE (const OrdinalType &m, const OrdinalType &n, const OrdinalType &p, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, ScalarType *C, ScalarType *D, ScalarType *X, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Solves the linear equality-constrained least squares (LSE) problem where A is an m by n matrix,B is a p by n matrix C is a given m-vector , and D is a given p-vector . More... | |
void | GEQRF (const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Computes a QR factorization of a general m by n matrix A . More... | |
void | GEQR2 (const OrdinalType &m, const OrdinalType &n, ScalarType A[], const OrdinalType &lda, ScalarType TAU[], ScalarType WORK[], OrdinalType *const info) const |
BLAS 2 version of GEQRF, with known workspace size. More... | |
void | GETRF (const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *IPIV, OrdinalType *info) const |
Computes an LU factorization of a general m by n matrix A using partial pivoting with row interchanges. More... | |
void | GETRS (const char &TRANS, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const OrdinalType *IPIV, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const |
Solves a system of linear equations A*X=B or A'*X=B with a general n by n matrix A using the LU factorization computed by GETRF. More... | |
void | LASCL (const char &TYPE, const OrdinalType &kl, const OrdinalType &ku, const MagnitudeType cfrom, const MagnitudeType cto, const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *info) const |
Multiplies the m by n matrix A by the real scalar cto/cfrom . More... | |
void | GEQP3 (const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *jpvt, ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const |
Computes a QR factorization with column pivoting of a matrix A: A*P = Q*R using Level 3 BLAS. More... | |
void | LASWP (const OrdinalType &N, ScalarType A[], const OrdinalType &LDA, const OrdinalType &K1, const OrdinalType &K2, const OrdinalType IPIV[], const OrdinalType &INCX) const |
Apply a series of row interchanges to the matrix A. More... | |
void | GBTRF (const OrdinalType &m, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, ScalarType *A, const OrdinalType &lda, OrdinalType *IPIV, OrdinalType *info) const |
Computes an LU factorization of a general banded m by n matrix A using partial pivoting with row interchanges. More... | |
void | GBTRS (const char &TRANS, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const OrdinalType *IPIV, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const |
Solves a system of linear equations A*X=B or A'*X=B with a general banded n by n matrix A using the LU factorization computed by GBTRF. More... | |
void | GTTRF (const OrdinalType &n, ScalarType *dl, ScalarType *d, ScalarType *du, ScalarType *du2, OrdinalType *IPIV, OrdinalType *info) const |
Computes an LU factorization of a n by n tridiagonal matrix A using partial pivoting with row interchanges. More... | |
void | GTTRS (const char &TRANS, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *dl, const ScalarType *d, const ScalarType *du, const ScalarType *du2, const OrdinalType *IPIV, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const |
Solves a system of linear equations A*X=B or A'*X=B or A^H*X=B with a tridiagonal matrix A using the LU factorization computed by GTTRF. More... | |
void | GETRI (const OrdinalType &n, ScalarType *A, const OrdinalType &lda, const OrdinalType *IPIV, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Computes the inverse of a matrix A using the LU factorization computed by GETRF. More... | |
void | LATRS (const char &UPLO, const char &TRANS, const char &DIAG, const char &NORMIN, const OrdinalType &N, ScalarType *A, const OrdinalType &LDA, ScalarType *X, MagnitudeType *SCALE, MagnitudeType *CNORM, OrdinalType *INFO) const |
Robustly solve a possibly singular triangular linear system. More... | |
void | GECON (const char &NORM, const OrdinalType &n, const ScalarType *A, const OrdinalType &lda, const ScalarType &anorm, ScalarType *rcond, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const |
Estimates the reciprocal of the condition number of a general real matrix A , in either the 1-norm or the infinity-norm, using the LU factorization computed by GETRF. More... | |
void | GBCON (const char &NORM, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const ScalarType *A, const OrdinalType &lda, OrdinalType *IPIV, const ScalarType &anorm, ScalarType *rcond, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const |
Estimates the reciprocal of the condition number of a general banded real matrix A , in either the 1-norm or the infinity-norm, using the LU factorization computed by GETRF. More... | |
ScalarTraits< ScalarType > ::magnitudeType | LANGB (const char &NORM, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const ScalarType *A, const OrdinalType &lda, MagnitudeType *WORK) const |
Returns the value of the one norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of an n by n band matrix A , with kl sub-diagonals and ku super-diagonals. More... | |
void | GESV (const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, OrdinalType *IPIV, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const |
Computes the solution to a real system of linear equations A*X=B , where A is factored through GETRF and the nrhs solutions are computed through GETRS. More... | |
void | GEEQU (const OrdinalType &m, const OrdinalType &n, const ScalarType *A, const OrdinalType &lda, ScalarType *R, ScalarType *C, ScalarType *rowcond, ScalarType *colcond, ScalarType *amax, OrdinalType *info) const |
Computes row and column scalings intended to equilibrate an m by n matrix A and reduce its condition number. More... | |
void | GERFS (const char &TRANS, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const ScalarType *AF, const OrdinalType &ldaf, const OrdinalType *IPIV, const ScalarType *B, const OrdinalType &ldb, ScalarType *X, const OrdinalType &ldx, ScalarType *FERR, ScalarType *BERR, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const |
Improves the computed solution to a system of linear equations and provides error bounds and backward error estimates for the solution. Use after GETRF/GETRS. More... | |
void | GBEQU (const OrdinalType &m, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const ScalarType *A, const OrdinalType &lda, MagnitudeType *R, MagnitudeType *C, MagnitudeType *rowcond, MagnitudeType *colcond, MagnitudeType *amax, OrdinalType *info) const |
Computes row and column scalings intended to equilibrate an m by n banded matrix A and reduce its condition number. More... | |
void | GBRFS (const char &TRANS, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const ScalarType *AF, const OrdinalType &ldaf, const OrdinalType *IPIV, const ScalarType *B, const OrdinalType &ldb, ScalarType *X, const OrdinalType &ldx, ScalarType *FERR, ScalarType *BERR, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const |
Improves the computed solution to a banded system of linear equations and provides error bounds and backward error estimates for the solution. Use after GBTRF/GBTRS. More... | |
void | GESVX (const char &FACT, const char &TRANS, const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, ScalarType *AF, const OrdinalType &ldaf, OrdinalType *IPIV, char *EQUED, ScalarType *R, ScalarType *C, ScalarType *B, const OrdinalType &ldb, ScalarType *X, const OrdinalType &ldx, ScalarType *rcond, ScalarType *FERR, ScalarType *BERR, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const |
Uses the LU factorization to compute the solution to a real system of linear equations A*X=B , returning error bounds on the solution and a condition estimate. More... | |
void | SYTRD (const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *D, ScalarType *E, ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Reduces a real symmetric matrix A to tridiagonal form by orthogonal similarity transformations. More... | |
void | GEHRD (const OrdinalType &n, const OrdinalType &ilo, const OrdinalType &ihi, ScalarType *A, const OrdinalType &lda, ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Reduces a real general matrix A to upper Hessenberg form by orthogonal similarity transformations. More... | |
void | TRTRS (const char &UPLO, const char &TRANS, const char &DIAG, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const |
Solves a triangular linear system of the form A*X=B or A**T*X=B , where A is a triangular matrix. More... | |
void | TRTRI (const char &UPLO, const char &DIAG, const OrdinalType &n, const ScalarType *A, const OrdinalType &lda, OrdinalType *info) const |
Computes the inverse of an upper or lower triangular matrix A . More... | |
Symmetric Eigenproblem Routines | |
void | SPEV (const char &JOBZ, const char &UPLO, const OrdinalType &n, ScalarType *AP, ScalarType *W, ScalarType *Z, const OrdinalType &ldz, ScalarType *WORK, OrdinalType *info) const |
Computes the eigenvalues and, optionally, eigenvectors of a symmetric n by n matrix A in packed storage. More... | |
void | SYEV (const char &JOBZ, const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *W, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Computes all the eigenvalues and, optionally, eigenvectors of a symmetric n by n matrix A. More... | |
void | SYGV (const OrdinalType &itype, const char &JOBZ, const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, ScalarType *W, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Computes all the eigenvalues and, optionally, eigenvectors of a symmetric n by n matrix pencil {A ,B}, where A is symmetric and B is symmetric positive-definite. More... | |
void | HEEV (const char &JOBZ, const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, MagnitudeType *W, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const |
Computes all the eigenvalues and, optionally, eigenvectors of a Hermitian n by n matrix A. More... | |
void | HEGV (const OrdinalType &itype, const char &JOBZ, const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, MagnitudeType *W, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const |
Computes all the eigenvalues and, optionally, eigenvectors of a generalized Hermitian-definite n by n matrix pencil {A ,B}, where A is Hermitian and B is Hermitian positive-definite. More... | |
void | STEQR (const char &COMPZ, const OrdinalType &n, MagnitudeType *D, MagnitudeType *E, ScalarType *Z, const OrdinalType &ldz, MagnitudeType *WORK, OrdinalType *info) const |
Computes the eigenvalues and, optionally, eigenvectors of a symmetric tridiagonal n by n matrix A using implicit QL/QR. The eigenvectors can only be computed if A was reduced to tridiagonal form by SYTRD. More... | |
void | PTEQR (const char &COMPZ, const OrdinalType &n, MagnitudeType *D, MagnitudeType *E, ScalarType *Z, const OrdinalType &ldz, MagnitudeType *WORK, OrdinalType *info) const |
Computes the eigenvalues and, optionally, eigenvectors of a symmetric positive-definite tridiagonal n by n matrix A using BDSQR, after factoring the matrix with PTTRF. More... | |
Non-Hermitian Eigenproblem Routines | |
void | HSEQR (const char &JOB, const char &COMPZ, const OrdinalType &n, const OrdinalType &ilo, const OrdinalType &ihi, ScalarType *H, const OrdinalType &ldh, ScalarType *WR, ScalarType *WI, ScalarType *Z, const OrdinalType &ldz, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Computes the eigenvalues of a real upper Hessenberg matrix H and, optionally, the matrices T and Z from the Schur decomposition, where T is an upper quasi-triangular matrix and Z contains the Schur vectors. More... | |
void | GEES (const char &JOBVS, const char &SORT, OrdinalType &(*ptr2func)(ScalarType *, ScalarType *), const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *sdim, ScalarType *WR, ScalarType *WI, ScalarType *VS, const OrdinalType &ldvs, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *BWORK, OrdinalType *info) const |
void | GEES (const char &JOBVS, const char &SORT, OrdinalType &(*ptr2func)(ScalarType *), const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *sdim, ScalarType *W, ScalarType *VS, const OrdinalType &ldvs, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *BWORK, OrdinalType *info) const |
void | GEES (const char &JOBVS, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *sdim, MagnitudeType *WR, MagnitudeType *WI, ScalarType *VS, const OrdinalType &ldvs, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *BWORK, OrdinalType *info) const |
void | GEEV (const char &JOBVL, const char &JOBVR, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, MagnitudeType *WR, MagnitudeType *WI, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const |
Computes for an n by n real nonsymmetric matrix A , the eigenvalues and, optionally, the left and/or right eigenvectors. More... | |
void | GEEVX (const char &BALANC, const char &JOBVL, const char &JOBVR, const char &SENSE, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *WR, ScalarType *WI, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, OrdinalType *ilo, OrdinalType *ihi, MagnitudeType *SCALE, MagnitudeType *abnrm, MagnitudeType *RCONDE, MagnitudeType *RCONDV, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *IWORK, OrdinalType *info) const |
void | GGEVX (const char &BALANC, const char &JOBVL, const char &JOBVR, const char &SENSE, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, ScalarType *BETA, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, OrdinalType *ilo, OrdinalType *ihi, MagnitudeType *lscale, MagnitudeType *rscale, MagnitudeType *abnrm, MagnitudeType *bbnrm, MagnitudeType *RCONDE, MagnitudeType *RCONDV, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *IWORK, OrdinalType *BWORK, OrdinalType *info) const |
void | GGEV (const char &JOBVL, const char &JOBVR, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, ScalarType *BETA, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
void | TRSEN (const char &JOB, const char &COMPQ, const OrdinalType *SELECT, const OrdinalType &n, ScalarType *T, const OrdinalType &ldt, ScalarType *Q, const OrdinalType &ldq, MagnitudeType *WR, MagnitudeType *WI, OrdinalType *M, ScalarType *S, MagnitudeType *SEP, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *IWORK, const OrdinalType &liwork, OrdinalType *info) const |
void | TGSEN (const OrdinalType &ijob, const OrdinalType &wantq, const OrdinalType &wantz, const OrdinalType *SELECT, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, MagnitudeType *BETA, ScalarType *Q, const OrdinalType &ldq, ScalarType *Z, const OrdinalType &ldz, OrdinalType *M, MagnitudeType *PL, MagnitudeType *PR, MagnitudeType *DIF, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *IWORK, const OrdinalType &liwork, OrdinalType *info) const |
void | GGES (const char &JOBVL, const char &JOBVR, const char &SORT, OrdinalType &(*ptr2func)(ScalarType *, ScalarType *, ScalarType *), const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, OrdinalType *sdim, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, MagnitudeType *BETA, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *BWORK, OrdinalType *info) const |
Singular Value Decompositon Routines | |
void | GESVD (const char &JOBU, const char &JOBVT, const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, MagnitudeType *S, ScalarType *U, const OrdinalType &ldu, ScalarType *V, const OrdinalType &ldv, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const |
Computes the singular values (and optionally, vectors) of a real matrix A . More... | |
Orthogonal matrix routines | |
void | ORMQR (const char &SIDE, const char &TRANS, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *C, const OrdinalType &ldc, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
void | ORM2R (const char &SIDE, const char &TRANS, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, const ScalarType A[], const OrdinalType &lda, const ScalarType TAU[], ScalarType C[], const OrdinalType &ldc, ScalarType WORK[], OrdinalType *const info) const |
BLAS 2 version of ORMQR; known workspace size. More... | |
void | UNMQR (const char &SIDE, const char &TRANS, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *C, const OrdinalType &ldc, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Apply Householder reflectors (complex case). More... | |
void | UNM2R (const char &SIDE, const char &TRANS, const OrdinalType &M, const OrdinalType &N, const OrdinalType &K, const ScalarType A[], const OrdinalType &LDA, const ScalarType TAU[], ScalarType C[], const OrdinalType &LDC, ScalarType WORK[], OrdinalType *const INFO) const |
BLAS 2 version of UNMQR; known workspace size. More... | |
void | ORGQR (const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Compute explicit Q factor from QR factorization (GEQRF) (real case). More... | |
void | UNGQR (const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Compute explicit QR factor from QR factorization (GEQRF) (complex case). More... | |
void | ORGHR (const OrdinalType &n, const OrdinalType &ilo, const OrdinalType &ihi, ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Generates a real orthogonal matrix Q which is the product of ihi-ilo elementary reflectors of order n , as returned by GEHRD. On return Q is stored in A . More... | |
void | ORMHR (const char &SIDE, const char &TRANS, const OrdinalType &m, const OrdinalType &n, const OrdinalType &ilo, const OrdinalType &ihi, const ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *C, const OrdinalType &ldc, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const |
Overwrites the general real m by n matrix C with the product of C and Q , which is a product of ihi-ilo elementary reflectors, as returned by GEHRD. More... | |
Triangular Matrix Routines | |
void | TREVC (const char &SIDE, const char &HOWMNY, OrdinalType *select, const OrdinalType &n, const ScalarType *T, const OrdinalType &ldt, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, const OrdinalType &mm, OrdinalType *m, ScalarType *WORK, OrdinalType *info) const |
void | TREVC (const char &SIDE, const OrdinalType &n, const ScalarType *T, const OrdinalType &ldt, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, const OrdinalType &mm, OrdinalType *m, ScalarType *WORK, MagnitudeType *RWORK, OrdinalType *info) const |
void | TREXC (const char &COMPQ, const OrdinalType &n, ScalarType *T, const OrdinalType &ldt, ScalarType *Q, const OrdinalType &ldq, OrdinalType *ifst, OrdinalType *ilst, ScalarType *WORK, OrdinalType *info) const |
void | TGEVC (const char &SIDE, const char &HOWMNY, const OrdinalType *SELECT, const OrdinalType &n, ScalarType *S, const OrdinalType &lds, ScalarType *P, const OrdinalType &ldp, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, const OrdinalType &mm, OrdinalType *M, ScalarType *WORK, OrdinalType *info) const |
Rotation/Reflection generators | |
void | LARTG (const ScalarType &f, const ScalarType &g, MagnitudeType *c, ScalarType *s, ScalarType *r) const |
Gnerates a plane rotation that zeros out the second component of the input vector. More... | |
void | LARFG (const OrdinalType &n, ScalarType *alpha, ScalarType *x, const OrdinalType &incx, ScalarType *tau) const |
Generates an elementary reflector of order n that zeros out the last n-1 components of the input vector. More... | |
Matrix Balancing Routines | |
void | GEBAL (const char &JOBZ, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *ilo, OrdinalType *ihi, MagnitudeType *scale, OrdinalType *info) const |
Balances a general matrix A, through similarity transformations to make the rows and columns as close in norm as possible. More... | |
void | GEBAK (const char &JOBZ, const char &SIDE, const OrdinalType &n, const OrdinalType &ilo, const OrdinalType &ihi, const MagnitudeType *scale, const OrdinalType &m, ScalarType *V, const OrdinalType &ldv, OrdinalType *info) const |
Forms the left or right eigenvectors of a general matrix that has been balanced by GEBAL by backward transformation of the computed eigenvectors V . More... | |
Random number generators | |
ScalarType | LARND (const OrdinalType &idist, OrdinalType *seed) const |
Returns a random number from a uniform or normal distribution. More... | |
void | LARNV (const OrdinalType &idist, OrdinalType *seed, const OrdinalType &n, ScalarType *v) const |
Returns a vector of random numbers from a chosen distribution. More... | |
Machine Characteristics Routines. | |
ScalarType | LAMCH (const char &CMACH) const |
Determines machine parameters for floating point characteristics. More... | |
OrdinalType | ILAENV (const OrdinalType &ispec, const std::string &NAME, const std::string &OPTS, const OrdinalType &N1=-1, const OrdinalType &N2=-1, const OrdinalType &N3=-1, const OrdinalType &N4=-1) const |
Chooses problem-dependent parameters for the local environment. More... | |
Miscellaneous Utilities. | |
ScalarType | LAPY2 (const ScalarType &x, const ScalarType &y) const |
Computes x^2 + y^2 safely, to avoid overflow. More... | |
The Templated LAPACK Wrapper Class.
The Teuchos::LAPACK class is a wrapper that encapsulates LAPACK (Linear Algebra Package). LAPACK provides portable, high- performance implementations of linear, eigen, SVD, etc solvers.
The standard LAPACK interface is Fortran-specific. Unfortunately, the interface between C++ and Fortran is not standard across all computer platforms. The Teuchos::LAPACK class provides C++ wrappers for the LAPACK kernels in order to insulate the rest of Teuchos from the details of C++ to Fortran translation. A Teuchos::LAPACK object is essentially nothing, but allows access to the LAPACK wrapper functions.
Teuchos::LAPACK is a serial interface only. This is appropriate since the standard LAPACK are only specified for serial execution (or shared memory parallel).
These templates are specialized to use the Fortran LAPACK routines for scalar types float
and double
.
If Teuchos is configured with -DTeuchos_ENABLE_COMPLEX
:BOOL=ON then these templates are specialized for scalar types std::complex<float>
and std::complex<double>
also.
http://www.netlib.org/lapack/
). Definition at line 96 of file Teuchos_LAPACK.hpp.
|
inline |
Default Constructor.
Definition at line 106 of file Teuchos_LAPACK.hpp.
|
inline |
Copy Constructor.
Definition at line 109 of file Teuchos_LAPACK.hpp.
|
inlinevirtual |
Destructor.
Definition at line 112 of file Teuchos_LAPACK.hpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::PTTRF | ( | const OrdinalType & | n, |
ScalarType * | d, | ||
ScalarType * | e, | ||
OrdinalType * | info | ||
) | const |
Computes the L*D*L'
factorization of a Hermitian/symmetric positive definite tridiagonal matrix A
.
Definition at line 114 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::PTTRS | ( | const OrdinalType & | n, |
const OrdinalType & | nrhs, | ||
const ScalarType * | d, | ||
const ScalarType * | e, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
OrdinalType * | info | ||
) | const |
Solves a tridiagonal system A*X=B
using the *D*L' factorization of A
computed by PTTRF.
Definition at line 118 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::POTRF | ( | const char & | UPLO, |
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | info | ||
) | const |
Computes Cholesky factorization of a real symmetric positive definite matrix A
.
Definition at line 122 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::POTRS | ( | const char & | UPLO, |
const OrdinalType & | n, | ||
const OrdinalType & | nrhs, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
OrdinalType * | info | ||
) | const |
Solves a system of linear equations A*X=B
, where A
is a symmetric positive definite matrix factored by POTRF and the nrhs
solutions are returned in B
.
Definition at line 126 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::POTRI | ( | const char & | UPLO, |
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | info | ||
) | const |
Computes the inverse of a real symmetric positive definite matrix A
using the Cholesky factorization A
from POTRF.
Definition at line 130 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::POCON | ( | const char & | UPLO, |
const OrdinalType & | n, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const ScalarType & | anorm, | ||
ScalarType * | rcond, | ||
ScalarType * | WORK, | ||
OrdinalType * | IWORK, | ||
OrdinalType * | info | ||
) | const |
Estimates the reciprocal of the condition number (1-norm) of a real symmetric positive definite matrix A
using the Cholesky factorization from POTRF.
Definition at line 134 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::POSV | ( | const char & | UPLO, |
const OrdinalType & | n, | ||
const OrdinalType & | nrhs, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
OrdinalType * | info | ||
) | const |
Computes the solution to a real system of linear equations A*X=B
, where A
is a symmetric positive definite matrix and the nrhs
solutions are returned in B
.
Definition at line 138 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::POEQU | ( | const OrdinalType & | n, |
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
MagnitudeType * | S, | ||
MagnitudeType * | scond, | ||
MagnitudeType * | amax, | ||
OrdinalType * | info | ||
) | const |
Computes row and column scalings intended to equilibrate a symmetric positive definite matrix A
and reduce its condition number (w.r.t. 2-norm).
Definition at line 142 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::PORFS | ( | const char & | UPLO, |
const OrdinalType & | n, | ||
const OrdinalType & | nrhs, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const ScalarType * | AF, | ||
const OrdinalType & | ldaf, | ||
const ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
ScalarType * | X, | ||
const OrdinalType & | ldx, | ||
ScalarType * | FERR, | ||
ScalarType * | BERR, | ||
ScalarType * | WORK, | ||
OrdinalType * | IWORK, | ||
OrdinalType * | info | ||
) | const |
Improves the computed solution to a system of linear equations when the coefficient matrix is symmetric positive definite, and provides error bounds and backward error estimates for the solution.
Definition at line 146 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::POSVX | ( | const char & | FACT, |
const char & | UPLO, | ||
const OrdinalType & | n, | ||
const OrdinalType & | nrhs, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | AF, | ||
const OrdinalType & | ldaf, | ||
char * | EQUED, | ||
ScalarType * | S, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
ScalarType * | X, | ||
const OrdinalType & | ldx, | ||
ScalarType * | rcond, | ||
ScalarType * | FERR, | ||
ScalarType * | BERR, | ||
ScalarType * | WORK, | ||
OrdinalType * | IWORK, | ||
OrdinalType * | info | ||
) | const |
Uses the Cholesky factorization to compute the solution to a real system of linear equations A*X=B
, where A
is symmetric positive definite. System can be equilibrated by POEQU and iteratively refined by PORFS, if requested.
Definition at line 149 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GELS | ( | const char & | TRANS, |
const OrdinalType & | m, | ||
const OrdinalType & | n, | ||
const OrdinalType & | nrhs, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Solves an over/underdetermined real m
by n
linear system A
using QR or LQ factorization of A.
Definition at line 153 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GELSS | ( | const OrdinalType & | m, |
const OrdinalType & | n, | ||
const OrdinalType & | nrhs, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
MagnitudeType * | S, | ||
const MagnitudeType | rcond, | ||
OrdinalType * | rank, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
MagnitudeType * | RWORK, | ||
OrdinalType * | info | ||
) | const |
Use the SVD to solve a possibly rank-deficient linear least-squares problem.
GELSS uses the singular value decomposition (SVD) to compute the minimum-norm solution to a possibly rank-deficient linear least-squares problem. The problem may be under- or overdetermined.
LAPACK's _GELSS routines take different arguments, depending on whether they are for real or complex arithmetic. This is because _GELSS imitates the interface of LAPACK's SVD routine. LAPACK's SVD routine takes an additional RWORK workspace array argument for COMPLEX*8 (CGELSS) and COMPLEX*16 (ZGELSS). LAPACK's real SVD routines (SGELSS and DGELSS) do not take the RWORK argument.
This class had already exposed GELSS for ScalarType = float and double that does not include an RWORK argument. Backwards compatibility requirements prevent us from simply changing that interface. We could provide a different interface for LAPACK specializations with ScalarType = std::complex<T>, but that would make the GELSS interface not generic at compile time. This would make using GELSS in generic code harder (for example, you would need to specialize code that uses GELSS on a Boolean, which specifies whether ScalarType is complex).
We fix this problem by providing an overloaded generic GELSS interface that does take an RWORK argument. This does not change the existing interface, but provides the additional capability to solve complex-valued least-squares problems. The RWORK argument is ignored when ScalarType is real, and may therefore be set to NULL in that case.
Definition at line 682 of file Teuchos_LAPACK.hpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GELSS | ( | const OrdinalType & | m, |
const OrdinalType & | n, | ||
const OrdinalType & | nrhs, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
ScalarType * | S, | ||
const ScalarType & | rcond, | ||
OrdinalType * | rank, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Legacy GELSS interface for real-valued ScalarType.
Definition at line 688 of file Teuchos_LAPACK.hpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GGLSE | ( | const OrdinalType & | m, |
const OrdinalType & | n, | ||
const OrdinalType & | p, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
ScalarType * | C, | ||
ScalarType * | D, | ||
ScalarType * | X, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Solves the linear equality-constrained least squares (LSE) problem where A
is an m
by n
matrix,B
is a p
by n
matrix C
is a given m-vector
, and D is a given p-vector
.
Definition at line 166 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GEQRF | ( | const OrdinalType & | m, |
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | TAU, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Computes a QR factorization of a general m
by n
matrix A
.
Definition at line 170 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GEQR2 | ( | const OrdinalType & | m, |
const OrdinalType & | n, | ||
ScalarType | A[], | ||
const OrdinalType & | lda, | ||
ScalarType | TAU[], | ||
ScalarType | WORK[], | ||
OrdinalType *const | info | ||
) | const |
BLAS 2 version of GEQRF, with known workspace size.
Definition at line 173 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GETRF | ( | const OrdinalType & | m, |
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | IPIV, | ||
OrdinalType * | info | ||
) | const |
Computes an LU factorization of a general m
by n
matrix A
using partial pivoting with row interchanges.
Definition at line 178 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GETRS | ( | const char & | TRANS, |
const OrdinalType & | n, | ||
const OrdinalType & | nrhs, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const OrdinalType * | IPIV, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
OrdinalType * | info | ||
) | const |
Solves a system of linear equations A*X=B
or A'*X=B
with a general n
by n
matrix A
using the LU factorization computed by GETRF.
Definition at line 182 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::LASCL | ( | const char & | TYPE, |
const OrdinalType & | kl, | ||
const OrdinalType & | ku, | ||
const MagnitudeType | cfrom, | ||
const MagnitudeType | cto, | ||
const OrdinalType & | m, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | info | ||
) | const |
Multiplies the m
by n
matrix A
by the real scalar cto/cfrom
.
Definition at line 186 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GEQP3 | ( | const OrdinalType & | m, |
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | jpvt, | ||
ScalarType * | TAU, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
MagnitudeType * | RWORK, | ||
OrdinalType * | info | ||
) | const |
Computes a QR factorization with column pivoting of a matrix A:
A*P = Q*R using Level 3 BLAS.
Definition at line 191 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::LASWP | ( | const OrdinalType & | N, |
ScalarType | A[], | ||
const OrdinalType & | LDA, | ||
const OrdinalType & | K1, | ||
const OrdinalType & | K2, | ||
const OrdinalType | IPIV[], | ||
const OrdinalType & | INCX | ||
) | const |
Apply a series of row interchanges to the matrix A.
Definition at line 207 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GBTRF | ( | const OrdinalType & | m, |
const OrdinalType & | n, | ||
const OrdinalType & | kl, | ||
const OrdinalType & | ku, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | IPIV, | ||
OrdinalType * | info | ||
) | const |
Computes an LU factorization of a general banded m
by n
matrix A
using partial pivoting with row interchanges.
Definition at line 218 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GBTRS | ( | const char & | TRANS, |
const OrdinalType & | n, | ||
const OrdinalType & | kl, | ||
const OrdinalType & | ku, | ||
const OrdinalType & | nrhs, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const OrdinalType * | IPIV, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
OrdinalType * | info | ||
) | const |
Solves a system of linear equations A*X=B
or A'*X=B
with a general banded n
by n
matrix A
using the LU factorization computed by GBTRF.
Definition at line 222 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GTTRF | ( | const OrdinalType & | n, |
ScalarType * | dl, | ||
ScalarType * | d, | ||
ScalarType * | du, | ||
ScalarType * | du2, | ||
OrdinalType * | IPIV, | ||
OrdinalType * | info | ||
) | const |
Computes an LU factorization of a n
by n
tridiagonal matrix A
using partial pivoting with row interchanges.
Definition at line 226 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GTTRS | ( | const char & | TRANS, |
const OrdinalType & | n, | ||
const OrdinalType & | nrhs, | ||
const ScalarType * | dl, | ||
const ScalarType * | d, | ||
const ScalarType * | du, | ||
const ScalarType * | du2, | ||
const OrdinalType * | IPIV, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
OrdinalType * | info | ||
) | const |
Solves a system of linear equations A*X=B
or A'*X=B
or A^H*X=B
with a tridiagonal matrix A
using the LU factorization computed by GTTRF.
Definition at line 230 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GETRI | ( | const OrdinalType & | n, |
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const OrdinalType * | IPIV, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Computes the inverse of a matrix A
using the LU factorization computed by GETRF.
Definition at line 234 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::LATRS | ( | const char & | UPLO, |
const char & | TRANS, | ||
const char & | DIAG, | ||
const char & | NORMIN, | ||
const OrdinalType & | N, | ||
ScalarType * | A, | ||
const OrdinalType & | LDA, | ||
ScalarType * | X, | ||
MagnitudeType * | SCALE, | ||
MagnitudeType * | CNORM, | ||
OrdinalType * | INFO | ||
) | const |
Robustly solve a possibly singular triangular linear system.
Definition at line 238 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GECON | ( | const char & | NORM, |
const OrdinalType & | n, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const ScalarType & | anorm, | ||
ScalarType * | rcond, | ||
ScalarType * | WORK, | ||
OrdinalType * | IWORK, | ||
OrdinalType * | info | ||
) | const |
Estimates the reciprocal of the condition number of a general real matrix A
, in either the 1-norm or the infinity-norm, using the LU factorization computed by GETRF.
Definition at line 264 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GBCON | ( | const char & | NORM, |
const OrdinalType & | n, | ||
const OrdinalType & | kl, | ||
const OrdinalType & | ku, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | IPIV, | ||
const ScalarType & | anorm, | ||
ScalarType * | rcond, | ||
ScalarType * | WORK, | ||
OrdinalType * | IWORK, | ||
OrdinalType * | info | ||
) | const |
Estimates the reciprocal of the condition number of a general banded real matrix A
, in either the 1-norm or the infinity-norm, using the LU factorization computed by GETRF.
Definition at line 268 of file Teuchos_LAPACK.cpp.
ScalarTraits< ScalarType >::magnitudeType Teuchos::LAPACK< OrdinalType, ScalarType >::LANGB | ( | const char & | NORM, |
const OrdinalType & | n, | ||
const OrdinalType & | kl, | ||
const OrdinalType & | ku, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
MagnitudeType * | WORK | ||
) | const |
Returns the value of the one norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of an n
by n
band matrix A
, with kl
sub-diagonals and ku
super-diagonals.
Definition at line 272 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GESV | ( | const OrdinalType & | n, |
const OrdinalType & | nrhs, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | IPIV, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
OrdinalType * | info | ||
) | const |
Computes the solution to a real system of linear equations A*X=B
, where A
is factored through GETRF and the nrhs
solutions are computed through GETRS.
Definition at line 276 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GEEQU | ( | const OrdinalType & | m, |
const OrdinalType & | n, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | R, | ||
ScalarType * | C, | ||
ScalarType * | rowcond, | ||
ScalarType * | colcond, | ||
ScalarType * | amax, | ||
OrdinalType * | info | ||
) | const |
Computes row and column scalings intended to equilibrate an m
by n
matrix A
and reduce its condition number.
Definition at line 280 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GERFS | ( | const char & | TRANS, |
const OrdinalType & | n, | ||
const OrdinalType & | nrhs, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const ScalarType * | AF, | ||
const OrdinalType & | ldaf, | ||
const OrdinalType * | IPIV, | ||
const ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
ScalarType * | X, | ||
const OrdinalType & | ldx, | ||
ScalarType * | FERR, | ||
ScalarType * | BERR, | ||
ScalarType * | WORK, | ||
OrdinalType * | IWORK, | ||
OrdinalType * | info | ||
) | const |
Improves the computed solution to a system of linear equations and provides error bounds and backward error estimates for the solution. Use after GETRF/GETRS.
Definition at line 284 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GBEQU | ( | const OrdinalType & | m, |
const OrdinalType & | n, | ||
const OrdinalType & | kl, | ||
const OrdinalType & | ku, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
MagnitudeType * | R, | ||
MagnitudeType * | C, | ||
MagnitudeType * | rowcond, | ||
MagnitudeType * | colcond, | ||
MagnitudeType * | amax, | ||
OrdinalType * | info | ||
) | const |
Computes row and column scalings intended to equilibrate an m
by n
banded matrix A
and reduce its condition number.
Definition at line 288 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GBRFS | ( | const char & | TRANS, |
const OrdinalType & | n, | ||
const OrdinalType & | kl, | ||
const OrdinalType & | ku, | ||
const OrdinalType & | nrhs, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const ScalarType * | AF, | ||
const OrdinalType & | ldaf, | ||
const OrdinalType * | IPIV, | ||
const ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
ScalarType * | X, | ||
const OrdinalType & | ldx, | ||
ScalarType * | FERR, | ||
ScalarType * | BERR, | ||
ScalarType * | WORK, | ||
OrdinalType * | IWORK, | ||
OrdinalType * | info | ||
) | const |
Improves the computed solution to a banded system of linear equations and provides error bounds and backward error estimates for the solution. Use after GBTRF/GBTRS.
Definition at line 292 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GESVX | ( | const char & | FACT, |
const char & | TRANS, | ||
const OrdinalType & | n, | ||
const OrdinalType & | nrhs, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | AF, | ||
const OrdinalType & | ldaf, | ||
OrdinalType * | IPIV, | ||
char * | EQUED, | ||
ScalarType * | R, | ||
ScalarType * | C, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
ScalarType * | X, | ||
const OrdinalType & | ldx, | ||
ScalarType * | rcond, | ||
ScalarType * | FERR, | ||
ScalarType * | BERR, | ||
ScalarType * | WORK, | ||
OrdinalType * | IWORK, | ||
OrdinalType * | info | ||
) | const |
Uses the LU factorization to compute the solution to a real system of linear equations A*X=B
, returning error bounds on the solution and a condition estimate.
Definition at line 295 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::SYTRD | ( | const char & | UPLO, |
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | D, | ||
ScalarType * | E, | ||
ScalarType * | TAU, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Reduces a real symmetric matrix A
to tridiagonal form by orthogonal similarity transformations.
std::complex<float>
or std::complex<double>
. Definition at line 299 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GEHRD | ( | const OrdinalType & | n, |
const OrdinalType & | ilo, | ||
const OrdinalType & | ihi, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | TAU, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Reduces a real general matrix A
to upper Hessenberg form by orthogonal similarity transformations.
Definition at line 303 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::TRTRS | ( | const char & | UPLO, |
const char & | TRANS, | ||
const char & | DIAG, | ||
const OrdinalType & | n, | ||
const OrdinalType & | nrhs, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
OrdinalType * | info | ||
) | const |
Solves a triangular linear system of the form A*X=B
or A**T*X=B
, where A
is a triangular matrix.
Definition at line 307 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::TRTRI | ( | const char & | UPLO, |
const char & | DIAG, | ||
const OrdinalType & | n, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | info | ||
) | const |
Computes the inverse of an upper or lower triangular matrix A
.
Definition at line 311 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::SPEV | ( | const char & | JOBZ, |
const char & | UPLO, | ||
const OrdinalType & | n, | ||
ScalarType * | AP, | ||
ScalarType * | W, | ||
ScalarType * | Z, | ||
const OrdinalType & | ldz, | ||
ScalarType * | WORK, | ||
OrdinalType * | info | ||
) | const |
Computes the eigenvalues and, optionally, eigenvectors of a symmetric n
by n
matrix A
in packed storage.
std::complex<float>
or std::complex<double>
. Definition at line 315 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::SYEV | ( | const char & | JOBZ, |
const char & | UPLO, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | W, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Computes all the eigenvalues and, optionally, eigenvectors of a symmetric n
by n
matrix A.
std::complex<float>
or std::complex<double>
. Definition at line 319 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::SYGV | ( | const OrdinalType & | itype, |
const char & | JOBZ, | ||
const char & | UPLO, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
ScalarType * | W, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Computes all the eigenvalues and, optionally, eigenvectors of a symmetric n
by n
matrix pencil {A
,B}, where A
is symmetric and B
is symmetric positive-definite.
std::complex<float>
or std::complex<double>
. Definition at line 323 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::HEEV | ( | const char & | JOBZ, |
const char & | UPLO, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
MagnitudeType * | W, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
MagnitudeType * | RWORK, | ||
OrdinalType * | info | ||
) | const |
Computes all the eigenvalues and, optionally, eigenvectors of a Hermitian n
by n
matrix A.
float
or double
. Definition at line 327 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::HEGV | ( | const OrdinalType & | itype, |
const char & | JOBZ, | ||
const char & | UPLO, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
MagnitudeType * | W, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
MagnitudeType * | RWORK, | ||
OrdinalType * | info | ||
) | const |
Computes all the eigenvalues and, optionally, eigenvectors of a generalized Hermitian-definite n
by n
matrix pencil {A
,B}, where A
is Hermitian and B
is Hermitian positive-definite.
float
or double
. Definition at line 331 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::STEQR | ( | const char & | COMPZ, |
const OrdinalType & | n, | ||
MagnitudeType * | D, | ||
MagnitudeType * | E, | ||
ScalarType * | Z, | ||
const OrdinalType & | ldz, | ||
MagnitudeType * | WORK, | ||
OrdinalType * | info | ||
) | const |
Computes the eigenvalues and, optionally, eigenvectors of a symmetric tridiagonal n
by n
matrix A
using implicit QL/QR. The eigenvectors can only be computed if A
was reduced to tridiagonal form by SYTRD.
Definition at line 335 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::PTEQR | ( | const char & | COMPZ, |
const OrdinalType & | n, | ||
MagnitudeType * | D, | ||
MagnitudeType * | E, | ||
ScalarType * | Z, | ||
const OrdinalType & | ldz, | ||
MagnitudeType * | WORK, | ||
OrdinalType * | info | ||
) | const |
Computes the eigenvalues and, optionally, eigenvectors of a symmetric positive-definite tridiagonal n
by n
matrix A
using BDSQR, after factoring the matrix with PTTRF.
Definition at line 339 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::HSEQR | ( | const char & | JOB, |
const char & | COMPZ, | ||
const OrdinalType & | n, | ||
const OrdinalType & | ilo, | ||
const OrdinalType & | ihi, | ||
ScalarType * | H, | ||
const OrdinalType & | ldh, | ||
ScalarType * | WR, | ||
ScalarType * | WI, | ||
ScalarType * | Z, | ||
const OrdinalType & | ldz, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Computes the eigenvalues of a real upper Hessenberg matrix H
and, optionally, the matrices T
and Z
from the Schur decomposition, where T is an upper quasi-triangular matrix and Z contains the Schur vectors.
Definition at line 343 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GEES | ( | const char & | JOBVS, |
const char & | SORT, | ||
OrdinalType &(*)(ScalarType *, ScalarType *) | ptr2func, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | sdim, | ||
ScalarType * | WR, | ||
ScalarType * | WI, | ||
ScalarType * | VS, | ||
const OrdinalType & | ldvs, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | BWORK, | ||
OrdinalType * | info | ||
) | const |
Computes for an n
by n
nonsymmetric matrix A
, the eigenvalues, the Schur form T
, and, optionally, the matrix of Schur vectors Z
. When ScalarType
is float
or double
, the real Schur form is computed.
float&
and double
, where select
requires two arguments to represent a complex eigenvalue.) Definition at line 1192 of file Teuchos_LAPACK.hpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GEES | ( | const char & | JOBVS, |
const char & | SORT, | ||
OrdinalType &(*)(ScalarType *) | ptr2func, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | sdim, | ||
ScalarType * | W, | ||
ScalarType * | VS, | ||
const OrdinalType & | ldvs, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
MagnitudeType * | RWORK, | ||
OrdinalType * | BWORK, | ||
OrdinalType * | info | ||
) | const |
Computes for an n
by n
nonsymmetric matrix A
, the eigenvalues, the Schur form T
, and, optionally, the matrix of Schur vectors Z
. When ScalarType
is float
or double
, the real Schur form is computed.
std::complex<float>
and std::complex<double>
, where select
requires one arguments to represent a complex eigenvalue.) Definition at line 1198 of file Teuchos_LAPACK.hpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GEES | ( | const char & | JOBVS, |
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | sdim, | ||
MagnitudeType * | WR, | ||
MagnitudeType * | WI, | ||
ScalarType * | VS, | ||
const OrdinalType & | ldvs, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
MagnitudeType * | RWORK, | ||
OrdinalType * | BWORK, | ||
OrdinalType * | info | ||
) | const |
Computes for an n
by n
nonsymmetric matrix A
, the eigenvalues, the Schur form T
, and, optionally, the matrix of Schur vectors Z
. When ScalarType
is float
or double
, the real Schur form is computed.
ScalarType
, when the user doesn't want to enable the sorting functionality.) Definition at line 1204 of file Teuchos_LAPACK.hpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GEEV | ( | const char & | JOBVL, |
const char & | JOBVR, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
MagnitudeType * | WR, | ||
MagnitudeType * | WI, | ||
ScalarType * | VL, | ||
const OrdinalType & | ldvl, | ||
ScalarType * | VR, | ||
const OrdinalType & | ldvr, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
MagnitudeType * | RWORK, | ||
OrdinalType * | info | ||
) | const |
Computes for an n
by n
real nonsymmetric matrix A
, the eigenvalues and, optionally, the left and/or right eigenvectors.
Real and imaginary parts of the eigenvalues are returned in separate arrays, WR for real and WI for complex. The RWORK array is only referenced if ScalarType is complex.
Definition at line 359 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GEEVX | ( | const char & | BALANC, |
const char & | JOBVL, | ||
const char & | JOBVR, | ||
const char & | SENSE, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | WR, | ||
ScalarType * | WI, | ||
ScalarType * | VL, | ||
const OrdinalType & | ldvl, | ||
ScalarType * | VR, | ||
const OrdinalType & | ldvr, | ||
OrdinalType * | ilo, | ||
OrdinalType * | ihi, | ||
MagnitudeType * | SCALE, | ||
MagnitudeType * | abnrm, | ||
MagnitudeType * | RCONDE, | ||
MagnitudeType * | RCONDV, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | IWORK, | ||
OrdinalType * | info | ||
) | const |
Computes for an n
by n
real nonsymmetric matrix A
, the eigenvalues and, optionally, the left and/or right eigenvectors. Optionally, it can compute a balancing transformation to improve the conditioning of the eigenvalues and eigenvectors.
ScalarType&
= float&
or double
.) Definition at line 372 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GGEVX | ( | const char & | BALANC, |
const char & | JOBVL, | ||
const char & | JOBVR, | ||
const char & | SENSE, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
MagnitudeType * | ALPHAR, | ||
MagnitudeType * | ALPHAI, | ||
ScalarType * | BETA, | ||
ScalarType * | VL, | ||
const OrdinalType & | ldvl, | ||
ScalarType * | VR, | ||
const OrdinalType & | ldvr, | ||
OrdinalType * | ilo, | ||
OrdinalType * | ihi, | ||
MagnitudeType * | lscale, | ||
MagnitudeType * | rscale, | ||
MagnitudeType * | abnrm, | ||
MagnitudeType * | bbnrm, | ||
MagnitudeType * | RCONDE, | ||
MagnitudeType * | RCONDV, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | IWORK, | ||
OrdinalType * | BWORK, | ||
OrdinalType * | info | ||
) | const |
Computes for a pair of n
by n
nonsymmetric matrices (A
,B
) the generalized eigenvalues, and optionally, the left and/or right generalized eigenvectors. Optionally, it can compute a balancing transformation to improve the conditioning of the eigenvalues and eigenvectors.
ScalarType
= float
or double
.) Definition at line 376 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GGEV | ( | const char & | JOBVL, |
const char & | JOBVR, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
MagnitudeType * | ALPHAR, | ||
MagnitudeType * | ALPHAI, | ||
ScalarType * | BETA, | ||
ScalarType * | VL, | ||
const OrdinalType & | ldvl, | ||
ScalarType * | VR, | ||
const OrdinalType & | ldvr, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Computes for a pair of n
by n
nonsymmetric matrices (A
,B
) the generalized eigenvalues, and optionally, the left and/or right generalized eigenvectors.
ScalarType
= float
or double
.) Definition at line 384 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::TRSEN | ( | const char & | JOB, |
const char & | COMPQ, | ||
const OrdinalType * | SELECT, | ||
const OrdinalType & | n, | ||
ScalarType * | T, | ||
const OrdinalType & | ldt, | ||
ScalarType * | Q, | ||
const OrdinalType & | ldq, | ||
MagnitudeType * | WR, | ||
MagnitudeType * | WI, | ||
OrdinalType * | M, | ||
ScalarType * | S, | ||
MagnitudeType * | SEP, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | IWORK, | ||
const OrdinalType & | liwork, | ||
OrdinalType * | info | ||
) | const |
Reorders the real Schur factorization of a real matrix so that a selected cluster of eigenvalues appears in the leading diagonal blocks of the upper quasi-triangular matrix T
, and the leading columns of Q
form an orthonormal basis of the corresponding right invariant subspace.
ScalarType
= float
or double
.) Definition at line 388 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::TGSEN | ( | const OrdinalType & | ijob, |
const OrdinalType & | wantq, | ||
const OrdinalType & | wantz, | ||
const OrdinalType * | SELECT, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
MagnitudeType * | ALPHAR, | ||
MagnitudeType * | ALPHAI, | ||
MagnitudeType * | BETA, | ||
ScalarType * | Q, | ||
const OrdinalType & | ldq, | ||
ScalarType * | Z, | ||
const OrdinalType & | ldz, | ||
OrdinalType * | M, | ||
MagnitudeType * | PL, | ||
MagnitudeType * | PR, | ||
MagnitudeType * | DIF, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | IWORK, | ||
const OrdinalType & | liwork, | ||
OrdinalType * | info | ||
) | const |
Reorders the generalized real Schur decomposition of a real matrix pair (A
, B
), so that a selected cluster of eigenvalues appears in the leading diagonal blocks of the upper quasi-triangular matrix A
and the upper triangular B
.
ScalarType
= float
or double
.) Definition at line 392 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GGES | ( | const char & | JOBVL, |
const char & | JOBVR, | ||
const char & | SORT, | ||
OrdinalType &(*)(ScalarType *, ScalarType *, ScalarType *) | ptr2func, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
ScalarType * | B, | ||
const OrdinalType & | ldb, | ||
OrdinalType * | sdim, | ||
MagnitudeType * | ALPHAR, | ||
MagnitudeType * | ALPHAI, | ||
MagnitudeType * | BETA, | ||
ScalarType * | VL, | ||
const OrdinalType & | ldvl, | ||
ScalarType * | VR, | ||
const OrdinalType & | ldvr, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | BWORK, | ||
OrdinalType * | info | ||
) | const |
Computes for a pair of n
by n
nonsymmetric matrices (A
,B
) the generalized eigenvalues, the generalized real Schur form (S
,T
), optionally, the left and/or right matrices of Schur vectors.
ScalarType
= float
or double
.) Definition at line 396 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GESVD | ( | const char & | JOBU, |
const char & | JOBVT, | ||
const OrdinalType & | m, | ||
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
MagnitudeType * | S, | ||
ScalarType * | U, | ||
const OrdinalType & | ldu, | ||
ScalarType * | V, | ||
const OrdinalType & | ldv, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
MagnitudeType * | RWORK, | ||
OrdinalType * | info | ||
) | const |
Computes the singular values (and optionally, vectors) of a real matrix A
.
Definition at line 368 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::ORMQR | ( | const char & | SIDE, |
const char & | TRANS, | ||
const OrdinalType & | m, | ||
const OrdinalType & | n, | ||
const OrdinalType & | k, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const ScalarType * | TAU, | ||
ScalarType * | C, | ||
const OrdinalType & | ldc, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Apply Householder reflectors (real case).
Overwrite the general real m
by n
matrix C
with the product of Q
and C
, whiere Q is the product of k
elementary (Householder) reflectors as returned by GEQRF.
Definition at line 400 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::ORM2R | ( | const char & | SIDE, |
const char & | TRANS, | ||
const OrdinalType & | m, | ||
const OrdinalType & | n, | ||
const OrdinalType & | k, | ||
const ScalarType | A[], | ||
const OrdinalType & | lda, | ||
const ScalarType | TAU[], | ||
ScalarType | C[], | ||
const OrdinalType & | ldc, | ||
ScalarType | WORK[], | ||
OrdinalType *const | info | ||
) | const |
BLAS 2 version of ORMQR; known workspace size.
Definition at line 404 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::UNMQR | ( | const char & | SIDE, |
const char & | TRANS, | ||
const OrdinalType & | m, | ||
const OrdinalType & | n, | ||
const OrdinalType & | k, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const ScalarType * | TAU, | ||
ScalarType * | C, | ||
const OrdinalType & | ldc, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Apply Householder reflectors (complex case).
Overwrite the general complex m
by n
matrix C
with the product of Q
and C
, where Q is the product of k
elementary (Householder) reflectors as returned by GEQRF.
Definition at line 408 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::UNM2R | ( | const char & | SIDE, |
const char & | TRANS, | ||
const OrdinalType & | M, | ||
const OrdinalType & | N, | ||
const OrdinalType & | K, | ||
const ScalarType | A[], | ||
const OrdinalType & | LDA, | ||
const ScalarType | TAU[], | ||
ScalarType | C[], | ||
const OrdinalType & | LDC, | ||
ScalarType | WORK[], | ||
OrdinalType *const | INFO | ||
) | const |
BLAS 2 version of UNMQR; known workspace size.
Definition at line 416 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::ORGQR | ( | const OrdinalType & | m, |
const OrdinalType & | n, | ||
const OrdinalType & | k, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const ScalarType * | TAU, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Compute explicit Q factor from QR factorization (GEQRF) (real case).
Generate the m
by n
matrix Q with orthonormal columns corresponding to the first n
columns of a product of k
elementary reflectors of order m
, as returned by GEQRF
.
UNGQR
in that case. ("OR" stands for "orthogonal"; "UN" stands for "unitary.") Definition at line 425 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::UNGQR | ( | const OrdinalType & | m, |
const OrdinalType & | n, | ||
const OrdinalType & | k, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const ScalarType * | TAU, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Compute explicit QR factor from QR factorization (GEQRF) (complex case).
Generate the m
by n
matrix Q with orthonormal columns corresponding tothe first n
columns of a product of k
elementary reflectors of order m
, as returned by GEQRF
.
ORGQR
when ScalarType is real. (Unitary real matrices are orthogonal.) Definition at line 429 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::ORGHR | ( | const OrdinalType & | n, |
const OrdinalType & | ilo, | ||
const OrdinalType & | ihi, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const ScalarType * | TAU, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Generates a real orthogonal matrix Q
which is the product of ihi-ilo
elementary reflectors of order n
, as returned by GEHRD. On return Q
is stored in A
.
Definition at line 433 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::ORMHR | ( | const char & | SIDE, |
const char & | TRANS, | ||
const OrdinalType & | m, | ||
const OrdinalType & | n, | ||
const OrdinalType & | ilo, | ||
const OrdinalType & | ihi, | ||
const ScalarType * | A, | ||
const OrdinalType & | lda, | ||
const ScalarType * | TAU, | ||
ScalarType * | C, | ||
const OrdinalType & | ldc, | ||
ScalarType * | WORK, | ||
const OrdinalType & | lwork, | ||
OrdinalType * | info | ||
) | const |
Overwrites the general real m
by n
matrix C
with the product of C
and Q
, which is a product of ihi-ilo
elementary reflectors, as returned by GEHRD.
Definition at line 437 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::TREVC | ( | const char & | SIDE, |
const char & | HOWMNY, | ||
OrdinalType * | select, | ||
const OrdinalType & | n, | ||
const ScalarType * | T, | ||
const OrdinalType & | ldt, | ||
ScalarType * | VL, | ||
const OrdinalType & | ldvl, | ||
ScalarType * | VR, | ||
const OrdinalType & | ldvr, | ||
const OrdinalType & | mm, | ||
OrdinalType * | m, | ||
ScalarType * | WORK, | ||
OrdinalType * | info | ||
) | const |
Computes some or all of the right and/or left eigenvectors of an upper triangular matrix T
. If ScalarType is float
or double
, then the matrix is quasi-triangular and arugments RWORK
is ignored.
Definition at line 1309 of file Teuchos_LAPACK.hpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::TREVC | ( | const char & | SIDE, |
const OrdinalType & | n, | ||
const ScalarType * | T, | ||
const OrdinalType & | ldt, | ||
ScalarType * | VL, | ||
const OrdinalType & | ldvl, | ||
ScalarType * | VR, | ||
const OrdinalType & | ldvr, | ||
const OrdinalType & | mm, | ||
OrdinalType * | m, | ||
ScalarType * | WORK, | ||
MagnitudeType * | RWORK, | ||
OrdinalType * | info | ||
) | const |
Computes some or all of the right and/or left eigenvectors of an upper triangular matrix T
. If ScalarType& is float
or double
, then the matrix is quasi-triangular and arugments RWORK
is ignored.
ScalarType
, when the user doesn't want to enable the selecting functionality, with HOWMNY='A'.) Definition at line 1315 of file Teuchos_LAPACK.hpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::TREXC | ( | const char & | COMPQ, |
const OrdinalType & | n, | ||
ScalarType * | T, | ||
const OrdinalType & | ldt, | ||
ScalarType * | Q, | ||
const OrdinalType & | ldq, | ||
OrdinalType * | ifst, | ||
OrdinalType * | ilst, | ||
ScalarType * | WORK, | ||
OrdinalType * | info | ||
) | const |
Reorders the Schur factorization of a matrix T
via unitary similarity transformations so that the diagonal element of T
with row index ifst
is moved to row ilst
. If ScalarType
is float
or double
, then T
should be in real Schur form and the operation affects the diagonal block referenced by ifst
.
std::complex<float>
or std::complex<double>
. Definition at line 452 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::TGEVC | ( | const char & | SIDE, |
const char & | HOWMNY, | ||
const OrdinalType * | SELECT, | ||
const OrdinalType & | n, | ||
ScalarType * | S, | ||
const OrdinalType & | lds, | ||
ScalarType * | P, | ||
const OrdinalType & | ldp, | ||
ScalarType * | VL, | ||
const OrdinalType & | ldvl, | ||
ScalarType * | VR, | ||
const OrdinalType & | ldvr, | ||
const OrdinalType & | mm, | ||
OrdinalType * | M, | ||
ScalarType * | WORK, | ||
OrdinalType * | info | ||
) | const |
Computes some or all of the right and/or left eigenvectors of a pair of real matrices ( S
, P
), where S
is a quasi-triangular matrix and P
is upper triangular.
ScalarType
= float
or double
. Definition at line 456 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::LARTG | ( | const ScalarType & | f, |
const ScalarType & | g, | ||
MagnitudeType * | c, | ||
ScalarType * | s, | ||
ScalarType * | r | ||
) | const |
Gnerates a plane rotation that zeros out the second component of the input vector.
Definition at line 460 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::LARFG | ( | const OrdinalType & | n, |
ScalarType * | alpha, | ||
ScalarType * | x, | ||
const OrdinalType & | incx, | ||
ScalarType * | tau | ||
) | const |
Generates an elementary reflector of order n
that zeros out the last n-1
components of the input vector.
Definition at line 464 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GEBAL | ( | const char & | JOBZ, |
const OrdinalType & | n, | ||
ScalarType * | A, | ||
const OrdinalType & | lda, | ||
OrdinalType * | ilo, | ||
OrdinalType * | ihi, | ||
MagnitudeType * | scale, | ||
OrdinalType * | info | ||
) | const |
Balances a general matrix A, through similarity transformations to make the rows and columns as close in norm as possible.
Definition at line 467 of file Teuchos_LAPACK.cpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::GEBAK | ( | const char & | JOBZ, |
const char & | SIDE, | ||
const OrdinalType & | n, | ||
const OrdinalType & | ilo, | ||
const OrdinalType & | ihi, | ||
const MagnitudeType * | scale, | ||
const OrdinalType & | m, | ||
ScalarType * | V, | ||
const OrdinalType & | ldv, | ||
OrdinalType * | info | ||
) | const |
Forms the left or right eigenvectors of a general matrix that has been balanced by GEBAL by backward transformation of the computed eigenvectors V
.
Definition at line 471 of file Teuchos_LAPACK.cpp.
ScalarType Teuchos::LAPACK< OrdinalType, ScalarType >::LARND | ( | const OrdinalType & | idist, |
OrdinalType * | seed | ||
) | const |
Returns a random number from a uniform or normal distribution.
Definition at line 1378 of file Teuchos_LAPACK.hpp.
void Teuchos::LAPACK< OrdinalType, ScalarType >::LARNV | ( | const OrdinalType & | idist, |
OrdinalType * | seed, | ||
const OrdinalType & | n, | ||
ScalarType * | v | ||
) | const |
Returns a vector of random numbers from a chosen distribution.
Definition at line 479 of file Teuchos_LAPACK.cpp.
ScalarType Teuchos::LAPACK< OrdinalType, ScalarType >::LAMCH | ( | const char & | CMACH | ) | const |
Determines machine parameters for floating point characteristics.
std::complex<float>
or std::complex<double>
. Definition at line 483 of file Teuchos_LAPACK.cpp.
OrdinalType Teuchos::LAPACK< OrdinalType, ScalarType >::ILAENV | ( | const OrdinalType & | ispec, |
const std::string & | NAME, | ||
const std::string & | OPTS, | ||
const OrdinalType & | N1 = -1 , |
||
const OrdinalType & | N2 = -1 , |
||
const OrdinalType & | N3 = -1 , |
||
const OrdinalType & | N4 = -1 |
||
) | const |
Chooses problem-dependent parameters for the local environment.
Definition at line 487 of file Teuchos_LAPACK.cpp.
ScalarType Teuchos::LAPACK< OrdinalType, ScalarType >::LAPY2 | ( | const ScalarType & | x, |
const ScalarType & | y | ||
) | const |
Computes x^2 + y^2 safely, to avoid overflow.
std::complex<float>
or std::complex<double>
. Definition at line 500 of file Teuchos_LAPACK.cpp.