Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Private Attributes | List of all members
Ifpack_Amesos Class Reference

Ifpack_Amesos: a class to use Amesos' factorizations as preconditioners. More...

#include <Ifpack_Amesos.h>

Inheritance diagram for Ifpack_Amesos:
Inheritance graph
[legend]

Private Attributes

Teuchos::RefCountPtr< const
Epetra_RowMatrix
Matrix_
 Pointers to the matrix to be preconditioned. More...
 
Teuchos::RefCountPtr
< Epetra_LinearProblem
Problem_
 Linear problem required by Solver_. More...
 
Teuchos::RefCountPtr
< Amesos_BaseSolver
Solver_
 Amesos solver, use to apply the inverse of the local matrix. More...
 
Teuchos::ParameterList List_
 Contains a copy of the input parameter list. More...
 
std::string Label_
 Contains the label of this object. More...
 
bool IsEmpty_
 If true, the linear system on this processor is empty, thus the preconditioner is null operation. More...
 
bool IsInitialized_
 If true, the preconditioner has been successfully initialized. More...
 
bool IsComputed_
 If true, the preconditioner has been successfully computed. More...
 
bool UseTranspose_
 If true, the preconditioner solves for the transpose of the matrix. More...
 
int NumInitialize_
 Contains the number of successful calls to Initialize(). More...
 
int NumCompute_
 Contains the number of successful call to Compute(). More...
 
int NumApplyInverse_
 Contains the number of successful call to ApplyInverse(). More...
 
double InitializeTime_
 Contains the time for all successful calls to Initialize(). More...
 
double ComputeTime_
 Contains the time for all successful calls to Compute(). More...
 
double ApplyInverseTime_
 Contains the time for all successful calls to ApplyInverse(). More...
 
Teuchos::RefCountPtr< Epetra_TimeTime_
 Time object. More...
 
double ComputeFlops_
 Contains the number of flops for Compute(). More...
 
double ApplyInverseFlops_
 Contain sthe number of flops for ApplyInverse(). More...
 
double Condest_
 Contains the estimated condition number. More...
 
 Ifpack_Amesos (Epetra_RowMatrix *Matrix)
 Constructor. More...
 
 Ifpack_Amesos (const Ifpack_Amesos &rhs)
 Copy constructor. More...
 
Ifpack_Amesosoperator= (const Ifpack_Amesos &rhs)
 Operator=. More...
 
virtual ~Ifpack_Amesos ()
 Destructor. More...
 
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) More...
 
virtual const char * Label () const
 Returns a character string describing the operator. More...
 
virtual bool UseTranspose () const
 Returns the current UseTranspose setting. More...
 
virtual bool HasNormInf () const
 Returns true if the this object can provide an approximate Inf-norm, false otherwise. More...
 
virtual const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this operator. More...
 
virtual const Epetra_MapOperatorDomainMap () const
 Returns the Epetra_Map object associated with the domain of this operator. More...
 
virtual const Epetra_MapOperatorRangeMap () const
 Returns the Epetra_Map object associated with the range of this operator. More...
 
virtual bool IsInitialized () const
 Returns true is the preconditioner has been successfully initialized. More...
 
virtual int Initialize ()
 Initializes the preconditioners. More...
 
virtual bool IsComputed () const
 Returns true if the preconditioner has been successfully computed. More...
 
virtual int Compute ()
 Computes the preconditioners. More...
 
virtual int SetParameters (Teuchos::ParameterList &List)
 Sets all the parameters for the preconditioner. More...
 
virtual const Epetra_RowMatrixMatrix () const
 Returns a const reference to the internally stored matrix. More...
 
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. More...
 
virtual double Condest () const
 Returns the estimated condition number, never computes it. More...
 
virtual int NumInitialize () const
 Returns the number of calls to Initialize(). More...
 
virtual int NumCompute () const
 Returns the number of calls to Compute(). More...
 
virtual int NumApplyInverse () const
 Returns the number of calls to ApplyInverse(). More...
 
