Teuchos - Trilinos Tools Package
Version of the Day
|
A class for solving dense linear problems. More...
#include <Teuchos_SerialQRDenseSolver.hpp>
Public Member Functions | |
Constructor/Destructor Methods | |
SerialQRDenseSolver () | |
Default constructor; matrix should be set using setMatrix(), LHS and RHS set with setVectors(). More... | |
virtual | ~SerialQRDenseSolver () |
SerialQRDenseSolver destructor. More... | |
Set Methods | |
int | setMatrix (const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &A) |
Sets the pointers for coefficient matrix. More... | |
int | setVectors (const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &X, const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &B) |
Sets the pointers for left and right hand side vector(s). More... | |
Strategy Modifying Methods | |
void | factorWithEquilibration (bool flag) |
Causes equilibration to be called just before the matrix factorization as part of the call to factor . More... | |
void | solveWithTranspose (bool flag) |
If flag is true, causes all subsequent function calls to work with the adjoint of this matrix, otherwise not. More... | |
void | solveWithTransposeFlag (Teuchos::ETransp trans) |
All subsequent function calls will work with the transpose-type set by this method (Teuchos::NO_TRANS or Teuchos::CONJ_TRANS). More... | |
Factor/Solve/Invert Methods | |
int | factor () |
Computes the in-place QR factorization of the matrix using the LAPACK routine _GETRF or the Eigen class HouseholderQR. More... | |
int | solve () |
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors().. More... | |
int | computeEquilibrateScaling () |
Determines if this matrix should be scaled. More... | |
int | equilibrateMatrix () |
Equilibrates the this matrix. More... | |
int | equilibrateRHS () |
Equilibrates the current RHS. More... | |
int | unequilibrateLHS () |
Unscales the solution vectors if equilibration was used to solve the system. More... | |
int | formQ () |
Explicitly forms the unitary matrix Q. More... | |
int | formR () |
Explicitly forms the upper triangular matrix R. More... | |
int | multiplyQ (ETransp transq, SerialDenseMatrix< OrdinalType, ScalarType > &C) |
Left multiply the input matrix by the unitary matrix Q or its adjoint. More... | |
int | solveR (ETransp transr, SerialDenseMatrix< OrdinalType, ScalarType > &C) |
Solve input matrix on the left with the upper triangular matrix R or its adjoint. More... | |
Query methods | |
bool | transpose () |
Returns true if adjoint of this matrix has and will be used. More... | |
bool | factored () |
Returns true if matrix is factored (factor available via getFactoredMatrix()). More... | |
bool | equilibratedA () |
Returns true if factor is equilibrated (factor available via getFactoredMatrix()). More... | |
bool | equilibratedB () |
Returns true if RHS is equilibrated (RHS available via getRHS()). More... | |
bool | shouldEquilibrate () |
Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system. More... | |
bool | solved () |
Returns true if the current set of vectors has been solved. More... | |
bool | formedQ () |
Returns true if Q has been formed explicitly. More... | |
bool | formedR () |
Returns true if R has been formed explicitly. More... | |
Data Accessor methods | |
RCP< SerialDenseMatrix < OrdinalType, ScalarType > > | getMatrix () const |
Returns pointer to current matrix. More... | |
RCP< SerialDenseMatrix < OrdinalType, ScalarType > > | getFactoredMatrix () const |
Returns pointer to factored matrix (assuming factorization has been performed). More... | |
RCP< SerialDenseMatrix < OrdinalType, ScalarType > > | getQ () const |
Returns pointer to Q (assuming factorization has been performed). More... | |
RCP< SerialDenseMatrix < OrdinalType, ScalarType > > | getR () const |
Returns pointer to R (assuming factorization has been performed). More... | |
RCP< SerialDenseMatrix < OrdinalType, ScalarType > > | getLHS () const |
Returns pointer to current LHS. More... | |
RCP< SerialDenseMatrix < OrdinalType, ScalarType > > | getRHS () const |
Returns pointer to current RHS. More... | |
OrdinalType | numRows () const |
Returns row dimension of system. More... | |
OrdinalType | numCols () const |
Returns column dimension of system. More... | |
std::vector< ScalarType > | tau () const |
Returns pointer to pivot vector (if factorization has been computed), zero otherwise. More... | |
MagnitudeType | ANORM () const |
Returns the absolute value of the largest element of this matrix (returns -1 if not yet computed). More... | |
I/O methods | |
void | Print (std::ostream &os) const |
Print service methods; defines behavior of ostream << operator. More... | |
Public Member Functions inherited from Teuchos::CompObject | |
CompObject () | |
Default constructor. More... | |
CompObject (const CompObject &source) | |
Copy Constructor. More... | |
virtual | ~CompObject () |
Destructor. More... | |
void | setFlopCounter (const Flops &FlopCounter) |
Set the internal Teuchos::Flops() pointer. More... | |
void | setFlopCounter (const CompObject &compObject) |
Set the internal Teuchos::Flops() pointer to the flop counter of another Teuchos::CompObject. More... | |
void | unsetFlopCounter () |
Set the internal Teuchos::Flops() pointer to 0 (no flops counted). More... | |
Flops * | getFlopCounter () const |
Get the pointer to the Teuchos::Flops() object associated with this object, returns 0 if none. More... | |
void | resetFlops () const |
Resets the number of floating point operations to zero for this multi-std::vector. More... | |
double | getFlops () const |
Returns the number of floating point operations with this multi-std::vector. More... | |
void | updateFlops (int addflops) const |
Increment Flop count for this object. More... | |
void | updateFlops (long int addflops) const |
Increment Flop count for this object. More... | |
void | updateFlops (double addflops) const |
Increment Flop count for this object. More... | |
void | updateFlops (float addflops) const |
Increment Flop count for this object. More... | |
Public Member Functions inherited from Teuchos::Object | |
Object (int tracebackModeIn=-1) | |
Default Constructor. More... | |
Object (const char *label, int tracebackModeIn=-1) | |
Labeling Constructor. More... | |
Object (const std::string &label, int tracebackModeIn=-1) | |
Create an Object with the given label, and optionally, with the given traceback mode. More... | |
virtual | ~Object () |
Destructor (virtual, for safety of derived classes). More... | |
virtual void | print (std::ostream &os) const |
Print the object to the given output stream. More... | |
virtual int | reportError (const std::string message, int errorCode) const |
Report an error with this Object. More... | |
virtual void | setLabel (const char *theLabel) |
virtual const char * | label () const |
Access the object's label (LEGACY; return std::string instead). More... | |
Public Member Functions inherited from Teuchos::BLAS< OrdinalType, ScalarType > | |
BLAS (void) | |
Default constructor. More... | |
BLAS (const BLAS< OrdinalType, ScalarType > &) | |
Copy constructor. More... | |
virtual | ~BLAS (void) |
Destructor. More... | |
Public Member Functions inherited from Teuchos::DefaultBLASImpl< OrdinalType, ScalarType > | |
DefaultBLASImpl (void) | |
Default constructor. More... | |
DefaultBLASImpl (const DefaultBLASImpl< OrdinalType, ScalarType > &) | |
Copy constructor. More... | |
virtual | ~DefaultBLASImpl (void) |
Destructor. More... | |
template<typename alpha_type , typename A_type , typename x_type , typename beta_type > | |
void | GEMV (ETransp trans, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const x_type *x, const OrdinalType &incx, const beta_type beta, ScalarType *y, const OrdinalType &incy) const |
Performs the matrix-vector operation: y <- alpha*A*x+beta*y or y <- alpha*A'*x+beta*y where A is a general m by n matrix. More... | |
template<typename A_type > | |
void | TRMV (EUplo uplo, ETransp trans, EDiag diag, const OrdinalType &n, const A_type *A, const OrdinalType &lda, ScalarType *x, const OrdinalType &incx) const |
Performs the matrix-vector operation: x <- A*x or x <- A'*x where A is a unit/non-unit n by n upper/lower triangular matrix. More... | |
template<typename alpha_type , typename x_type , typename y_type > | |
void | GER (const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const x_type *x, const OrdinalType &incx, const y_type *y, const OrdinalType &incy, ScalarType *A, const OrdinalType &lda) const |
Performs the rank 1 operation: A <- alpha*x*y'+A . More... | |
template<typename alpha_type , typename A_type , typename B_type , typename beta_type > | |
void | GEMM (ETransp transa, ETransp transb, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const |
General matrix-matrix multiply. More... | |
void | SWAP (const OrdinalType &n, ScalarType *const x, const OrdinalType &incx, ScalarType *const y, const OrdinalType &incy) const |
Swap the entries of x and y. More... | |
template<typename alpha_type , typename A_type , typename B_type , typename beta_type > | |
void | SYMM (ESide side, EUplo uplo, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const |
Performs the matrix-matrix operation: C <- alpha*A*B+beta*C or C <- alpha*B*A+beta*C where A is an m by m or n by n symmetric matrix and B is a general matrix. More... | |
template<typename alpha_type , typename A_type , typename beta_type > | |
void | SYRK (EUplo uplo, ETransp trans, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const |
Performs the symmetric rank k operation: C <- alpha*A*A'+beta*C or C <- alpha*A'*A+beta*C , where alpha and beta are scalars, C is an n by n symmetric matrix and A is an n by k matrix in the first case or k by n matrix in the second case. More... | |
template<typename alpha_type , typename A_type > | |
void | TRMM (ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const |
Performs the matrix-matrix operation: B <- alpha*op (A)*B or B <- alpha*B*op (A) where op(A) is an unit/non-unit, upper/lower triangular matrix and B is a general matrix. More... | |
template<typename alpha_type , typename A_type > | |
void | TRSM (ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const |
Solves the matrix equations: op(A)*X=alpha*B or X*op (A)=alpha*B where X and B are m by n matrices, A is a unit/non-unit, upper/lower triangular matrix and op(A) is A or A' . The matrix X is overwritten on B . More... | |
void | ROTG (ScalarType *da, ScalarType *db, rotg_c_type *c, ScalarType *s) const |
Computes a Givens plane rotation. More... | |
void | ROT (const OrdinalType &n, ScalarType *dx, const OrdinalType &incx, ScalarType *dy, const OrdinalType &incy, MagnitudeType *c, ScalarType *s) const |
Applies a Givens plane rotation. More... | |
void | SCAL (const OrdinalType &n, const ScalarType &alpha, ScalarType *x, const OrdinalType &incx) const |
Scale the vector x by the constant alpha . More... | |
void | COPY (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const |
Copy the vector x to the vector y . More... | |
template<typename alpha_type , typename x_type > | |
void | AXPY (const OrdinalType &n, const alpha_type alpha, const x_type *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const |
Perform the operation: y <- y+alpha*x . More... | |
ScalarTraits< ScalarType > ::magnitudeType | ASUM (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const |
Sum the absolute values of the entries of x . More... | |
template<typename x_type , typename y_type > | |
ScalarType | DOT (const OrdinalType &n, const x_type *x, const OrdinalType &incx, const y_type *y, const OrdinalType &incy) const |
Form the dot product of the vectors x and y . More... | |
ScalarTraits< ScalarType > ::magnitudeType | NRM2 (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const |
Compute the 2-norm of the vector x . More... | |
OrdinalType | IAMAX (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const |
Return the index of the element of x with the maximum magnitude. More... | |
Public Member Functions inherited from Teuchos::LAPACK< OrdinalType, ScalarType > | |
LAPACK (void) | |
Default Constructor. More... | |
LAPACK (const LAPACK< OrdinalType, ScalarType > &lapack) | |
Copy Constructor. More... | |
virtual | ~LAPACK (void) |
Destructor. More... | |
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... | |
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... | |
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... | |
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 |
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... | |
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... | |
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 |
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... | |
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... | |
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... | |
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... | |
ScalarType | LAPY2 (const ScalarType &x, const ScalarType &y) const |
Computes x^2 + y^2 safely, to avoid overflow. More... | |
Additional Inherited Members | |
Public Types inherited from Teuchos::DefaultBLASImpl< OrdinalType, ScalarType > | |
typedef details::GivensRotator < ScalarType >::c_type | rotg_c_type |
The type used for c in ROTG. More... | |
Static Public Member Functions inherited from Teuchos::Object | |
static void | setTracebackMode (int tracebackModeValue) |
Set the value of the Object error traceback report mode. More... | |
static int | getTracebackMode () |
Get the value of the Object error traceback report mode. More... | |
Related Functions inherited from Teuchos::Object | |
std::ostream & | operator<< (std::ostream &os, const Teuchos::Object &obj) |
Print the given Object to the given output stream. More... | |
A class for solving dense linear problems.
The Teuchos::SerialQRDenseSolver class enables the definition, in terms of Teuchos::SerialDenseMatrix and Teuchos::SerialDenseVector objects, of a dense linear problem, followed by the solution of that problem via the most sophisticated techniques available in LAPACK.
The Teuchos::SerialQRDenseSolver class is intended to provide full-featured support for solving linear problems for general dense rectangular (or square) matrices. It is written on top of BLAS and LAPACK and thus has excellent performance and numerical capabilities. Using this class, one can either perform simple factorizations and solves or apply all the tricks available in LAPACK to get the best possible solution for very ill-conditioned problems.
Teuchos::SerialQRDenseSolver vs. Teuchos::LAPACK
The Teuchos::LAPACK class provides access to most of the same functionality as Teuchos::SerialQRDenseSolver. The primary difference is that Teuchos::LAPACK is a "thin" layer on top of LAPACK and Teuchos::SerialQRDenseSolver attempts to provide easy access to the more sophisticated aspects of solving dense linear and eigensystems.
Constructing Teuchos::SerialQRDenseSolver Objects
There is a single Teuchos::SerialQRDenseSolver constructor. However, the matrix, right hand side and solution vectors must be set prior to executing most methods in this class.
Setting vectors used for linear solves
The matrix A, the left hand side X and the right hand side B (when solving AX = B, for X), can be set by appropriate set methods. Each of these three objects must be an Teuchos::SerialDenseMatrix or and Teuchos::SerialDenseVector object. The set methods are as follows:
Vector and Utility Functions
Once a Teuchos::SerialQRDenseSolver is constructed, several mathematical functions can be applied to the object. Specifically:
Strategies for Solving Linear Systems In many cases, linear least squares systems can be accurately solved by simply computing the QR factorization of the matrix and then performing a forward back solve with a given set of right hand side vectors. However, in some instances, the factorization may be very poorly conditioned and this simple approach may not work. In these situations, equilibration and iterative refinement may improve the accuracy, or prevent a breakdown in the factorization.
Teuchos::SerialQRDenseSolver will use equilibration with the factorization if, once the object is constructed and before it is factored, you call the function factorWithEquilibration(true) to force equilibration to be used. If you are uncertain if equilibration should be used, you may call the function shouldEquilibrate() which will return true if equilibration could possibly help. shouldEquilibrate() uses guidelines specified in the LAPACK User Guide to determine if equilibration might be useful.
Examples using Teuchos::SerialQRDenseSolver can be found in the Teuchos test directories.
Definition at line 132 of file Teuchos_SerialQRDenseSolver.hpp.
Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::SerialQRDenseSolver | ( | ) |
Default constructor; matrix should be set using setMatrix(), LHS and RHS set with setVectors().
Definition at line 400 of file Teuchos_SerialQRDenseSolver.hpp.
|
virtual |
SerialQRDenseSolver destructor.
Definition at line 436 of file Teuchos_SerialQRDenseSolver.hpp.
int Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::setMatrix | ( | const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > & | A | ) |
Sets the pointers for coefficient matrix.
Row dimension of A must be greater than or equal to the column dimension of A.
Definition at line 482 of file Teuchos_SerialQRDenseSolver.hpp.
int Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::setVectors | ( | const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > & | X, |
const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > & | B | ||
) |
Sets the pointers for left and right hand side vector(s).
Row dimension of X must match column dimension of matrix A, row dimension of B must match row dimension of A.
Definition at line 509 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Causes equilibration to be called just before the matrix factorization as part of the call to factor
.
Definition at line 186 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
If flag
is true, causes all subsequent function calls to work with the adjoint of this matrix, otherwise not.
Definition at line 189 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
All subsequent function calls will work with the transpose-type set by this method (Teuchos::NO_TRANS
or Teuchos::CONJ_TRANS).
Definition at line 192 of file Teuchos_SerialQRDenseSolver.hpp.
int Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::factor | ( | ) |
Computes the in-place QR factorization of the matrix using the LAPACK routine _GETRF or the Eigen class HouseholderQR.
Definition at line 533 of file Teuchos_SerialQRDenseSolver.hpp.
int Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::solve | ( | ) |
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()..
Definition at line 575 of file Teuchos_SerialQRDenseSolver.hpp.
int Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::computeEquilibrateScaling | ( | ) |
Determines if this matrix should be scaled.
Definition at line 661 of file Teuchos_SerialQRDenseSolver.hpp.
int Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::equilibrateMatrix | ( | ) |
Equilibrates the this matrix.
Definition at line 700 of file Teuchos_SerialQRDenseSolver.hpp.
int Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::equilibrateRHS | ( | ) |
Equilibrates the current RHS.
Definition at line 746 of file Teuchos_SerialQRDenseSolver.hpp.
int Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::unequilibrateLHS | ( | ) |
Unscales the solution vectors if equilibration was used to solve the system.
Definition at line 792 of file Teuchos_SerialQRDenseSolver.hpp.
int Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::formQ | ( | ) |
Explicitly forms the unitary matrix Q.
Definition at line 829 of file Teuchos_SerialQRDenseSolver.hpp.
int Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::formR | ( | ) |
Explicitly forms the upper triangular matrix R.
Definition at line 866 of file Teuchos_SerialQRDenseSolver.hpp.
int Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::multiplyQ | ( | ETransp | transq, |
SerialDenseMatrix< OrdinalType, ScalarType > & | C | ||
) |
Left multiply the input matrix by the unitary matrix Q or its adjoint.
Definition at line 893 of file Teuchos_SerialQRDenseSolver.hpp.
int Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::solveR | ( | ETransp | transr, |
SerialDenseMatrix< OrdinalType, ScalarType > & | C | ||
) |
Solve input matrix on the left with the upper triangular matrix R or its adjoint.
Definition at line 954 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns true if adjoint of this matrix has and will be used.
Definition at line 267 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns true if matrix is factored (factor available via getFactoredMatrix()).
Definition at line 270 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns true if factor is equilibrated (factor available via getFactoredMatrix()).
Definition at line 273 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns true if RHS is equilibrated (RHS available via getRHS()).
Definition at line 276 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system.
Definition at line 279 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns true if the current set of vectors has been solved.
Definition at line 282 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns true if Q has been formed explicitly.
Definition at line 285 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns true if R has been formed explicitly.
Definition at line 288 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns pointer to current matrix.
Definition at line 296 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns pointer to factored matrix (assuming factorization has been performed).
Definition at line 299 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns pointer to Q (assuming factorization has been performed).
Definition at line 302 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns pointer to R (assuming factorization has been performed).
Definition at line 305 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns pointer to current LHS.
Definition at line 308 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns pointer to current RHS.
Definition at line 311 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns row dimension of system.
Definition at line 314 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns column dimension of system.
Definition at line 317 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
Definition at line 320 of file Teuchos_SerialQRDenseSolver.hpp.
|
inline |
Returns the absolute value of the largest element of this matrix (returns -1 if not yet computed).
Definition at line 323 of file Teuchos_SerialQRDenseSolver.hpp.
void Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType >::Print | ( | std::ostream & | os | ) | const |
Print service methods; defines behavior of ostream << operator.
Definition at line 1016 of file Teuchos_SerialQRDenseSolver.hpp.