IFPACK  Development
 All Classes Namespaces Files Functions Variables Enumerations Friends Pages
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | List of all members
Ifpack Class Reference

Ifpack: a function class to define Ifpack preconditioners. More...

#include <Ifpack.h>

Public Types

enum  parameter {
  absolute_threshold, relative_threshold, drop_tolerance, fill_tolerance,
  relax_value, level_fill, level_overlap, num_steps,
  use_reciprocal, overlap_mode
}
 
enum  EPrecType {
  POINT_RELAXATION, POINT_RELAXATION_STAND_ALONE, BLOCK_RELAXATION, BLOCK_RELAXATION_STAND_ALONE,
  BLOCK_RELAXATION_STAND_ALONE_ILU, BLOCK_RELAXATION_STAND_ALONE_ILUT, BLOCK_RELAXATION_STAND_ALONE_IC, BLOCK_RELAXATION_STAND_ALONE_SILU,
  BLOCK_RELAXATION_STAND_ALONE_AMESOS, BLOCK_RELAXATION_AMESOS, AMESOS, AMESOS_STAND_ALONE,
  IC, IC_STAND_ALONE, ICT, ICT_STAND_ALONE,
  ILU, ILU_STAND_ALONE, ILUT, ILUT_STAND_ALONE,
  SILU, CHEBYSHEV, POLYNOMIAL, KRYLOV,
  IHSS, SORA, TRIDI_RELAXATION, TRIDI_RELAXATION_STAND_ALONE
}
 Enum for the type of preconditioner.
 

Public Member Functions

Ifpack_PreconditionerCreate (const std::string PrecType, Epetra_RowMatrix *Matrix, const int overlap=0, bool overrideSerialDefault=false)
 Creates an instance of Ifpack_Preconditioner given the std::string name of the preconditioner type (can fail with bad input). More...
 
int SetParameters (int argc, char *argv[], Teuchos::ParameterList &List, std::string &PrecType, int &Overlap)
 Sets the options in List from the command line. More...
 

Static Public Member Functions

static const char * toString (const EPrecType precType)
 Function that gives the std::string name for preconditioner given its enumerication value.
 
static Ifpack_PreconditionerCreate (EPrecType PrecType, Epetra_RowMatrix *Matrix, const int overlap=0, bool overrideSerialDefault=false)
 Creates an instance of Ifpack_Preconditioner given the enum value of the preconditioner type (can not fail, no bad input possible). More...
 

Static Public Attributes

static const int numPrecTypes
 
static const EPrecType precTypeValues [numPrecTypes]
 List of the preconditioner types as enum values .
 
static const char * precTypeNames [numPrecTypes]
 List of preconditioner types as std::string values.
 
static const bool supportsUnsymmetric [numPrecTypes]
 List of bools that determines if the preconditioner type supports unsymmetric matrices.
 

Detailed Description

Ifpack: a function class to define Ifpack preconditioners.

Class Ifpack is a function class, that contains just one method: Create(). Using Create(), users can easily define a variety of IFPACK preconditioners.

Create requires 3 arguments:

The first argument can assume the following values:

Note
Objects in stand-alone mode cannot use reordering, variable overlap, and singleton filters. However, their construction can be slightly faster than the non stand-alone counterpart.

The following fragment of code shows the basic usage of this class.

#include "Ifpack.h"
...
Ifpack Factory;
Epetra_RowMatrix* A; // A is FillComplete()'d.
std::string PrecType = "ILU"; // use incomplete LU on each process
int OverlapLevel = 1; // one row of overlap among the processes
Ifpack_Preconditioner* Prec = Factory.Create(PrecType, A, OverlapLevel);
assert (Prec != 0);
Teuchos::ParameterList List;
List.set("fact: level-of-fill", 5); // use ILU(5)
IFPACK_CHK_ERR(Prec->SetParameters(List));
IFPACK_CHK_ERR(Prec->Initialize());
IFPACK_CHK_ERR(Prec->Compute());
// now Prec can be used as AztecOO preconditioner
// like for instance
AztecOO AztecOOSolver(*Problem);
// specify solver
AztecOOSolver.SetAztecOption(AZ_solver,AZ_gmres);
AztecOOSolver.SetAztecOption(AZ_output,32);
// Set Prec as preconditioning operator
AztecOOSolver.SetPrecOperator(Prec);
// Call the solver
AztecOOSolver.Iterate(1550,1e-8);
// print information on stdout
cout << *Prec;
// delete the preconditioner
delete Prec;
Author
Marzio Sala, (formally) SNL org. 1414
Date
Last updated on 25-Jan-05.

Definition at line 138 of file Ifpack.h.

Member Function Documentation

Ifpack_Preconditioner * Ifpack::Create ( EPrecType  PrecType,
Epetra_RowMatrix Matrix,
const int  overlap = 0,
bool  overrideSerialDefault = false 
)
static

Creates an instance of Ifpack_Preconditioner given the enum value of the preconditioner type (can not fail, no bad input possible).

Parameters
PrecType(In) - Enum value of preconditioner type to be created.
Matrix(In) - Matrix used to define the preconditioner
overlap(In) - specified overlap, defaulted to 0.

Definition at line 253 of file Ifpack.cpp.

References Epetra_Operator::Comm(), and Epetra_Comm::NumProc().

Referenced by Create().

Ifpack_Preconditioner * Ifpack::Create ( const std::string  PrecType,
Epetra_RowMatrix Matrix,
const int  overlap = 0,
bool  overrideSerialDefault = false 
)

Creates an instance of Ifpack_Preconditioner given the std::string name of the preconditioner type (can fail with bad input).

Parameters
PrecType(In) - String name of preconditioner type to be created.
Matrix(In) - Matrix used to define the preconditioner
overlap(In) - specified overlap, defaulted to 0.

Returns 0 if the preconditioner with that input name does not exist. Otherwise, return a newly created preconditioner object. Note that the client is responsible for calling delete on the returned object once it is finished using it!

Definition at line 393 of file Ifpack.cpp.

References Create().

int Ifpack::SetParameters ( int  argc,
char *  argv[],
Teuchos::ParameterList &  List,
std::string &  PrecType,
int &  Overlap 
)

Sets the options in List from the command line.

Note: If you want full support for all parameters, consider reading in a parameter list from an XML file as supported by the Teuchos helper function Teuchos::updateParametersFromXmlFile() or Teuchos::updateParametersFromXmlStream().

Definition at line 413 of file Ifpack.cpp.

Member Data Documentation

const int Ifpack::numPrecTypes
static
Initial value:
=
+7
+4
+8
+2
+7

Definition at line 195 of file Ifpack.h.


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