virtual double InitializeTime () const
 Returns the total time spent in Initialize(). More...
 
virtual double ComputeTime () const
 Returns the total time spent in Compute(). More...
 
virtual double ApplyInverseTime () const
 Returns the total time spent in ApplyInverse(). More...
 
virtual double InitializeFlops () const
 Returns the number of flops in the initialization phase. More...
 
virtual double ComputeFlops () const
 Returns the total number of flops to computate the preconditioner. More...
 
virtual double ApplyInverseFlops () const
 Returns the total number of flops to apply the preconditioner. More...
 
virtual const
Teuchos::ParameterList
List () const
 
virtual std::ostream & Print (std::ostream &os) const
 Prints on ostream basic information about this object. More...
 
void SetLabel (const char *Label_in)
 Sets the label. More...
 
void SetIsInitialized (const bool IsInitialized_in)
 Sets IsInitialized_. More...
 
void SetIsComputed (const int IsComputed_in)
 Sets IsComputed_. More...
 
void SetNumInitialize (const int NumInitialize_in)
 Sets NumInitialize_. More...
 
void SetNumCompute (const int NumCompute_in)
 Sets NumCompute_. More...
 
void SetNumApplyInverse (const int NumApplyInverse_in)
 Sets NumApplyInverse_. More...
 
void SetInitializeTime (const double InitializeTime_in)
 Sets InitializeTime_. More...
 
void SetComputeTime (const double ComputeTime_in)
 Sets ComputeTime_. More...
 
void SetApplyInverseTime (const double ApplyInverseTime_in)
 Sets ApplyInverseTime_. More...
 
void SetComputeFlops (const double ComputeFlops_in)
 Sets ComputeFlops_. More...
 
void SetApplyInverseFlops (const double ApplyInverseFlops_in)
 Sets ComputeFlops_. More...
 
void SetList (const Teuchos::ParameterList &List_in)
 Set List_. More...
 

Additional Inherited Members

Detailed Description

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.

Warning
The number of flops is NOT updated.
Author
Marzio Sala, SNL 9214.
Date
Last update Oct-04.

Definition at line 81 of file Ifpack_Amesos.h.

Constructor & Destructor Documentation

Ifpack_Amesos::Ifpack_Amesos ( Epetra_RowMatrix Matrix)

Constructor.

Definition at line 59 of file Ifpack_Amesos.cpp.

Ifpack_Amesos::Ifpack_Amesos ( const Ifpack_Amesos rhs)

Copy constructor.

Definition at line 80 of file Ifpack_Amesos.cpp.

virtual Ifpack_Amesos::~Ifpack_Amesos ( )
inlinevirtual

Destructor.

Definition at line 98 of file Ifpack_Amesos.h.

Member Function Documentation

Ifpack_Amesos& Ifpack_Amesos::operator= ( const Ifpack_Amesos rhs)

Operator=.

int Ifpack_Amesos::SetUseTranspose ( bool  UseTranspose_in)
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.

Parameters
UseTranspose_in- (In) If true, multiply by the transpose of operator, otherwise just use operator.
Returns
Integer error code, set to 0 if successful. Set to -1 if this implementation does not support transpose.

Implements Epetra_Operator.

Definition at line 239 of file Ifpack_Amesos.cpp.

int Ifpack_Amesos::Apply ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Applies the matrix to an Epetra_MultiVector.

Parameters
X- (In) A Epetra_MultiVector of dimension NumVectors to multiply with matrix.
Y- (Out) A Epetra_MultiVector of dimension NumVectors containing the result.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_Operator.

Definition at line 252 of file Ifpack_Amesos.cpp.

int Ifpack_Amesos::ApplyInverse ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Applies the preconditioner to X, returns the result in Y.

Parameters
X- (In) A Epetra_MultiVector of dimension NumVectors to be preconditioned.
Y- (Out) A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.
Warning
In order to work with AztecOO, any implementation of this method must support the case where X and Y are the same object.

