43 #ifndef IFPACK_AMESOS_H
44 #define IFPACK_AMESOS_H
48 #include "Epetra_Operator.h"
49 #include "Teuchos_ParameterList.hpp"
50 #include "Teuchos_RefCountPtr.hpp"
145 virtual double NormInf()
const;
151 virtual const char *
Label()
const;
218 const int MaxIters = 1550,
219 const double Tol = 1e-9,
289 virtual std::ostream&
Print(std::ostream& os)
const;
373 Teuchos::RefCountPtr<const Epetra_RowMatrix>
Matrix_;
376 Teuchos::RefCountPtr<Epetra_LinearProblem>
Problem_;
378 Teuchos::RefCountPtr<Amesos_BaseSolver>
Solver_;
407 Teuchos::RefCountPtr<Epetra_Time>
Time_;
418 #endif // IFPACK_AMESOS_H
virtual const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
void SetComputeFlops(const double ComputeFlops_in)
Sets ComputeFlops_.
void SetComputeTime(const double ComputeTime_in)
Sets ComputeTime_.
void SetNumInitialize(const int NumInitialize_in)
Sets NumInitialize_.
virtual double ComputeTime() const
Returns the total time spent in Compute().
virtual int NumInitialize() const
Returns the number of calls to Initialize().
double Condest_
Contains the estimated condition number.
virtual double ApplyInverseTime() const
Returns the total time spent in ApplyInverse().
virtual int SetParameters(Teuchos::ParameterList &List)
Sets all the parameters for the preconditioner.
virtual double ComputeFlops() const
Returns the total number of flops to computate the preconditioner.
bool UseTranspose_
If true, the preconditioner solves for the transpose of the matrix.
double ComputeTime_
Contains the time for all successful calls to Compute().
virtual int NumCompute() const
Returns the number of calls to Compute().
int NumApplyInverse_
Contains the number of successful call to ApplyInverse().
std::string Label_
Contains the label of this object.
int NumInitialize_
Contains the number of successful calls to Initialize().
Ifpack_Amesos(Epetra_RowMatrix *Matrix)
Constructor.
virtual double ApplyInverseFlops() const
Returns the total number of flops to apply the preconditioner.
Teuchos::ParameterList List_
Contains a copy of the input parameter list.
virtual int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Applies the matrix to an Epetra_MultiVector.
void SetIsInitialized(const bool IsInitialized_in)
Sets IsInitialized_.
virtual const Epetra_RowMatrix & Matrix() const
Returns a const reference to the internally stored matrix.
void SetList(const Teuchos::ParameterList &List_in)
Set List_.
Teuchos::RefCountPtr< Amesos_BaseSolver > Solver_
Amesos solver, use to apply the inverse of the local matrix.
double InitializeTime_
Contains the time for all successful calls to Initialize().
virtual double InitializeFlops() const
Returns the number of flops in the initialization phase.
virtual bool IsComputed() const
Returns true if the preconditioner has been successfully computed.
virtual const Teuchos::ParameterList & List() const
Teuchos::RefCountPtr< const Epetra_RowMatrix > Matrix_
Pointers to the matrix to be preconditioned.
Ifpack_CondestType
Ifpack_CondestType: enum to define the type of condition number estimate.
Ifpack_Amesos: a class to use Amesos' factorizations as preconditioners.
virtual int Compute()
Computes the preconditioners.
double ComputeFlops_
Contains the number of flops for Compute().
void SetApplyInverseTime(const double ApplyInverseTime_in)
Sets ApplyInverseTime_.
virtual int Initialize()
Initializes the preconditioners.
void SetIsComputed(const int IsComputed_in)
Sets IsComputed_.
virtual int SetUseTranspose(bool UseTranspose_in)
If set true, transpose of this operator will be applied (not implemented).
virtual const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
Ifpack_Preconditioner: basic class for preconditioning in Ifpack.
virtual const char * Label() const
Returns a character string describing the operator.
virtual double Condest() const
Returns the estimated condition number, never computes it.
Teuchos::RefCountPtr< Epetra_LinearProblem > Problem_
Linear problem required by Solver_.
virtual int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Applies the preconditioner to X, returns the result in Y.
virtual double InitializeTime() const
Returns the total time spent in Initialize().
virtual int NumApplyInverse() const
Returns the number of calls to ApplyInverse().
virtual bool UseTranspose() const
Returns the current UseTranspose setting.
bool IsComputed_
If true, the preconditioner has been successfully computed.
void SetNumApplyInverse(const int NumApplyInverse_in)
Sets NumApplyInverse_.
bool IsInitialized_
If true, the preconditioner has been successfully initialized.
void SetNumCompute(const int NumCompute_in)
Sets NumCompute_.
virtual bool HasNormInf() const
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
void SetLabel(const char *Label_in)
Sets the label.
virtual ~Ifpack_Amesos()
Destructor.
int NumCompute_
Contains the number of successful call to Compute().
Ifpack_Amesos & operator=(const Ifpack_Amesos &rhs)
Operator=.
virtual bool IsInitialized() const
Returns true is the preconditioner has been successfully initialized.
double ApplyInverseTime_
Contains the time for all successful calls to ApplyInverse().
bool IsEmpty_
If true, the linear system on this processor is empty, thus the preconditioner is null operation...
virtual double NormInf() const
Returns the infinity norm of the global matrix (not implemented)
void SetInitializeTime(const double InitializeTime_in)
Sets InitializeTime_.
double ApplyInverseFlops_
Contain sthe number of flops for ApplyInverse().
virtual const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this operator.
virtual std::ostream & Print(std::ostream &os) const
Prints on ostream basic information about this object.
Teuchos::RefCountPtr< Epetra_Time > Time_
Time object.
void SetApplyInverseFlops(const double ApplyInverseFlops_in)
Sets ComputeFlops_.