Ifpack_Amesos: a class to use Amesos' factorizations as preconditioners. More...
#include <Ifpack_Amesos.h>
Public Member Functions | |
Ifpack_Amesos (Epetra_RowMatrix *Matrix) | |
Constructor. | |
Ifpack_Amesos (const Ifpack_Amesos &rhs) | |
Copy constructor. | |
Ifpack_Amesos & | operator= (const Ifpack_Amesos &rhs) |
Operator=. | |
virtual | ~Ifpack_Amesos () |
Destructor. | |
virtual int | SetUseTranspose (bool UseTranspose_in) |
If set true, transpose of this operator will be applied (not implemented). More... | |
virtual int | Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
Applies the matrix to an Epetra_MultiVector. More... | |
virtual int | ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
Applies the preconditioner to X, returns the result in Y. More... | |
virtual double | NormInf () const |
Returns the infinity norm of the global matrix (not implemented) | |
virtual const char * | Label () const |
Returns a character string describing the operator. | |
virtual bool | UseTranspose () const |
Returns the current UseTranspose setting. | |
virtual bool | HasNormInf () const |
Returns true if the this object can provide an approximate Inf-norm, false otherwise. | |
virtual const Epetra_Comm & | Comm () const |
Returns a pointer to the Epetra_Comm communicator associated with this operator. | |
virtual const Epetra_Map & | OperatorDomainMap () const |
Returns the Epetra_Map object associated with the domain of this operator. | |
virtual const Epetra_Map & | OperatorRangeMap () const |
Returns the Epetra_Map object associated with the range of this operator. | |
virtual bool | IsInitialized () const |
Returns true is the preconditioner has been successfully initialized. | |
virtual int | Initialize () |
Initializes the preconditioners. More... | |
virtual bool | IsComputed () const |
Returns true if the preconditioner has been successfully computed. | |
virtual int | Compute () |
Computes the preconditioners. More... | |
virtual int | SetParameters (Teuchos::ParameterList &List) |
Sets all the parameters for the preconditioner. More... | |
virtual const Epetra_RowMatrix & | Matrix () const |
Returns a const reference to the internally stored matrix. | |
virtual double | Condest (const Ifpack_CondestType CT=Ifpack_Cheap, const int MaxIters=1550, const double Tol=1e-9, Epetra_RowMatrix *Matrix_in=0) |
Returns the estimated condition number, computes it if necessary. | |
virtual double | Condest () const |
Returns the estimated condition number, never computes it. | |
virtual int | NumInitialize () const |
Returns the number of calls to Initialize(). | |
virtual int | NumCompute () const |
Returns the number of calls to Compute(). | |
virtual int | NumApplyInverse () const |
Returns the number of calls to ApplyInverse(). | |
virtual double | InitializeTime () const |
Returns the total time spent in Initialize(). | |
virtual double | ComputeTime () const |
Returns the total time spent in Compute(). | |
virtual double | ApplyInverseTime () const |
Returns the total time spent in ApplyInverse(). | |
virtual double | InitializeFlops () const |
Returns the number of flops in the initialization phase. | |
virtual double | ComputeFlops () const |
Returns the total number of flops to computate the preconditioner. | |
virtual double | ApplyInverseFlops () const |
Returns the total number of flops to apply the preconditioner. | |
virtual const Teuchos::ParameterList & | List () const |
virtual std::ostream & | Print (std::ostream &os) const |
Prints on ostream basic information about this object. | |
Protected Member Functions | |
void | SetLabel (const char *Label_in) |
Sets the label. | |
void | SetIsInitialized (const bool IsInitialized_in) |
Sets IsInitialized_ . | |
void | SetIsComputed (const int IsComputed_in) |
Sets IsComputed_ . | |
void | SetNumInitialize (const int NumInitialize_in) |
Sets NumInitialize_ . | |
void | SetNumCompute (const int NumCompute_in) |
Sets NumCompute_ . | |
void | SetNumApplyInverse (const int NumApplyInverse_in) |
Sets NumApplyInverse_ . | |
void | SetInitializeTime (const double InitializeTime_in) |
Sets InitializeTime_ . | |
void | SetComputeTime (const double ComputeTime_in) |
Sets ComputeTime_ . | |
void | SetApplyInverseTime (const double ApplyInverseTime_in) |
Sets ApplyInverseTime_ . | |
void | SetComputeFlops (const double ComputeFlops_in) |
Sets ComputeFlops_ . | |
void | SetApplyInverseFlops (const double ApplyInverseFlops_in) |
Sets ComputeFlops_ . | |
void | SetList (const Teuchos::ParameterList &List_in) |
Set List_ . | |
Ifpack_Amesos: a class to use Amesos' factorizations as preconditioners.
Class Ifpack_Amesos enables the use of Amesos' factorizations as Ifpack_Preconditioners.
Ifpack_Amesos is just a bare-bone wrap to Amesos. Currently, the only parameter required recognized by SetParameters() is "amesos: solver type"
(defaulted to "Amesos_Klu"
), which defined the Amesos solver. The Teuchos list in input to SetParameters() is copied, then the copied list is used to set the parameters of the Amesos object.
This class works with matrices whose communicator contains only one process, that is, either serial matrices, or Ifpack_LocalFilter'd matrices.
Definition at line 81 of file Ifpack_Amesos.h.
|
virtual |
Applies the matrix to an Epetra_MultiVector.
X | - (In) A Epetra_MultiVector of dimension NumVectors to multiply with matrix. |
Y | - (Out) A Epetra_MultiVector of dimension NumVectors containing the result. |
Implements Epetra_Operator.
Definition at line 252 of file Ifpack_Amesos.cpp.
|
virtual |
Applies the preconditioner to X, returns the result in Y.
X | - (In) A Epetra_MultiVector of dimension NumVectors to be preconditioned. |
Y | - (Out) A Epetra_MultiVector of dimension NumVectors containing result. |
Implements Ifpack_Preconditioner.
Definition at line 261 of file Ifpack_Amesos.cpp.
References IsComputed().
|
virtual |
Computes the preconditioners.
Implements Ifpack_Preconditioner.
Definition at line 212 of file Ifpack_Amesos.cpp.
References Initialize(), and IsInitialized().
Referenced by Ifpack_Amesos().
|
virtual |
Initializes the preconditioners.
Implements Ifpack_Preconditioner.
Definition at line 126 of file Ifpack_Amesos.cpp.
References Comm().
Referenced by Compute(), and Ifpack_Amesos().
|
virtual |
Sets all the parameters for the preconditioner.
Parameters currently supported:
"amesos: solver type"
: Specifies the solver type for Amesos. Default: Amesos_Klu
.The input list will be copied, then passed to the Amesos object through Amesos::SetParameters().
Implements Ifpack_Preconditioner.
Definition at line 117 of file Ifpack_Amesos.cpp.
|
virtual |
If set true, transpose of this operator will be applied (not implemented).
This flag allows the transpose of the given operator to be used
implicitly.
UseTranspose_in | - (In) If true, multiply by the transpose of operator, otherwise just use operator. |
Implements Epetra_Operator.
Definition at line 239 of file Ifpack_Amesos.cpp.