Implements Ifpack_Preconditioner.

Definition at line 261 of file Ifpack_Amesos.cpp.

double Ifpack_Amesos::NormInf ( ) const
virtual

Returns the infinity norm of the global matrix (not implemented)

Implements Epetra_Operator.

Definition at line 295 of file Ifpack_Amesos.cpp.

const char * Ifpack_Amesos::Label ( ) const
virtual

Returns a character string describing the operator.

Implements Epetra_Operator.

Definition at line 301 of file Ifpack_Amesos.cpp.

bool Ifpack_Amesos::UseTranspose ( ) const
virtual

Returns the current UseTranspose setting.

Implements Epetra_Operator.

Definition at line 307 of file Ifpack_Amesos.cpp.

bool Ifpack_Amesos::HasNormInf ( ) const
virtual

Returns true if the this object can provide an approximate Inf-norm, false otherwise.

Implements Epetra_Operator.

Definition at line 313 of file Ifpack_Amesos.cpp.

const Epetra_Comm & Ifpack_Amesos::Comm ( ) const
virtual

Returns a pointer to the Epetra_Comm communicator associated with this operator.

Implements Epetra_Operator.

Definition at line 319 of file Ifpack_Amesos.cpp.

const Epetra_Map & Ifpack_Amesos::OperatorDomainMap ( ) const
virtual

Returns the Epetra_Map object associated with the domain of this operator.

Implements Epetra_Operator.

Definition at line 325 of file Ifpack_Amesos.cpp.

const Epetra_Map & Ifpack_Amesos::OperatorRangeMap ( ) const
virtual

Returns the Epetra_Map object associated with the range of this operator.

Implements Epetra_Operator.

Definition at line 331 of file Ifpack_Amesos.cpp.

virtual bool Ifpack_Amesos::IsInitialized ( ) const
inlinevirtual

Returns true is the preconditioner has been successfully initialized.

Implements Ifpack_Preconditioner.

Definition at line 173 of file Ifpack_Amesos.h.

int Ifpack_Amesos::Initialize ( )
virtual

Initializes the preconditioners.

Returns
0 if successful, 1 if problems occurred.

Implements Ifpack_Preconditioner.

Definition at line 126 of file Ifpack_Amesos.cpp.

virtual bool Ifpack_Amesos::IsComputed ( ) const
inlinevirtual

Returns true if the preconditioner has been successfully computed.

Implements Ifpack_Preconditioner.

Definition at line 185 of file Ifpack_Amesos.h.

int Ifpack_Amesos::Compute ( )
virtual

Computes the preconditioners.

Returns
0 if successful, 1 if problems occurred.

Implements Ifpack_Preconditioner.

Definition at line 212 of file Ifpack_Amesos.cpp.

int Ifpack_Amesos::SetParameters ( Teuchos::ParameterList List)
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 const Epetra_RowMatrix& Ifpack_Amesos::Matrix ( ) const
inlinevirtual

Returns a const reference to the internally stored matrix.

Implements Ifpack_Preconditioner.

Definition at line 211 of file Ifpack_Amesos.h.

double Ifpack_Amesos::Condest ( const Ifpack_CondestType  CT = Ifpack_Cheap,
const int  MaxIters = 1550,
const double  Tol = 1e-9,
Epetra_RowMatrix Matrix_in = 0 
)
virtual

Returns the estimated condition number, computes it if necessary.

Implements Ifpack_Preconditioner.

Definition at line 337 of file Ifpack_Amesos.cpp.

virtual double Ifpack_Amesos::Condest ( ) const
inlinevirtual

Returns the estimated condition number, never computes it.

Implements Ifpack_Preconditioner.

Definition at line 223 of file Ifpack_Amesos.h.

virtual int Ifpack_Amesos::NumInitialize ( ) const
inlinevirtual

Returns the number of calls to Initialize().

Implements Ifpack_Preconditioner.

Definition at line 229 of file Ifpack_Amesos.h.

