44 #ifndef EPETRA_SERIALDENSESOLVER_H
45 #define EPETRA_SERIALDENSESOLVER_H
47 #if defined(Epetra_SHOW_DEPRECATED_WARNINGS)
49 #warning "The Epetra package is deprecated"
177 void SolveWithTranspose(
bool Flag) {Transpose_ = Flag;
if (Flag) TRANS_ =
'T';
else TRANS_ =
'N';
return;};
186 void EstimateSolutionErrors(
bool Flag) ;
196 virtual int Factor(
void);
202 virtual int Solve(
void);
208 virtual int Invert(
void);
214 virtual int ComputeEquilibrateScaling(
void);
220 virtual int EquilibrateMatrix(
void);
226 int EquilibrateRHS(
void);
233 virtual int ApplyRefinement(
void);
239 int UnequilibrateLHS(
void);
248 virtual int ReciprocalConditionEstimate(
double & Value);
301 int M()
const {
return(M_);};
304 int N()
const {
return(N_);};
307 double *
A()
const {
return(A_);};
310 int LDA()
const {
return(LDA_);};
313 double *
B()
const {
return(B_);};
316 int LDB()
const {
return(LDB_);};
319 int NRHS()
const {
return(NRHS_);};
322 double *
X()
const {
return(X_);};
325 int LDX()
const {
return(LDX_);};
328 double *
AF()
const {
return(AF_);};
331 int LDAF()
const {
return(LDAF_);};
334 int *
IPIV()
const {
return(IPIV_);};
337 double ANORM()
const {
return(ANORM_);};
340 double RCOND()
const {
return(RCOND_);};
345 double ROWCND()
const {
return(ROWCND_);};
350 double COLCND()
const {
return(COLCND_);};
353 double AMAX()
const {
return(AMAX_);};
356 double *
FERR()
const {
return(FERR_);};
359 double *
BERR()
const {
return(BERR_);};
362 double *
R()
const {
return(R_);};
365 double *
C()
const {
return(C_);};
370 virtual void Print(std::ostream& os)
const;
375 void AllocateWORK() {
if (WORK_==0) {LWORK_ = 4*N_; WORK_ =
new double[LWORK_];}
return;};
bool Transpose()
Returns true if transpose of this matrix has and will be used.
double * B() const
Returns pointer to current RHS.
virtual void Print(std::ostream &os) const
Print object to an output stream Print method.
double ROWCND() const
Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed)...
bool B_Equilibrated()
Returns true if RHS is equilibrated (RHS available via B() and LDB()).
bool SolutionRefined()
Returns true if the current set of vectors has been refined.
double COLCND() const
Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed...
double * C() const
Returns a pointer to the column scale vector used for equilibration.
double * BERR() const
Returns a pointer to the backward error estimates computed by LAPACK.
int M() const
Returns row dimension of system.
void SolveWithTranspose(bool Flag)
If Flag is true, causes all subsequent function calls to work with the transpose of this matrix...
double * R() const
Returns a pointer to the row scaling vector used for equilibration.
bool SolutionErrorsEstimated()
Returns true if forward and backward error estimated have been computed (available via FERR() and BER...
Epetra_CompObject & operator=(const Epetra_CompObject &src)
Epetra_SerialDenseMatrix: A class for constructing and using real double precision general dense matr...
Epetra_SerialDenseMatrix * LHS() const
Returns pointer to current LHS.
double ANORM() const
Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
int * IPIV() const
Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
Epetra_SerialDenseMatrix * RHS_
double RCOND() const
Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed)...
Epetra_BLAS: The Epetra BLAS Wrapper Class.
int N() const
Returns column dimension of system.
Epetra_SerialDenseMatrix * RHS() const
Returns pointer to current RHS.
Epetra_Object: The base Epetra class.
virtual bool ShouldEquilibrate()
Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system...
Epetra_SerialDenseMatrix * Matrix_
Epetra_CompObject: Functionality and data that is common to all computational classes.
int LDB() const
Returns the leading dimension of the RHS.
bool Inverted()
Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()).
double * X() const
Returns pointer to current solution.
int LDA() const
Returns the leading dimension of the this matrix.
bool SolutionErrorsEstimated_
bool ReciprocalConditionEstimated()
Returns true if the condition number of the this matrix has been computed (value available via Recipr...
Epetra_LAPACK: The Epetra LAPACK Wrapper Class.
int LDAF() const
Returns the leading dimension of the factored matrix.
void SolveToRefinedSolution(bool Flag)
Causes all solves to compute solution to best ability using iterative refinement. ...
int LDX() const
Returns the leading dimension of the solution.
Epetra_SerialDenseMatrix * FactoredMatrix() const
Returns pointer to factored matrix (assuming factorization has been performed).
bool ReciprocalConditionEstimated_
bool A_Equilibrated()
Returns true if factor is equilibrated (factor available via AF() and LDAF()).
Epetra_SerialDenseMatrix * Factor_
void FactorWithEquilibration(bool Flag)
Causes equilibration to be called just before the matrix factorization as part of the call to Factor...
bool Factored()
Returns true if matrix is factored (factor available via AF() and LDAF()).
Epetra_SerialDenseMatrix * Matrix() const
Returns pointer to current matrix.
double AMAX() const
Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed)...
double * FERR() const
Returns a pointer to the forward error estimates computed by LAPACK.
bool Solved()
Returns true if the current set of vectors has been solved.
int Solve(int, TYPE *, TYPE *, TYPE *)
double * AF() const
Returns pointer to the factored matrix (may be the same as A() if factorization done in place)...
int NRHS() const
Returns the number of current right hand sides and solution vectors.
Epetra_SerialDenseMatrix * LHS_
Epetra_SerialDenseSolver: A class for solving dense linear problems.
bool EstimateSolutionErrors_
double * A() const
Returns pointer to the this matrix.