57 #ifdef HAVE_IFPACK_AMESOS 
   60 #ifdef HAVE_IFPACK_HIPS 
   63 #ifdef HAVE_IFPACK_SUPERLU 
   66 #ifdef HAVE_IFPACK_SUPPORTGRAPH 
   76 #include "Teuchos_CommandLineProcessor.hpp" 
   77 #include "Teuchos_StringToIntMap.hpp" 
   78 #include "Epetra_CrsMatrix.h" 
   94   ,POINT_RELAXATION_STAND_ALONE
 
   96   ,BLOCK_RELAXATION_STAND_ALONE
 
   97   ,BLOCK_RELAXATION_STAND_ALONE_ILU
 
   98   ,BLOCK_RELAXATION_STAND_ALONE_ILUT
 
   99   ,BLOCK_RELAXATION_STAND_ALONE_IC
 
  100 #ifdef HAVE_IFPACK_SUPERLU 
  101   ,BLOCK_RELAXATION_STAND_ALONE_SILU
 
  103 #ifdef HAVE_IFPACK_AMESOS 
  104   ,BLOCK_RELAXATION_STAND_ALONE_AMESOS
 
  105   ,BLOCK_RELAXATION_AMESOS
 
  108 #endif // HAVE_IFPACK_AMESOS 
  117 #ifdef HAVE_IFPACK_SPARSKIT 
  119 #endif // HAVE_IFPACK_SPARSKIT 
  120 #ifdef HAVE_IFPACK_HIPS 
  123 #ifdef HAVE_IFPACK_HYPRE 
  126 #ifdef HAVE_IFPACK_SUPERLU 
  129 #if defined (HAVE_IFPACK_SUPPORTGRAPH) && defined (HAVE_IFPACK_AMESOS) 
  132 #ifdef HAVE_IFPACK_SUPPORTGRAPH 
  141   ,TRIDI_RELAXATION_STAND_ALONE
 
  148   ,
"point relaxation stand-alone" 
  150   ,
"block relaxation stand-alone" 
  151   ,
"block relaxation stand-alone (ILU)" 
  152   ,
"block relaxation stand-alone (ILUT)" 
  153   ,
"block relaxation stand-alone (IC)" 
  154 #ifdef HAVE_IFPACK_SUPERLU 
  155   ,
"block relaxation stand-alone (SILU)" 
  157 #ifdef HAVE_IFPACK_AMESOS 
  158   ,
"block relaxation stand-alone (Amesos)" 
  159   ,
"block relaxation (Amesos)" 
  161   ,
"Amesos stand-alone" 
  171 #ifdef HAVE_IFPACK_SPARSKIT 
  174 #ifdef HAVE_IFPACK_HIPS 
  177 #ifdef HAVE_IFPACK_HYPRE 
  180 #ifdef HAVE_IFPACK_SUPERLU 
  183 #if defined (HAVE_IFPACK_SUPPORTGRAPH) && defined (HAVE_IFPACK_AMESOS) 
  186 #ifdef HAVE_IFPACK_SUPPORTGRAPH 
  195   ,
"tridi relaxation stand-alone" 
  208 #ifdef HAVE_IFPACK_SUPERLU 
  211 #ifdef HAVE_IFPACK_AMESOS 
  225 #ifdef HAVE_IFPACK_SPARSKIT 
  228 #ifdef HAVE_IFPACK_HIPS 
  231 #ifdef HAVE_IFPACK_HYPRE 
  234 #ifdef HAVE_IFPACK_SUPERLU 
  237 #if defined (HAVE_IFPACK_SUPPORTGRAPH) && defined (HAVE_IFPACK_AMESOS) 
  240 #ifdef HAVE_IFPACK_SUPPORTGRAPH 
  256                                       bool overrideSerialDefault)
 
  258   const bool serial = (Matrix->
Comm().
NumProc() == 1);
 
  262       if (serial && !overrideSerialDefault)
 
  269       if (serial && !overrideSerialDefault)
 
  282 #ifdef HAVE_IFPACK_SUPERLU 
  283     case BLOCK_RELAXATION_STAND_ALONE_SILU:
 
  286 #ifdef HAVE_IFPACK_AMESOS 
  287     case BLOCK_RELAXATION_STAND_ALONE_AMESOS:
 
  289     case BLOCK_RELAXATION_AMESOS:
 
  293       if (serial && !overrideSerialDefault)
 
  297     case AMESOS_STAND_ALONE:
 
  301       if (serial && !overrideSerialDefault)
 
  308       if (serial && !overrideSerialDefault)
 
  315       if (serial && !overrideSerialDefault)
 
  322       if (serial && !overrideSerialDefault)
 
  328 #ifdef HAVE_IFPACK_SPARSKIT 
  330       return(
new Ifpack_SPARSKIT(Matrix));
 
  332 #ifdef HAVE_IFPACK_HIPS 
  334       return(
new Ifpack_HIPS(Matrix));
 
  336 #ifdef HAVE_IFPACK_HYPRE 
  338       return(
new Ifpack_Hypre(Matrix));
 
  340 #ifdef HAVE_IFPACK_SUPERLU 
  342       return(
new Ifpack_SILU(Matrix));
 
  344 #if defined (HAVE_IFPACK_SUPPORTGRAPH) && defined (HAVE_IFPACK_AMESOS) 
  346       if (serial && !overrideSerialDefault)
 
  351 #ifdef HAVE_IFPACK_SUPPORTGRAPH 
  353       if (serial && !overrideSerialDefault)
 
  369 #ifdef HAVE_IFPACK_EPETRAEXT 
  371       return(
new Ifpack_IHSS(Matrix));
 
  373       return(
new Ifpack_SORa(Matrix));
 
  376      if (serial && !overrideSerialDefault)
 
  396                                       bool overrideSerialDefault)
 
  399     return Ifpack::Create(Teuchos::get<EPrecType>(::precTypeNameToIntMap,PrecType),Matrix,Overlap,overrideSerialDefault);
 
  422   std::string ifp_prec_type = 