virtual int Ifpack_Amesos::NumCompute ( ) const
inlinevirtual

Returns the number of calls to Compute().

Implements Ifpack_Preconditioner.

Definition at line 235 of file Ifpack_Amesos.h.

virtual int Ifpack_Amesos::NumApplyInverse ( ) const
inlinevirtual

Returns the number of calls to ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 241 of file Ifpack_Amesos.h.

virtual double Ifpack_Amesos::InitializeTime ( ) const
inlinevirtual

Returns the total time spent in Initialize().

Implements Ifpack_Preconditioner.

Definition at line 247 of file Ifpack_Amesos.h.

virtual double Ifpack_Amesos::ComputeTime ( ) const
inlinevirtual

Returns the total time spent in Compute().

Implements Ifpack_Preconditioner.

Definition at line 253 of file Ifpack_Amesos.h.

virtual double Ifpack_Amesos::ApplyInverseTime ( ) const
inlinevirtual

Returns the total time spent in ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 259 of file Ifpack_Amesos.h.

virtual double Ifpack_Amesos::InitializeFlops ( ) const
inlinevirtual

Returns the number of flops in the initialization phase.

Implements Ifpack_Preconditioner.

Definition at line 265 of file Ifpack_Amesos.h.

virtual double Ifpack_Amesos::ComputeFlops ( ) const
inlinevirtual

Returns the total number of flops to computate the preconditioner.

Implements Ifpack_Preconditioner.

Definition at line 271 of file Ifpack_Amesos.h.

virtual double Ifpack_Amesos::ApplyInverseFlops ( ) const
inlinevirtual

Returns the total number of flops to apply the preconditioner.

Implements Ifpack_Preconditioner.

Definition at line 277 of file Ifpack_Amesos.h.

virtual const Teuchos::ParameterList& Ifpack_Amesos::List ( ) const
inlinevirtual

Definition at line 283 of file Ifpack_Amesos.h.

std::ostream & Ifpack_Amesos::Print ( std::ostream &  os) const
virtual

Prints on ostream basic information about this object.

Implements Ifpack_Preconditioner.

Definition at line 352 of file Ifpack_Amesos.cpp.

void Ifpack_Amesos::SetLabel ( const char *  Label_in)
inlineprotected

Sets the label.

Definition at line 298 of file Ifpack_Amesos.h.

void Ifpack_Amesos::SetIsInitialized ( const bool  IsInitialized_in)
inlineprotected

Sets IsInitialized_.

Definition at line 304 of file Ifpack_Amesos.h.

void Ifpack_Amesos::SetIsComputed ( const int  IsComputed_in)
inlineprotected

Sets IsComputed_.

Definition at line 310 of file Ifpack_Amesos.h.

void Ifpack_Amesos::SetNumInitialize ( const int  NumInitialize_in)
inlineprotected

Sets NumInitialize_.

Definition at line 316 of file Ifpack_Amesos.h.

void Ifpack_Amesos::SetNumCompute ( const int  NumCompute_in)
inlineprotected

Sets NumCompute_.

Definition at line 322 of file Ifpack_Amesos.h.

void Ifpack_Amesos::SetNumApplyInverse ( const int  NumApplyInverse_in)
inlineprotected

Sets NumApplyInverse_.

Definition at line 328 of file Ifpack_Amesos.h.

void Ifpack_Amesos::SetInitializeTime ( const double  InitializeTime_in)
inlineprotected

Sets InitializeTime_.

Definition at line 334 of file Ifpack_Amesos.h.

void Ifpack_Amesos::SetComputeTime ( const double  ComputeTime_in)
inlineprotected

Sets ComputeTime_.

Definition at line 340 of file Ifpack_Amesos.h.

void Ifpack_Amesos::SetApplyInverseTime ( const double  ApplyInverseTime_in)
inlineprotected

Sets ApplyInverseTime_.

Definition at line 346 of file Ifpack_Amesos.h.

