44 #ifndef IFPACK_SERIALTRIDISOLVER_H 
   45 #define IFPACK_SERIALTRIDISOLVER_H 
   50 #include "Epetra_Object.h" 
   51 #include "Epetra_CompObject.h" 
   52 #include "Epetra_BLAS.h" 
   53 #include "Teuchos_LAPACK.hpp" 
  168   void SolveWithTranspose(
bool Flag) {Transpose_ = Flag; 
if (Flag) TRANS_ = 
'T'; 
else TRANS_ = 
'N'; 
return;};
 
  177   void EstimateSolutionErrors(
bool Flag) ;
 
  187   virtual int Factor(
void);
 
  193   virtual int Solve(
void);
 
  199   virtual int Invert(
void);
 
  205   virtual int ApplyRefinement(
void);
 
  220   virtual int ReciprocalConditionEstimate(
double & Value);
 
  264   int N()
  const {
return(N_);};
 
  267   double * 
A()
  const {
return(A_);};
 
  270   int LDA()
  const {
return(LDA_);};
 
  273   double * 
B()
  const {
return(B_);};
 
  276   int LDB()
  const {
return(LDB_);};
 
  279   int NRHS()
  const {
return(NRHS_);};
 
  282   double * 
X()
  const {
return(X_);};
 
  285   int LDX()
  const {
return(LDX_);};
 
  288   double * 
AF()
  const {
return(AF_);};
 
  291   int LDAF()
  const {
return(LDAF_);};
 
  294   int * 
IPIV()
  const {
return(IPIV_);};
 
  297   double ANORM()
  const {
return(ANORM_);};
 
  300   double RCOND()
  const {
return(RCOND_);};
 
  305   double ROWCND()
  const {
return(ROWCND_);};
 
  310   double COLCND()
  const {
return(COLCND_);};
 
  313   double AMAX()
  const {
return(AMAX_);};
 
  316   double * 
FERR()
  const {
return(FERR_);};
 
  319   double * 
BERR()
  const {
return(BERR_);};
 
  325   virtual void Print(std::ostream& os) 
const;
 
  330   void AllocateWORK() {
if (WORK_==0) {LWORK_ = 4*N_; WORK_ = 
new double[LWORK_];} 
return;};
 
  331   void AllocateIWORK() {
if (IWORK_==0) IWORK_ = 
new int[N_]; 
return;};
 
  339   bool EstimateSolutionErrors_;
 
  340   bool SolutionErrorsEstimated_;
 
  343   bool ReciprocalConditionEstimated_;
 
  344   bool RefineSolution_;
 
  345   bool SolutionRefined_;
 
  384   Teuchos::LAPACK<int,double> lapack;
 
double COLCND() const 
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
 
int LDB() const 
Returns the leading dimension of the RHS. 
 
double AMAX() const 
Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed)...
 
bool Solved()
Returns true if the current set of vectors has been solved. 
 
bool SolutionRefined()
Returns true if the current set of vectors has been refined. 
 
void SolveWithTranspose(bool Flag)
Causes equilibration to be called just before the matrix factorization as part of the call to Factor...
 
void SolveToRefinedSolution(bool Flag)
Causes all solves to compute solution to best ability using iterative refinement. ...
 
bool SolutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
 
int N() const 
Returns column dimension of system. 
 
double ROWCND() const 
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed)...
 
Ifpack_SerialTriDiSolver: A class for solving TriDi linear problems. 
 
bool Inverted()
Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()). 
 
double * X() const 
Returns pointer to current solution. 
 
bool Factored()
Returns true if matrix is factored (factor available via AF() and LDAF()). 
 
Ifpack_SerialTriDiMatrix * FactoredMatrix() const 
Returns pointer to factored matrix (assuming factorization has been performed). 
 
Epetra_SerialDenseMatrix * RHS() const 
Returns pointer to current RHS. 
 
double * B() const 
Returns pointer to current RHS. 
 
double * A() const 
Returns pointer to the this matrix. 
 
double ANORM() const 
Returns the 1-Norm of the this matrix (returns -1 if not yet computed). 
 
double * AF() const 
Returns pointer to the factored matrix (may be the same as A() if factorization done in place)...
 
double RCOND() const 
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed)...
 
int LDA() const 
Returns the leading dimension of the this matrix. 
 
bool ReciprocalConditionEstimated()
Returns true if the condition number of the this matrix has been computed (value available via Recipr...
 
int LDAF() const 
Returns the leading dimension of the factored matrix. 
 
bool Transpose()
Returns true if transpose of this matrix has and will be used. 
 
Ifpack_SerialTriDiMatrix: A class for constructing and using real double precision general TriDi matr...
 
double * FERR() const 
Returns a pointer to the forward error estimates computed by LAPACK. 
 
int NRHS() const 
Returns the number of current right hand sides and solution vectors. 
 
Epetra_SerialDenseMatrix * LHS() const 
Returns pointer to current LHS. 
 
Ifpack_SerialTriDiMatrix * Matrix() const 
Returns pointer to current matrix. 
 
int LDX() const 
Returns the leading dimension of the solution. 
 
int * IPIV() const 
Returns pointer to pivot vector (if factorization has been computed), zero otherwise. 
 
double * BERR() const 
Returns a pointer to the backward error estimates computed by LAPACK.