Teuchos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType > Class Template Reference

A class for representing and solving banded dense linear systems. More...

#include <Teuchos_SerialBandDenseSolver.hpp>

Inheritance diagram for Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >:
Inheritance graph
[legend]

Public Types

typedef Teuchos::ScalarTraits
< ScalarType >::magnitudeType 
MagnitudeType
 
- 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...
 
- Public Types inherited from Teuchos::LAPACK< OrdinalType, ScalarType >
typedef Teuchos::ScalarTraits
< ScalarType >::magnitudeType 
MagnitudeType
 

Protected Member Functions

void allocateWORK ()
 
void resetMatrix ()
 
void resetVectors ()
 

Protected Attributes

bool equilibrate_
 
bool shouldEquilibrate_
 
bool equilibratedA_
 
bool equilibratedB_
 
bool transpose_
 
bool factored_
 
bool estimateSolutionErrors_
 
bool solutionErrorsEstimated_
 
bool solved_
 
bool reciprocalConditionEstimated_
 
bool refineSolution_
 
bool solutionRefined_
 
Teuchos::ETransp TRANS_
 
OrdinalType M_
 
OrdinalType N_
 
OrdinalType KL_
 
OrdinalType KU_
 
OrdinalType Min_MN_
 
OrdinalType LDA_
 
OrdinalType LDAF_
 
OrdinalType INFO_
 
OrdinalType LWORK_
 
std::vector< OrdinalType > IPIV_
 
std::vector< int > IWORK_
 
MagnitudeType ANORM_
 
MagnitudeType RCOND_
 
MagnitudeType ROWCND_
 
MagnitudeType COLCND_
 
MagnitudeType AMAX_
 
RCP< SerialBandDenseMatrix
< OrdinalType, ScalarType > > 
Matrix_
 
RCP< SerialDenseMatrix
< OrdinalType, ScalarType > > 
LHS_
 
RCP< SerialDenseMatrix
< OrdinalType, ScalarType > > 
RHS_
 
RCP< SerialBandDenseMatrix
< OrdinalType, ScalarType > > 
Factor_
 
ScalarType * A_
 
ScalarType * AF_
 
std::vector< MagnitudeTypeFERR_
 
std::vector< MagnitudeTypeBERR_
 
std::vector< ScalarType > WORK_
 
std::vector< MagnitudeTypeR_
 
std::vector< MagnitudeTypeC_
 
- Protected Attributes inherited from Teuchos::CompObject
FlopsflopCounter_
 

Private Member Functions

 SerialBandDenseSolver (const SerialBandDenseSolver< OrdinalType, ScalarType > &Source)
 
SerialBandDenseSolveroperator= (const SerialBandDenseSolver< OrdinalType, ScalarType > &Source)
 

Constructor/Destructor Methods

 SerialBandDenseSolver ()
 Default constructor; matrix should be set using setMatrix(), LHS and RHS set with setVectors(). More...
 
virtual ~SerialBandDenseSolver ()
 SerialBandDenseSolver destructor. More...
 

Set Methods

int setMatrix (const RCP< SerialBandDenseMatrix< OrdinalType, ScalarType > > &AB)
 Sets the pointer 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)
 
void solveWithTranspose (bool flag)
 
void solveWithTransposeFlag (Teuchos::ETransp trans)
 
void solveToRefinedSolution (bool flag)
 Set whether or not to use iterative refinement to improve solutions to linear systems. More...
 
void estimateSolutionErrors (bool flag)
 Causes all solves to estimate the forward and backward solution error. More...
 

Factor/Solve Methods

int factor ()
 Computes the in-place LU factorization of the matrix using the LAPACK routine _GBTRF. More...
 
int solve ()
 Computes the solution X to AX = B for the this matrix and the B provided to SetVectors().. More...
 
int computeEquilibrateScaling ()
 Computes the scaling vector S(i) = 1/sqrt(A(i,i)) of the this matrix. More...
 
int equilibrateMatrix ()
 Equilibrates the this matrix. More...
 
int equilibrateRHS ()
 Equilibrates the current RHS. More...
 
int applyRefinement ()
 Apply Iterative Refinement. More...
 
int unequilibrateLHS ()
 Unscales the solution vectors if equilibration was used to solve the system. More...
 