void Ifpack_Amesos::SetComputeFlops ( const double  ComputeFlops_in)
inlineprotected

Sets ComputeFlops_.

Definition at line 352 of file Ifpack_Amesos.h.

void Ifpack_Amesos::SetApplyInverseFlops ( const double  ApplyInverseFlops_in)
inlineprotected

Sets ComputeFlops_.

Definition at line 358 of file Ifpack_Amesos.h.

void Ifpack_Amesos::SetList ( const Teuchos::ParameterList List_in)
inlineprotected

Set List_.

Definition at line 364 of file Ifpack_Amesos.h.

Member Data Documentation

Teuchos::RefCountPtr<const Epetra_RowMatrix> Ifpack_Amesos::Matrix_
private

Pointers to the matrix to be preconditioned.

Definition at line 373 of file Ifpack_Amesos.h.

Teuchos::RefCountPtr<Epetra_LinearProblem> Ifpack_Amesos::Problem_
private

Linear problem required by Solver_.

Definition at line 376 of file Ifpack_Amesos.h.

Teuchos::RefCountPtr<Amesos_BaseSolver> Ifpack_Amesos::Solver_
private

Amesos solver, use to apply the inverse of the local matrix.

Definition at line 378 of file Ifpack_Amesos.h.

Teuchos::ParameterList Ifpack_Amesos::List_
private

Contains a copy of the input parameter list.

Definition at line 380 of file Ifpack_Amesos.h.

std::string Ifpack_Amesos::Label_
private

Contains the label of this object.

Definition at line 383 of file Ifpack_Amesos.h.

bool Ifpack_Amesos::IsEmpty_
private

If true, the linear system on this processor is empty, thus the preconditioner is null operation.

Definition at line 385 of file Ifpack_Amesos.h.

bool Ifpack_Amesos::IsInitialized_
private

If true, the preconditioner has been successfully initialized.

Definition at line 387 of file Ifpack_Amesos.h.

bool Ifpack_Amesos::IsComputed_
private

If true, the preconditioner has been successfully computed.

Definition at line 389 of file Ifpack_Amesos.h.

bool Ifpack_Amesos::UseTranspose_
private

If true, the preconditioner solves for the transpose of the matrix.

Definition at line 391 of file Ifpack_Amesos.h.

int Ifpack_Amesos::NumInitialize_
private

Contains the number of successful calls to Initialize().

Definition at line 394 of file Ifpack_Amesos.h.

int Ifpack_Amesos::NumCompute_
private

Contains the number of successful call to Compute().

Definition at line 396 of file Ifpack_Amesos.h.

int Ifpack_Amesos::NumApplyInverse_
mutableprivate

Contains the number of successful call to ApplyInverse().

Definition at line 398 of file Ifpack_Amesos.h.

double Ifpack_Amesos::InitializeTime_
private

Contains the time for all successful calls to Initialize().

Definition at line 401 of file Ifpack_Amesos.h.

double Ifpack_Amesos::ComputeTime_
private

Contains the time for all successful calls to Compute().

Definition at line 403 of file Ifpack_Amesos.h.

double Ifpack_Amesos::ApplyInverseTime_
mutableprivate

Contains the time for all successful calls to ApplyInverse().

Definition at line 405 of file Ifpack_Amesos.h.

Teuchos::RefCountPtr<Epetra_Time> Ifpack_Amesos::Time_
private

Time object.

Definition at line 407 of file Ifpack_Amesos.h.

double Ifpack_Amesos::ComputeFlops_
private

Contains the number of flops for Compute().

Definition at line 410 of file Ifpack_Amesos.h.

double Ifpack_Amesos::ApplyInverseFlops_
private

Contain sthe number of flops for ApplyInverse().

Definition at line 412 of file Ifpack_Amesos.h.

double Ifpack_Amesos::Condest_
private

Contains the estimated condition number.

Definition at line 415 of file Ifpack_Amesos.h.


The documentation for this class was generated from the following files: