44 #ifndef IFPACK_SERIALTRIDISOLVER_H 
   45 #define IFPACK_SERIALTRIDISOLVER_H 
   47 #if defined(Ifpack_SHOW_DEPRECATED_WARNINGS) 
   49 #warning "The Ifpack package is deprecated" 
   56 #include "Epetra_Object.h" 
   57 #include "Epetra_CompObject.h" 
   58 #include "Epetra_BLAS.h" 
   59 #include "Teuchos_LAPACK.hpp" 
  199   virtual int Solve(
void);
 
  270   int N()
  const {
return(
N_);};
 
  273   double * 
A()
  const {
return(
A_);};
 
  279   double * 
B()
  const {
return(
B_);};
 
  288   double * 
X()
  const {
return(
X_);};
 
  331   virtual void Print(std::ostream& os) 
const;
 
double COLCND() const 
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
virtual int Invert(void)
Inverts the this matrix. 
Epetra_SerialDenseMatrix * RHS_
Ifpack_SerialTriDiSolver()
Default constructor; matrix should be set using SetMatrix(), LHS and RHS set with SetVectors()...
bool ReciprocalConditionEstimated_
virtual int Solve(void)
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors().. 
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)...
virtual int ApplyRefinement(void)
Apply Iterative Refinement. 
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. 
virtual int Factor(void)
Computes the in-place LU factorization of the matrix using the LAPACK routine DGETRF. 
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. 
int SetMatrix(Ifpack_SerialTriDiMatrix &A)
Sets the pointers for coefficient matrix. 
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. 
virtual ~Ifpack_SerialTriDiSolver()
Ifpack_SerialTriDiSolver destructor. 
bool Inverted()
Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()). 
void EstimateSolutionErrors(bool Flag)
Causes all solves to estimate the forward and backward solution error. 
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)...
bool SolutionErrorsEstimated_
Ifpack_SerialTriDiMatrix * Matrix_
double RCOND() const 
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed)...
virtual int ReciprocalConditionEstimate(double &Value)
Unscales the solution vectors if equilibration was used to solve the system. 
Ifpack_SerialTriDiSolver & operator=(const Ifpack_SerialTriDiSolver &Source)
bool EstimateSolutionErrors_
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...
Teuchos::LAPACK< int, double > lapack
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. 
Ifpack_SerialTriDiMatrix * Factor_
int SetVectors(Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &B)
Sets the pointers for left and right hand side vector(s). 
Epetra_SerialDenseMatrix * LHS() const 
Returns pointer to current LHS. 
Ifpack_SerialTriDiMatrix * Matrix() const 
Returns pointer to current matrix. 
virtual void Print(std::ostream &os) const 
Print service methods; defines behavior of ostream << operator. 
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. 
Epetra_SerialDenseMatrix * LHS_
double * BERR() const 
Returns a pointer to the backward error estimates computed by LAPACK.