Stratimikos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AztecOOParameterList.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Stratimikos: Thyra-based strategies for linear solvers
5 // Copyright (2006) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #include "AztecOOParameterList.hpp"
46 
47 namespace {
48 
49 //
50 // Define the names of the different parameters. Since the name of a
51 // parameter is used several times, it is a good idea to define a variable that
52 // stores the std::string name so that typing errors get caught at compile-time.
53 //
54 
55 const std::string AztecSolver_name = "Aztec Solver";
56 
57 const std::string AztecPreconditioner_name = "Aztec Preconditioner";
58 
60  AZTEC_PREC_NONE, AZTEC_PREC_ILU, AZTEC_PREC_ILUT, AZTEC_PREC_JACOBI,
61  AZTEC_PREC_SYMMGS, AZTEC_PREC_POLY, AZTEC_PREC_LSPOLY
62 };
63 
65 inline std::istream& operator>>(std::istream& is, EAztecPreconditioner& prec){
66  int intval;
67  is >> intval;
68  prec = (EAztecPreconditioner)intval;
69  return is;
70 }
71 
72 
73 const std::string Overlap_name = "Overlap";
74 
75 const std::string GraphFill_name = "Graph Fill";
76 
77 const std::string DropTolerance_name = "Drop Tolerance";
78 
79 const std::string FillFactor_name = "Fill Factor";
80 
81 const std::string Steps_name = "Steps";
82 
83 const std::string PolynomialOrder_name = "Polynomial Order";
84 
85 const std::string RCMReordering_name = "RCM Reordering";
86 
87 const std::string Orthogonalization_name = "Orthogonalization";
88 
89 const std::string SizeOfKrylovSubspace_name = "Size of Krylov Subspace";
90 
91 const std::string ConvergenceTest_name = "Convergence Test";
92 
93 const std::string IllConditioningThreshold_name = "Ill-Conditioning Threshold";
94 
95 const std::string OutputFrequency_name = "Output Frequency";
96 
97 Teuchos::RCP<Teuchos::ParameterList> validAztecOOParams;
98 
99 } // namespace
100 
103  ,AztecOO *solver
104  )
105 {
106  using Teuchos::getIntegralValue;
107  using Teuchos::getParameter;
108  TEUCHOS_TEST_FOR_EXCEPT(pl==NULL);
109  TEUCHOS_TEST_FOR_EXCEPT(solver==NULL);
110  // Validate the parameters and set their defaults! This also sets the
111  // validators needed to read in the parameters in an integral form.
113  // Aztec Solver
114  solver->SetAztecOption(
115  AZ_solver
116  ,getIntegralValue<int>(*pl,AztecSolver_name)
117  );
118  // Aztec Preconditioner
119  switch(
120  getIntegralValue<EAztecPreconditioner>(
121  *pl,AztecPreconditioner_name
122  )
123  )
124  {
125  case AZTEC_PREC_NONE:
126  solver->SetAztecOption(AZ_precond,AZ_none);
127  break;
128  case AZTEC_PREC_ILU:
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));
133  break;
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));
140  break;
141  case AZTEC_PREC_JACOBI:
142  solver->SetAztecOption(AZ_precond,AZ_Jacobi);
143  solver->SetAztecOption(AZ_poly_ord,getParameter<int>(*pl,Steps_name));
144  break;
145  case AZTEC_PREC_SYMMGS:
146  solver->SetAztecOption(AZ_precond,AZ_sym_GS);
147  solver->SetAztecOption(AZ_poly_ord,getParameter<int>(*pl,Steps_name));
148  break;
149  case AZTEC_PREC_POLY:
150  solver->SetAztecOption(AZ_precond,AZ_Neumann);
151  solver->SetAztecOption(AZ_poly_ord,getParameter<int>(*pl,PolynomialOrder_name));
152  break;
153  case AZTEC_PREC_LSPOLY:
154  solver->SetAztecOption(AZ_precond,AZ_ls);
155  solver->SetAztecOption(AZ_poly_ord,getParameter<int>(*pl,PolynomialOrder_name));
156  break;
157  default:
158  TEUCHOS_TEST_FOR_EXCEPT(true); // Should never get here!
159  }
160  // RCM Reordering (in conjunction with domain decomp preconditioning)
161  solver->SetAztecOption(
162  AZ_reorder
163  ,getIntegralValue<int>(*pl,RCMReordering_name)
164  );
165  // Gram-Schmidt orthogonalization procedure (GMRES only)
166  solver->SetAztecOption(
167  AZ_orthog
168  ,getIntegralValue<int>(*pl,Orthogonalization_name)
169  );
170  // Size of the krylov subspace
171  solver->SetAztecOption(AZ_kspace,getParameter<int>(*pl,SizeOfKrylovSubspace_name));
172  // Convergence criteria to use in the linear solver
173  solver->SetAztecOption(
174  AZ_conv
175  ,getIntegralValue<int>(*pl,ConvergenceTest_name)
176  );
177  // Set the ill-conditioning threshold for the upper hessenberg matrix
178  solver->SetAztecParam(
179  AZ_ill_cond_thresh, getParameter<double>(*pl,IllConditioningThreshold_name)
180  );
181  // Frequency of linear solve residual output
182  solver->SetAztecOption(
183  AZ_output, getParameter<int>(*pl,OutputFrequency_name)
184  );
185 #ifdef TEUCHOS_DEBUG
186  // Check to make sure that I did not use the PL incorrectly!
188 #endif // TEUCHOS_DEBUG
189 }
190 
193 {
194  //
195  // This function defines the valid parameter list complete with validators
196  // and default values. The default values never need to be repeated because
197  // if the use of the function validateParametersAndSetDefaults(...) used
198  // above in setAztecOOParameters(...). Also, the validators do not need to
199  // be kept track of since they will be set in the input list also.
200  //
201  using Teuchos::RCP;
202  using Teuchos::rcp;
203  using Teuchos::tuple;
204  using Teuchos::setStringToIntegralParameter;
205  using Teuchos::setIntParameter;
206  using Teuchos::setDoubleParameter;
217  >::getDummyObject(),
219  EAztecPreconditioner> >::getDummyObject());
220  //
221  RCP<ParameterList> pl = validAztecOOParams;
222  if(pl.get()) return pl;
223  pl = validAztecOOParams = rcp(new ParameterList());
224  //
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),
230  &*pl
231  );
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!",
238  tuple<std::string>(
239  "none","ilu","ilut","Jacobi",
240  "Symmetric Gauss-Seidel","Polynomial","Least-squares Polynomial"
241  ),
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
245  ),
246  &*pl
247  );
248  setIntParameter(
249  Overlap_name, 0,
250  "The amount of overlap used for the internal \"ilu\" and \"ilut\" preconditioners.",
251  &*pl
252  );
253  setIntParameter(
254  GraphFill_name, 0,
255  "The amount of fill allowed for the internal \"ilu\" preconditioner.",
256  &*pl
257  );
258  setDoubleParameter(
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.",
262  &*pl
263  );
264  setDoubleParameter(
265  FillFactor_name, 1.0,
266  "The amount of fill allowed for an internal \"ilut\" preconditioner.",
267  &*pl
268  );
269  setIntParameter(
270  Steps_name, 3,
271  "Number of steps taken for the \"Jacobi\" or the \"Symmetric Gauss-Seidel\"\n"
272  "internal preconditioners for each preconditioner application.",
273  &*pl
274  );
275  setIntParameter(
276  PolynomialOrder_name, 3,
277  "The order for of the polynomials used for the \"Polynomial\" and\n"
278  "\"Least-squares Polynomial\" internal preconditioners.",
279  &*pl
280  );
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"),
286  tuple<int>(1,0),
287  &*pl
288  );
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),
294  &*pl
295  );
296  setIntParameter(
297  SizeOfKrylovSubspace_name, 300,
298  "The maximum size of the Krylov subspace used with \"GMRES\" before\n"
299  "a restart is performed.",
300  &*pl
301  );
302  setStringToIntegralParameter<int>(
303  ConvergenceTest_name, "r0", // Same as "rhs" when x=0
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),
307  &*pl
308  );
309  setDoubleParameter(
310  IllConditioningThreshold_name, 1e+11,
311  "The threshold tolerance above which a system is considered\n"
312  "ill conditioned.",
313  &*pl
314  );
315  setIntParameter(
316  OutputFrequency_name, 0, // By default, no output from Aztec!
317  "The number of iterations between each output of the solver's progress.",
318  &*pl
319  );
320  //
321  return pl;
322 }
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)
EAztecPreconditioner
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)