int reciprocalConditionEstimate (MagnitudeType &Value)
 Returns the reciprocal of the 1-norm condition number of the this matrix. More...
 

Query methods

bool transpose ()
 Returns true if transpose of this matrix has and will be used. More...
 
bool factored ()
 Returns true if matrix is factored (factor available via AF() and LDAF()). More...
 
bool equilibratedA ()
 Returns true if factor is equilibrated (factor available via AF() and LDAF()). More...
 
bool equilibratedB ()
 Returns true if RHS is equilibrated (RHS available via B() and LDB()). More...
 
bool shouldEquilibrate ()
 Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system. More...
 
bool solutionErrorsEstimated ()
 Returns true if forward and backward error estimated have been computed (available via FERR() and BERR()). More...
 
bool reciprocalConditionEstimated ()
 Returns true if the condition number of the this matrix has been computed (value available via ReciprocalConditionEstimate()). More...
 
bool solved ()
 Returns true if the current set of vectors has been solved. More...
 
bool solutionRefined ()
 Returns true if the current set of vectors has been refined. More...
 

Data Accessor methods

RCP< SerialBandDenseMatrix
< OrdinalType, ScalarType > > 
getMatrix () const
 Returns pointer to current matrix. More...
 
RCP< SerialBandDenseMatrix
< OrdinalType, ScalarType > > 
getFactoredMatrix () const
 Returns pointer to factored matrix (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...
 
OrdinalType numLower () const
 Returns lower bandwidth of system. More...
 
OrdinalType numUpper () const
 Returns upper bandwidth of system. More...
 
std::vector< OrdinalType > IPIV () const
 Returns pointer to pivot vector (if factorization has been computed), zero otherwise. More...
 
MagnitudeType ANORM () const
 Returns the 1-Norm of the this matrix (returns -1 if not yet computed). More...
 
MagnitudeType RCOND () const
 Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed). More...
 
MagnitudeType ROWCND () const
 Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed). More...
 
MagnitudeType COLCND () const
 Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed). More...
 
MagnitudeType AMAX () const
 Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed). More...
 
std::vector< MagnitudeTypeFERR () const
 Returns a pointer to the forward error estimates computed by LAPACK. More...
 
std::vector< MagnitudeTypeBERR () const
 Returns a pointer to the backward error estimates computed by LAPACK. More...
 
std::vector< MagnitudeTypeR () const
 Returns a pointer to the row scaling vector used for equilibration. More...
 
std::vector< MagnitudeTypeC () const
 Returns a pointer to the column scale vector used for equilibration. More...
 

I/O methods

void Print (std::ostream &os) const
 Print service methods; defines behavior of ostream << operator. More...
 

Additional Inherited Members

- 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...
 
FlopsgetFlopCounter () 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 setLabel (const char *theLabel)
 
virtual const char * label () const
 Access the object's label (LEGACY; return std::string instead). 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...
 
- 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...
 
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...
 
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...
 
- 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, MagnitudeType *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 MagnitudeType *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, const 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, const 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, 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, const 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, const 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, const ScalarType *S, const OrdinalType &lds, const 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...
 
- 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...
 
- Static Public Attributes inherited from Teuchos::Object
static int tracebackMode = -1
 

Detailed Description

template<typename OrdinalType, typename ScalarType>
class Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >

A class for representing and solving banded dense linear systems.

Template Parameters
OrdinalTypeThe index type used by the linear algebra implementation. This should always be int.
ScalarTypeThe type of entries in the matrix.

Introduction

This class defines a banded dense matrix, which may have any number of rows or columns (not necessarily equal). It's called "serial" because the matrix lives in a single memory space. Thus, it's the kind of matrix that one might give to the BLAS or LAPACK, not a distributed matrix like one would give to ScaLAPACK.

This class also has methods for computing the (banded) LU factorization of the matrix, and solving linear systems with the matrix. We use instances of SerialDenseVector to represent the right-hand side b or the solution x in the linear system $Ax=b$. The instance of this class used store the banded matrix must contain KL extra superdiagonals to store the L and U factors (see details below).

Users have the option to do equilibration before factoring the matrix. This may improve accuracy when solving ill-conditioned problems.

SerialBandDenseSolver and LAPACK

