55 const std::string AztecSolver_name = 
"Aztec Solver";
 
   57 const std::string AztecPreconditioner_name = 
"Aztec Preconditioner";
 
   60   AZTEC_PREC_NONE, AZTEC_PREC_ILU, AZTEC_PREC_ILUT, AZTEC_PREC_JACOBI,
 
   61   AZTEC_PREC_SYMMGS, AZTEC_PREC_POLY, AZTEC_PREC_LSPOLY
 
   73 const std::string Overlap_name = 
"Overlap";
 
   75 const std::string GraphFill_name = 
"Graph Fill";
 
   77 const std::string DropTolerance_name = 
"Drop Tolerance";
 
   79 const std::string FillFactor_name = 
"Fill Factor";
 
   81 const std::string Steps_name = 
"Steps";
 
   83 const std::string PolynomialOrder_name = 
"Polynomial Order";
 
   85 const std::string RCMReordering_name = 
"RCM Reordering";
 
   87 const std::string Orthogonalization_name = 
"Orthogonalization";
 
   89 const std::string SizeOfKrylovSubspace_name = 
"Size of Krylov Subspace";
 
   91 const std::string ConvergenceTest_name = 
"Convergence Test";
 
   93 const std::string IllConditioningThreshold_name = 
"Ill-Conditioning Threshold";
 
   95 const std::string OutputFrequency_name = 
