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"
174 void SolveWithTranspose(
bool Flag) {Transpose_ = Flag;
if (Flag) TRANS_ =
'T';
else TRANS_ =
'N';
return;};
199 virtual int Solve(
void);
270 int N()
const {
return(N_);};
273 double *
A()
const {
return(A_);};
276 int LDA()
const {
return(LDA_);};
279 double *
B()
const {
return(B_);};
282 int LDB()
const {
return(LDB_);};
285 int NRHS()
const {
return(NRHS_);};
288 double *
X()
const {
return(X_);};
291 int LDX()
const {
return(LDX_);};
294 double *
AF()
const {
return(AF_);};
297 int LDAF()
const {
return(LDAF_);};
300 int *
IPIV()
const {
return(IPIV_);};
303 double ANORM()
const {
return(ANORM_);};
306 double RCOND()
const {
return(RCOND_);};
311 double ROWCND()
const {
return(ROWCND_);};
316 double COLCND()
const {
return(COLCND_);};
319 double AMAX()
const {
return(AMAX_);};
322 double *
FERR()
const {
return(FERR_);};
325 double *
BERR()
const {
return(BERR_);};
332 virtual void Print(std::ostream& os)
const;
336 void AllocateWORK() {
if (WORK_==0) {LWORK_ = 4*N_; WORK_ =
new double[LWORK_];}
return;};
337 void AllocateIWORK() {
if (IWORK_==0) IWORK_ =
new int[N_];
return;};
345 bool EstimateSolutionErrors_;
346 bool SolutionErrorsEstimated_;
349 bool ReciprocalConditionEstimated_;
350 bool RefineSolution_;
351 bool SolutionRefined_;
390 Teuchos::LAPACK<int,double> lapack;
Ifpack_SerialTriDiMatrix: A class for constructing and using real double precision general TriDi matr...
Ifpack_SerialTriDiSolver: A class for solving TriDi linear problems.
bool ReciprocalConditionEstimated()
Returns true if the condition number of the this matrix has been computed (value available via Recipr...
double * BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.
bool Inverted()
Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()).
bool Solved()
Returns true if the current set of vectors has been solved.
void EstimateSolutionErrors(bool Flag)
Causes all solves to estimate the forward and backward solution error.
int N() const
Returns column dimension of system.
int LDX() const
Returns the leading dimension of the solution.
virtual int ApplyRefinement(void)
Apply Iterative Refinement.
int NRHS() const
Returns the number of current right hand sides and solution vectors.
double * X() const
Returns pointer to current solution.
double ROWCND() const
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed).
virtual void Print(std::ostream &os) const
Print service methods; defines behavior of ostream << operator.
double * A() const
Returns pointer to the this matrix.
bool Factored()
Returns true if matrix is factored (factor available via AF() and LDAF()).
virtual int Invert(void)
Inverts the this matrix.
int SetMatrix(Ifpack_SerialTriDiMatrix &A)
Sets the pointers for coefficient matrix.
double AMAX() const
Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed).
virtual int Solve(void)
Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()....
bool Transpose()
Returns true if transpose of this matrix has and will be used.
Epetra_SerialDenseMatrix * LHS() const
Returns pointer to current LHS.
Ifpack_SerialTriDiMatrix * FactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
virtual int Factor(void)
Computes the in-place LU factorization of the matrix using the LAPACK routine DGETRF.
int LDAF() const
Returns the leading dimension of the factored matrix.
double * FERR() const
Returns a pointer to the forward error estimates computed by LAPACK.
virtual int ReciprocalConditionEstimate(double &Value)
Unscales the solution vectors if equilibration was used to solve the system.
int * IPIV() const
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
int LDB() const
Returns the leading dimension of the RHS.
Epetra_SerialDenseMatrix * RHS() const
Returns pointer to current RHS.
double COLCND() const
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
int SetVectors(Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &B)
Sets the pointers for left and right hand side vector(s).
Ifpack_SerialTriDiSolver()
Default constructor; matrix should be set using SetMatrix(), LHS and RHS set with SetVectors().
virtual ~Ifpack_SerialTriDiSolver()
Ifpack_SerialTriDiSolver destructor.
void SolveToRefinedSolution(bool Flag)
Causes all solves to compute solution to best ability using iterative refinement.
double * AF() const
Returns pointer to the factored matrix (may be the same as A() if factorization done in place).
void SolveWithTranspose(bool Flag)
Causes equilibration to be called just before the matrix factorization as part of the call to Factor.
double RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed).
bool SolutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
double * B() const
Returns pointer to current RHS.
int LDA() const
Returns the leading dimension of the this matrix.
bool SolutionRefined()
Returns true if the current set of vectors has been refined.
double ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
Ifpack_SerialTriDiMatrix * Matrix() const
Returns pointer to current matrix.