Teuchos' LAPACK class wraps LAPACK's LU factorizations, including the banded factorizations. It thus gives access to most of the same functionality as SerialBandDenseSolver. The main difference is that SerialBandDenseSolver offers a higher level of abstraction. It hides the details of which LAPACK routines to call. Furthermore, if you have built Teuchos with support for the third-party library Eigen, SerialBandDenseSolver lets you solve linear systems for ScalarType other than the four supported by the LAPACK library.

Constructing SerialBandDenseSolver objects

There is a single Teuchos::SerialBandDenseSolver 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 a SerialDenseMatrix or a SerialDenseVector object. The set methods are as follows:

Format of the matrix A

The SerialBandDenseMatrix must contain KL extra superdiagonals to store the L and U factors, where KL is the upper bandwidth. Consider using the non-member conversion routines generalToBanded and bandedToGeneral if the full SerialDenseMatrix is already in storage. However, it is more efficient simply to construct the SerialBandDenseMatrix with the desired parameters and use the provided matrix access operators so that the full rectangular matrix need not be stored. The conversion routine generalToBanded has a flag to store the input Teuchos::SerialDenseMatrix in banded format with KL extra superdiagonals so this class can use it. Again, it is more efficient to simply construct a Teuchos::SerialBandDenseMatrix object with KL extra superdiagonals than are needed for the matrix data and fill the matrix using the matrix access operators.

See the documentation of Teuchos::SerialBandDenseMatrix for further details on the storage format.

Vector and Utility Functions

Once a Teuchos::SerialBandDenseSolver is constructed, several mathematical functions can be applied to the object. Specifically:

Strategies for Solving Linear Systems

In many cases, linear systems can be accurately solved by simply computing the LU 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.

SerialBandDenseSolver 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, namely if SCOND < 0.1 and AMAX < Underflow or AMAX > Overflow, to determine if equilibration might be useful.

SerialBandDenseSolver will use iterative refinement after a forward/back solve if you call solveToRefinedSolution(true). It will also compute forward and backward error estimates if you call estimateSolutionErrors(true). Access to the forward (back) error estimates is available via FERR() (BERR()).

Examples using SerialBandDenseSolver can be found in the Teuchos test directories.

Definition at line 166 of file Teuchos_SerialBandDenseSolver.hpp.

Member Typedef Documentation

template<typename OrdinalType, typename ScalarType>
typedef Teuchos::ScalarTraits<ScalarType>::magnitudeType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::MagnitudeType

Definition at line 172 of file Teuchos_SerialBandDenseSolver.hpp.

Constructor & Destructor Documentation

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::SerialBandDenseSolver ( )

Default constructor; matrix should be set using setMatrix(), LHS and RHS set with setVectors().

Definition at line 464 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::~SerialBandDenseSolver ( )
virtual

SerialBandDenseSolver destructor.

Definition at line 501 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::SerialBandDenseSolver ( const SerialBandDenseSolver< OrdinalType, ScalarType > &  Source)
private

Member Function Documentation

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::setMatrix ( const RCP< SerialBandDenseMatrix< OrdinalType, ScalarType > > &  AB)

Sets the pointer for coefficient matrix.

Definition at line 546 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialBandDenseSolver< 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. X and B must have the same dimensions.

Definition at line 576 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
void Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::factorWithEquilibration ( bool  flag)
inline

Set whether or not to equilibrate just before the matrix factorization. This function must be called before the factorization is performed.

Definition at line 219 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
void Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::solveWithTranspose ( bool  flag)
inline

If flag is true, causes all subsequent function calls to work with the transpose of this matrix, otherwise not.

Note
This interface will not work correctly for complex-valued linear systems, use solveWithTransposeFlag().

Definition at line 224 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
void Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::solveWithTransposeFlag ( Teuchos::ETransp  trans)
inline

All subsequent function calls will work with the transpose-type set by this method (Teuchos::NO_TRANS, Teuchos::TRANS, and Teuchos::CONJ_TRANS).

Note
This interface will allow correct behavior for complex-valued linear systems, solveWithTranspose() will not.

Definition at line 228 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
void Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::solveToRefinedSolution ( bool  flag)
inline

Set whether or not to use iterative refinement to improve solutions to linear systems.