"ILU";
 
  423   CLP.
setOption(
"ifp-prec-type",&ifp_prec_type,
"Preconditioner type");
 
  426   CLP.
setOption(
"ifp-overlap",&ifp_overlap,
"Overlap among processors");
 
  428   std::string ifp_relax_type = 
"Jacobi";
 
  429   CLP.
setOption(
"ifp-relax-type",&ifp_relax_type,
"Relaxation type");
 
  431   int ifp_relax_sweeps = 1;
 
  433                 &ifp_relax_sweeps,
"Number of sweeps for relaxation");
 
  435   double ifp_relax_damping = 1.0;
 
  437                 &ifp_relax_damping,
"Damping for relaxation");
 
  439   std::string ifp_part_type = 
"greedy";
 
  440   CLP.
setOption(
"ifp-part-type",&ifp_part_type,
"Partitioner type");
 
  442   int ifp_part_local = 1;
 
  443   CLP.
setOption(
"ifp-part-local",&ifp_part_local,
"number of local partitions");
 
  448   CLP.
parse(argc,argv);
 
  451   PrecType = ifp_prec_type;
 
  452   Overlap = ifp_overlap;
 
  455   List.
set(
"relaxation: type", ifp_relax_type);
 
  456   List.
set(
"relaxation: sweeps", ifp_relax_sweeps);
 
  457   List.
set(
"relaxation: damping factor", ifp_relax_damping);
 
  458   List.
set(
"partitioner: type", ifp_part_type);
 
  459   List.
set(
"partitioner: local parts", ifp_part_local);
 
static const int numPrecTypes
Ifpack_BlockRelaxation: a class to define block relaxation preconditioners of Epetra_RowMatrix's. 
void recogniseAllOptions(const bool &recogniseAllOptions)
static const EPrecType precTypeValues[numPrecTypes]
List of the preconditioner types as enum values . 
int SetParameters(int argc, char *argv[], Teuchos::ParameterList &List, std::string &PrecType, int &Overlap)
Sets the options in List from the command line. 
EPrecType
Enum for the type of preconditioner. 
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
static const char * precTypeNames[numPrecTypes]
List of preconditioner types as std::string values. 
Ifpack_ILU: A class for constructing and using an incomplete lower/upper (ILU) factorization of a giv...
Ifpack_Amesos: a class to use Amesos' factorizations as preconditioners. 
virtual const Epetra_Comm & Comm() const =0
Ifpack_AdditiveSchwarz: a class to define Additive Schwarz preconditioners of Epetra_RowMatrix's. 
void setOption(const char option_true[], const char option_false[], bool *option_val, const char documentation[]=NULL)
Ifpack_SparseContainer: a class for storing and solving linear systems using sparse matrices...
Ifpack_Chebyshev: class for preconditioning with Chebyshev polynomials in Ifpack. ...
Ifpack_Preconditioner: basic class for preconditioning in Ifpack. 
Ifpack_PointRelaxation: a class to define point relaxation preconditioners of for Epetra_RowMatrix's...
EParseCommandLineReturn parse(int argc, char *argv[], std::ostream *errout=&std::cerr) const 
void throwExceptions(const bool &throwExceptions)
Ifpack_Polynomial: class for preconditioning with least squares polynomials in Ifpack. 
Ifpack_Krylov: class for smoothing with Krylov solvers in Ifpack. 
virtual int NumProc() const =0
Ifpack_IC: A class for constructing and using an incomplete Cholesky factorization of a given Epetra_...
Ifpack_ICT: A class for constructing and using an incomplete Cholesky factorization of a given Epetra...
static Ifpack_Preconditioner * Create(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...
static const bool supportsUnsymmetric[numPrecTypes]
List of bools that determines if the preconditioner type supports unsymmetric matrices. 
Ifpack_ILUT: A class for constructing and using an incomplete LU factorization of a given Epetra_RowM...
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)