"Output Frequency";
 
  106   using Teuchos::getIntegralValue;
 
  107   using Teuchos::getParameter;
 
  114   solver->SetAztecOption(
 
  116     ,getIntegralValue<int>(*pl,AztecSolver_name)
 
  120     getIntegralValue<EAztecPreconditioner>(
 
  121       *pl,AztecPreconditioner_name
 
  125     case AZTEC_PREC_NONE:
 
  126       solver->SetAztecOption(AZ_precond,AZ_none);
 
  129       solver->SetAztecOption(AZ_precond,AZ_dom_decomp);
 
  130       solver->SetAztecOption(AZ_overlap,getParameter<int>(*pl,Overlap_name));
 
  131       solver->SetAztecOption(AZ_subdomain_solve,AZ_ilu);
 
  132       solver->SetAztecOption(AZ_graph_fill,getParameter<int>(*pl,GraphFill_name));
 
  134     case AZTEC_PREC_ILUT:
 
  135       solver->SetAztecOption(AZ_precond,AZ_dom_decomp);
 
  136       solver->SetAztecOption(AZ_overlap,getParameter<int>(*pl,Overlap_name));
 
  137       solver->SetAztecOption(AZ_subdomain_solve,AZ_ilut);
 
  138       solver->SetAztecParam(AZ_drop,getParameter<double>(*pl,DropTolerance_name));
 
  139       solver->SetAztecParam(AZ_ilut_fill,getParameter<double>(*pl,FillFactor_name));
 
  141     case AZTEC_PREC_JACOBI:
 
  142       solver->SetAztecOption(AZ_precond,AZ_Jacobi);
 
  143       solver->SetAztecOption(AZ_poly_ord,getParameter<int>(*pl,Steps_name));
 
  145     case AZTEC_PREC_SYMMGS:
 
  146       solver->SetAztecOption(AZ_precond,AZ_sym_GS);
 
  147       solver->SetAztecOption(AZ_poly_ord,getParameter<int>(*pl,Steps_name));
 
  149     case AZTEC_PREC_POLY:
 
  150       solver->SetAztecOption(AZ_precond,AZ_Neumann);
 
  151       solver->SetAztecOption(AZ_poly_ord,getParameter<int>(*pl,PolynomialOrder_name));
 
  153     case AZTEC_PREC_LSPOLY:
 
  154       solver->SetAztecOption(AZ_precond,AZ_ls);
 
  155       solver->SetAztecOption(AZ_poly_ord,getParameter<int>(*pl,PolynomialOrder_name));
 
  161   solver->SetAztecOption(
 
  163     ,getIntegralValue<int>(*pl,RCMReordering_name)
 
  166   solver->SetAztecOption(
 
  168     ,getIntegralValue<int>(*pl,Orthogonalization_name)
 
  171   solver->SetAztecOption(AZ_kspace,getParameter<int>(*pl,SizeOfKrylovSubspace_name));
 
  173   solver->SetAztecOption(
 
  175     ,getIntegralValue<int>(*pl,ConvergenceTest_name)
 
  178   solver->SetAztecParam(
 
  179     AZ_ill_cond_thresh, getParameter<double>(*pl,IllConditioningThreshold_name)
 
  182   solver->SetAztecOption(
 
  183     AZ_output, getParameter<int>(*pl,OutputFrequency_name)
 
  188 #endif // TEUCHOS_DEBUG 
  203   using Teuchos::tuple;
 
  204   using Teuchos::setStringToIntegralParameter;
 
  205   using Teuchos::setIntParameter;
 
  206   using Teuchos::setDoubleParameter;
 
  221   RCP<ParameterList> pl = validAztecOOParams;
 
  222   if(pl.get()) 
return pl;
 
  223   pl = validAztecOOParams = 
rcp(
new ParameterList());
 
  225   setStringToIntegralParameter<int>(
 
  226     AztecSolver_name, 
"GMRES",
 
  227     "Type of linear solver algorithm to use.",
 
  228     tuple<std::string>(
"CG",
"GMRES",
"CGS",
"TFQMR",
"BiCGStab",
"LU",
"GMRESR",
"FixedPoint"),
 
  229     tuple<int>(AZ_cg,AZ_gmres,AZ_cgs,AZ_tfqmr,AZ_bicgstab,AZ_lu,AZ_GMRESR,AZ_fixed_pt),
 
  232   setStringToIntegralParameter<EAztecPreconditioner>(
 
  233     AztecPreconditioner_name, 
"ilu",
 
  234     "Type of internal preconditioner to use.\n" 
  235     "Note! this preconditioner will only be used if the input operator\n" 
  236     "supports the Epetra_RowMatrix interface and the client does not pass\n" 
  237     "in an external preconditioner!",
 
  239       "none",
"ilu",
"ilut",
"Jacobi",
 
  240       "Symmetric Gauss-Seidel",
"Polynomial",
"Least-squares Polynomial" 
  242     tuple<EAztecPreconditioner>(
 
  243       AZTEC_PREC_NONE,AZTEC_PREC_ILU,AZTEC_PREC_ILUT,AZTEC_PREC_JACOBI,
 
  244       AZTEC_PREC_SYMMGS,AZTEC_PREC_POLY,AZTEC_PREC_LSPOLY
 
  250     "The amount of overlap used for the internal \"ilu\" and \"ilut\" preconditioners.",
 
  255     "The amount of fill allowed for the internal \"ilu\" preconditioner.",
 
  259     DropTolerance_name, 0.0,
 
  260     "The tolerance below which an entry from the factors of an internal \"ilut\"\n" 
  261     "preconditioner will be dropped.",
 
  265     FillFactor_name, 1.0,
 
  266     "The amount of fill allowed for an internal \"ilut\" preconditioner.",
 
  271     "Number of steps taken for the \"Jacobi\" or the \"Symmetric Gauss-Seidel\"\n" 
  272     "internal preconditioners for each preconditioner application.",
 
  276     PolynomialOrder_name, 3,
 
  277     "The order for of the polynomials used for the \"Polynomial\" and\n" 
  278     "\"Least-squares Polynomial\" internal preconditioners.",
 
  281   setStringToIntegralParameter<int>(
 
  282     RCMReordering_name, 
"Disabled",
 
  283     "Determines if RCM reordering is used with the internal\n" 
  284     "\"ilu\" or \"ilut\" preconditioners.",
 
  285     tuple<std::string>(
"Enabled",
"Disabled"),
 
  289   setStringToIntegralParameter<int>(
 
  290     Orthogonalization_name, 
"Classical",
 
  291     "The type of orthogonalization to use with the \"GMRES\" solver.",
 
  292     tuple<std::string>(
"Classical",
"Modified"),
 
  293     tuple<int>(AZ_classic,AZ_modified),
 
  297     SizeOfKrylovSubspace_name, 300,
 
  298     "The maximum size of the Krylov subspace used with \"GMRES\" before\n" 
  299     "a restart is performed.",
 
  302   setStringToIntegralParameter<int>(
 
  303     ConvergenceTest_name, 
"r0", 
 
  304     "The convergence test to use for terminating the iterative solver.",
 
  305     tuple<std::string>(
"r0",
"rhs",
"Anorm",
"no scaling",
"sol"),
 
  306     tuple<int>(AZ_r0,AZ_rhs,AZ_Anorm,AZ_noscaled,AZ_sol),
 
  310     IllConditioningThreshold_name, 1e+11,
 
  311     "The threshold tolerance above which a system is considered\n" 
  316     OutputFrequency_name, 0, 
 
  317     "The number of iterations between each output of the solver's progress.",
 
void setAztecOOParameters(Teuchos::ParameterList *pl, AztecOO *solver)
Setup an AztecOO solver object with a set of parameters. 
 
static void addConverter(RCP< const ParameterEntryValidator > validator, RCP< ValidatorXMLConverter > converterToAdd)
 
std::istringstream & operator>>(std::istringstream &in, TwoDArray< T > &array)
 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
 
void validateParametersAndSetDefaults(ParameterList const &validParamList, int const depth=1000)
 
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const 
 
Teuchos::RCP< const Teuchos::ParameterList > getValidAztecOOParameters()
Return the list of all valid AztecOO parameters (to validate against). 
 
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)