Definition at line 231 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
void Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::estimateSolutionErrors ( bool  flag)

Causes all solves to estimate the forward and backward solution error.

Error estimates will be in the arrays FERR and BERR, resp, after the solve step is complete. These arrays are accessible via the FERR() and BERR() access functions.

Definition at line 599 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::factor ( )

Computes the in-place LU factorization of the matrix using the LAPACK routine _GBTRF.

Returns
Integer error code, set to 0 if successful.

Definition at line 609 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::solve ( )

Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()..

Returns
Integer error code, set to 0 if successful.

Definition at line 667 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::computeEquilibrateScaling ( )

Computes the scaling vector S(i) = 1/sqrt(A(i,i)) of the this matrix.

Returns
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Definition at line 767 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::equilibrateMatrix ( )

Equilibrates the this matrix.

Returns
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Definition at line 788 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::equilibrateRHS ( )

Equilibrates the current RHS.

Returns
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Definition at line 844 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::applyRefinement ( )

Apply Iterative Refinement.

Returns
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Definition at line 733 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::unequilibrateLHS ( )

Unscales the solution vectors if equilibration was used to solve the system.

Returns
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Definition at line 876 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
int Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::reciprocalConditionEstimate ( MagnitudeType Value)

Returns the reciprocal of the 1-norm condition number of the this matrix.

Parameters
ValueOut On return contains the reciprocal of the 1-norm condition number of the this matrix.
Returns
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Definition at line 902 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::transpose ( )
inline

Returns true if transpose of this matrix has and will be used.

Definition at line 299 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::factored ( )
inline

Returns true if matrix is factored (factor available via AF() and LDAF()).

Definition at line 302 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::equilibratedA ( )
inline

Returns true if factor is equilibrated (factor available via AF() and LDAF()).

Definition at line 305 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::equilibratedB ( )
inline

Returns true if RHS is equilibrated (RHS available via B() and LDB()).

Definition at line 308 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::shouldEquilibrate ( )
inline

Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system.

Definition at line 311 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::solutionErrorsEstimated ( )
inline

Returns true if forward and backward error estimated have been computed (available via FERR() and BERR()).

Definition at line 314 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::reciprocalConditionEstimated ( )
inline

Returns true if the condition number of the this matrix has been computed (value available via ReciprocalConditionEstimate()).

Definition at line 317 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::solved ( )
inline

Returns true if the current set of vectors has been solved.

Definition at line 320 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::solutionRefined ( )
inline

Returns true if the current set of vectors has been refined.

Definition at line 323 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
RCP<SerialBandDenseMatrix<OrdinalType, ScalarType> > Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::getMatrix ( ) const
inline

Returns pointer to current matrix.

Definition at line 330 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
RCP<SerialBandDenseMatrix<OrdinalType, ScalarType> > Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::getFactoredMatrix ( ) const
inline

Returns pointer to factored matrix (assuming factorization has been performed).

Definition at line 333 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
RCP<SerialDenseMatrix<OrdinalType, ScalarType> > Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::getLHS ( ) const
inline

Returns pointer to current LHS.

Definition at line 336 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
RCP<SerialDenseMatrix<OrdinalType, ScalarType> > Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::getRHS ( ) const
inline

Returns pointer to current RHS.

Definition at line 339 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::numRows ( ) const
inline

Returns row dimension of system.

Definition at line 342 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::numCols ( ) const
inline

Returns column dimension of system.

Definition at line 345 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::numLower ( ) const
inline

Returns lower bandwidth of system.

Definition at line 348 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::numUpper ( ) const
inline

Returns upper bandwidth of system.

Definition at line 351 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
std::vector<OrdinalType> Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::IPIV ( ) const
inline

Returns pointer to pivot vector (if factorization has been computed), zero otherwise.

Definition at line 354 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
MagnitudeType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::ANORM ( ) const
inline

Returns the 1-Norm of the this matrix (returns -1 if not yet computed).

Definition at line 357 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
MagnitudeType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::RCOND ( ) const
inline

Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed).

Definition at line 360 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
MagnitudeType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::ROWCND ( ) const
inline

Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed).

If ROWCND() is >= 0.1 and AMAX() is not close to overflow or underflow, then equilibration is not needed.

Definition at line 365 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
MagnitudeType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::COLCND ( ) const
inline

Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed).

If COLCND() is >= 0.1 then equilibration is not needed.

Definition at line 370 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
MagnitudeType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::AMAX ( ) const
inline

Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed).

Definition at line 373 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
std::vector<MagnitudeType> Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::FERR ( ) const
inline

Returns a pointer to the forward error estimates computed by LAPACK.

Definition at line 376 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
std::vector<MagnitudeType> Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::BERR ( ) const
inline

Returns a pointer to the backward error estimates computed by LAPACK.

Definition at line 379 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
std::vector<MagnitudeType> Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::R ( ) const
inline

Returns a pointer to the row scaling vector used for equilibration.

Definition at line 382 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
std::vector<MagnitudeType> Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::C ( ) const
inline

Returns a pointer to the column scale vector used for equilibration.

Definition at line 385 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
void Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::Print ( std::ostream &  os) const

Print service methods; defines behavior of ostream << operator.

Definition at line 935 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
void Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::allocateWORK ( )
inlineprotected

Definition at line 395 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
void Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::resetMatrix ( )
protected

Definition at line 520 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType , typename ScalarType >
void Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::resetVectors ( )
protected

Definition at line 507 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
SerialBandDenseSolver& Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::operator= ( const SerialBandDenseSolver< OrdinalType, ScalarType > &  Source)
private

Member Data Documentation

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::equilibrate_
protected

Definition at line 400 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::shouldEquilibrate_
protected

Definition at line 401 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::equilibratedA_
protected

Definition at line 402 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::equilibratedB_
protected

Definition at line 403 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::transpose_
protected

Definition at line 404 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::factored_
protected

Definition at line 405 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::estimateSolutionErrors_
protected

Definition at line 406 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::solutionErrorsEstimated_
protected

Definition at line 407 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::solved_
protected

Definition at line 408 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::reciprocalConditionEstimated_
protected

Definition at line 409 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::refineSolution_
protected

Definition at line 410 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
bool Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::solutionRefined_
protected

Definition at line 411 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
Teuchos::ETransp Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::TRANS_
protected

Definition at line 413 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::M_
protected

Definition at line 415 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::N_
protected

Definition at line 416 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::KL_
protected

Definition at line 417 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::KU_
protected

Definition at line 418 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::Min_MN_
protected

Definition at line 419 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::LDA_
protected

Definition at line 420 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::LDAF_
protected

Definition at line 421 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::INFO_
protected

Definition at line 422 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::LWORK_
protected

Definition at line 423 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
std::vector<OrdinalType> Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::IPIV_
protected

Definition at line 425 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
std::vector<int> Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::IWORK_
protected

Definition at line 426 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
MagnitudeType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::ANORM_
protected

Definition at line 428 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
MagnitudeType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::RCOND_
protected

Definition at line 429 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
MagnitudeType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::ROWCND_
protected

Definition at line 430 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
MagnitudeType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::COLCND_
protected

Definition at line 431 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
MagnitudeType Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::AMAX_
protected

Definition at line 432 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
RCP<SerialBandDenseMatrix<OrdinalType, ScalarType> > Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::Matrix_
protected

Definition at line 434 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
RCP<SerialDenseMatrix<OrdinalType, ScalarType> > Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::LHS_
protected

Definition at line 435 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
RCP<SerialDenseMatrix<OrdinalType, ScalarType> > Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::RHS_
protected

Definition at line 436 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
RCP<SerialBandDenseMatrix<OrdinalType, ScalarType> > Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::Factor_
protected

Definition at line 437 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
ScalarType* Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::A_
protected

Definition at line 439 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
ScalarType* Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::AF_
protected

Definition at line 440 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
std::vector<MagnitudeType> Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::FERR_
protected

Definition at line 441 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
std::vector<MagnitudeType> Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::BERR_
protected

Definition at line 442 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
std::vector<ScalarType> Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::WORK_
protected

Definition at line 443 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
std::vector<MagnitudeType> Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::R_
protected

Definition at line 444 of file Teuchos_SerialBandDenseSolver.hpp.

template<typename OrdinalType, typename ScalarType>
std::vector<MagnitudeType> Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType >::C_
protected

Definition at line 445 of file Teuchos_SerialBandDenseSolver.hpp.


The documentation for this class was generated from the following file: