00001 #include "TSFAztecSolver.hpp"
00002 #include "TSFEpetraVector.hpp"
00003 #include "TSFEpetraMatrix.hpp"
00004 #include "Ifpack_Preconditioner.h"
00005 #include "Ifpack.h"
00006 #include "EpetraTSFOperator.hpp"
00007 #include "Teuchos_basic_oblackholestream.hpp"
00008
00009
00010
00011 #ifndef HAVE_TEUCHOS_EXPLICIT_INSTANTIATION
00012 #include "TSFVectorImpl.hpp"
00013 #include "TSFLinearOperatorImpl.hpp"
00014 #include "TSFLinearSolverImpl.hpp"
00015 #endif
00016
00017 #ifdef HAVE_ML
00018 #include "ml_include.h"
00019 #include "ml_epetra_utils.h"
00020 #include "ml_epetra_operator.h"
00021 #include "ml_aztec_utils.h"
00022 #include "ml_MultiLevelPreconditioner.h"
00023 using namespace ML_Epetra;
00024 #else
00025 #error blarf
00026 #endif
00027
00028 using namespace TSFExtended;
00029 using namespace Teuchos;
00030
00031
00032 AztecSolver::AztecSolver(const ParameterList& params)
00033 : LinearSolverBase<double>(params),
00034 options_(AZ_OPTIONS_SIZE),
00035 parameters_(AZ_PARAMS_SIZE),
00036 useML_(false),
00037 useIfpack_(false),
00038 useUserPrec_(false),
00039 aztec_recursive_iterate_(false),
00040 precParams_(),
00041 userPrec_(),
00042 aztec_status(AZ_STATUS_SIZE),
00043 aztec_proc_config(AZ_PROC_SIZE)
00044 {
00045 setName("AztecSolver");
00046 initParamMap();
00047
00048
00049 AZ_defaults((int*) &(options_[0]), (double*) &(parameters_[0]));
00050
00051
00052
00053 ParameterList::ConstIterator iter;
00054 for (iter=params.begin(); iter != params.end(); ++iter)
00055 {
00056 const std::string& name = params.name(iter);
00057 const ParameterEntry& entry = params.entry(iter);
00058
00059 if (entry.isList())
00060 {
00061 if (name=="Preconditioner")
00062 {
00063 precParams_ = params.sublist("Preconditioner");
00064 TEST_FOR_EXCEPTION(!precParams_.isParameter("Type"), std::runtime_error,
00065 "preconditioner type not specified in parameter list "
00066 << precParams_);
00067 if (precParams_.get<string>("Type") == "ML")
00068 {
00069 useML_ = true;
00070 }
00071 else if (precParams_.get<string>("Type") == "Ifpack")
00072 {
00073 useIfpack_ = true;
00074 }
00075 else if (precParams_.get<string>("Type") == "User")
00076 {
00077 useUserPrec_ = true;
00078 }
00079 continue;
00080 }
00081 }
00082
00083
00084 if (paramMap().find(name) == paramMap().end()) continue;
00085
00086
00087 int aztecCode = paramMap()[name];
00088
00089 if (name=="Recursive Iterate")
00090 {
00091 if (getValue<int>(entry))
00092 aztec_recursive_iterate_ = true;
00093 else
00094 aztec_recursive_iterate_ = false;
00095 continue;
00096 }
00097
00098
00099
00100
00101
00102
00103
00104
00105 if (entry.isType<string>())
00106 {
00107 std::string val = getValue<string>(entry);
00108 TEST_FOR_EXCEPTION(paramMap().find(val) == paramMap().end(),
00109 std::runtime_error,
00110 "Aztec solver ctor: [" << val << "] is not a "
00111 "valid Aztec option value");
00112 int optionVal = paramMap()[val];
00113 options_[aztecCode] = optionVal;
00114 }
00115 else if (entry.isType<int>())
00116 {
00117 int val = getValue<int>(entry);
00118 options_[aztecCode] = val;
00119 if (name=="Verbosity") setVerbosity(val);
00120 }
00121 else if (entry.isType<double>())
00122 {
00123 double val = getValue<double>(entry);
00124 parameters_[aztecCode] = val;
00125 }
00126 }
00127 }
00128
00129
00130 AztecSolver::AztecSolver(const Teuchos::map<int, int>& aztecOptions,
00131 const Teuchos::map<int, double>& aztecParameters)
00132 : LinearSolverBase<double>(ParameterList()),
00133 options_(AZ_OPTIONS_SIZE),
00134 parameters_(AZ_PARAMS_SIZE),
00135 useML_(false),
00136 useIfpack_(false),
00137 useUserPrec_(false),
00138 aztec_recursive_iterate_(false),
00139 precParams_(),
00140 userPrec_(),
00141 aztec_status(AZ_STATUS_SIZE),
00142 aztec_proc_config(AZ_PROC_SIZE)
00143 {
00144 setName("AztecSolver");
00145 if (aztecOptions.find(AZ_recursive_iterate) != aztecOptions.end())
00146 {
00147 aztec_recursive_iterate_ = true;
00148 }
00149
00150
00151 AZ_defaults((int*) &(options_[0]), (double*) &(parameters_[0]));
00152
00153
00154 map<int, int>::const_iterator opIter;
00155 for (opIter=aztecOptions.begin(); opIter!=aztecOptions.end(); opIter++)
00156 {
00157 int opKey = opIter->first;
00158 if (opKey==AZ_recursive_iterate) continue;
00159 int opValue = opIter->second;
00160 options_[opKey] = opValue;
00161 }
00162
00163
00164 map<int, double>::const_iterator parIter;
00165 for (parIter=aztecParameters.begin(); parIter!=aztecParameters.end();
00166 parIter++)
00167 {
00168 int parKey = parIter->first;
00169 double parValue = parIter->second;
00170 parameters_[parKey] = parValue;
00171 }
00172 }
00173
00174
00175 void AztecSolver::updateTolerance(const double& tol)
00176 {
00177 parameters_[AZ_tol] = tol;
00178 }
00179
00180 SolverState<double> AztecSolver::solve(const LinearOperator<double>& op,
00181 const Vector<double>& rhs,
00182 Vector<double>& soln) const
00183 {
00184 RCP<ostream> out;
00185 if (verb()==0)
00186 {
00187 out = rcp(new oblackholestream());
00188 }
00189 else
00190 {
00191 out = rcp(&Out::os(), false);
00192 }
00193 RCP<MultiLevelPreconditioner> mlPrec;
00194 RCP<Ifpack_Preconditioner> ifpackPrec;
00195 RCP<Epetra_Operator> prec;
00196
00197 TSFExtended::Vector<double> bCopy = rhs.copy();
00198 TSFExtended::Vector<double> xCopy = rhs.copy();
00199
00200 if (verb() > 4)
00201 {
00202 *out << "rhs=" << bCopy << std::endl;
00203 }
00204
00205 Epetra_Vector* b = EpetraVector::getConcretePtr(bCopy);
00206 Epetra_Vector* x = EpetraVector::getConcretePtr(xCopy);
00207
00208 Epetra_CrsMatrix& A = EpetraMatrix::getConcrete(op);
00209
00210 AztecOO aztec(&A, x, b);
00211
00212 aztec.SetAllAztecOptions((int*) &(options_[0]));
00213 aztec.SetAllAztecParams((double*) &(parameters_[0]));
00214 aztec.SetOutputStream(*out);
00215 aztec.SetErrorStream(*out);
00216
00217 int maxIters = options_[AZ_max_iter];
00218 double tol = parameters_[AZ_tol];
00219
00220 if (useML_)
00221 {
00222 std::string precType = precParams_.get<string>("Problem Type");
00223 ParameterList mlParams;
00224 ML_Epetra::SetDefaults(precType, mlParams);
00225
00226
00227
00228 ParameterList::ConstIterator iter;
00229 ParameterList mlSettings = precParams_.sublist("ML Settings");
00230 for (iter=mlSettings.begin(); iter!=mlSettings.end(); ++iter)
00231 {
00232 const std::string& name = mlSettings.name(iter);
00233 const ParameterEntry& entry = mlSettings.entry(iter);
00234 mlParams.setEntry(name, entry);
00235 }
00236
00237 mlPrec = rcp(new ML_Epetra::MultiLevelPreconditioner(A, mlParams));
00238 prec = rcp_dynamic_cast<Epetra_Operator>(mlPrec);
00239 }
00240 else if (useIfpack_)
00241 {
00242 Ifpack precFactory;
00243 int overlap = precParams_.get<int>("Overlap");
00244 std::string precType = precParams_.get<string>("Prec Type");
00245
00246 ParameterList ifpackParams = precParams_.sublist("Ifpack Settings");
00247
00248 ifpackPrec = rcp(precFactory.Create(precType, &A, overlap));
00249 prec = rcp_dynamic_cast<Epetra_Operator>(ifpackPrec);
00250 ifpackPrec->SetParameters(ifpackParams);
00251 ifpackPrec->Initialize();
00252 ifpackPrec->Compute();
00253 }
00254 else if (useUserPrec_)
00255 {
00256 TEST_FOR_EXCEPT(userPrec_.get() == 0);
00257 prec = userPrec_;
00258 }
00259
00260
00261 if (prec.get() != 0) aztec.SetPrecOperator(prec.get());
00262
00263 aztec.CheckInput();
00264
00265
00266
00267 if (aztec_recursive_iterate_)
00268 {
00269 aztec.recursiveIterate(maxIters, tol);
00270 }
00271 else
00272 {
00273 aztec.Iterate(maxIters, tol);
00274 }
00275
00276 soln = xCopy;
00277
00278 const double* status = aztec.GetAztecStatus();
00279 SolverStatusCode state = SolveCrashed;
00280
00281 std::string msg;
00282 switch((int) status[AZ_why])
00283 {
00284 case AZ_normal:
00285 state = SolveConverged;
00286 msg = "converged";
00287 break;
00288 case AZ_param:
00289 state = SolveCrashed;
00290 msg = "failed: parameter not available";
00291 break;
00292 case AZ_breakdown:
00293 state = SolveCrashed;
00294 msg = "failed: numerical breakdown";
00295 break;
00296 case AZ_loss:
00297 state = SolveCrashed;
00298 msg = "failed: numerical loss of precision";
00299 break;
00300 case AZ_ill_cond:
00301 state = SolveCrashed;
00302 msg = "failed: ill-conditioned Hessenberg matrix in GMRES";
00303 break;
00304 case AZ_maxits:
00305 state = SolveFailedToConverge;
00306 msg = "failed: maxiters reached without converged";
00307 break;
00308 }
00309 SolverState<double> rtn(state, "Aztec solver " + msg, (int) status[AZ_its],
00310 status[AZ_r]);
00311 return rtn;
00312 }
00313
00314
00315 void AztecSolver::setUserPrec(const LinearOperator<double>& P,
00316 const LinearSolver<double>& pSolver)
00317 {
00318 if (useUserPrec_)
00319 {
00320 userPrec_ = rcp(new Epetra::Epetra_TSFOperator(P, pSolver));
00321 }
00322 else
00323 {
00324 TEST_FOR_EXCEPTION(!useUserPrec_, std::runtime_error,
00325 "Attempt to set user-defined preconditioner "
00326 "after another preconditioner has been specified");
00327 }
00328 }
00329
00330 void AztecSolver::initParamMap()
00331 {
00332 static bool first = true;
00333 if (first)
00334 {
00335 paramMap()["Method"]=AZ_solver;
00336 paramMap()["CG"]=AZ_cg;
00337 paramMap()["GMRES"]=AZ_gmres;
00338 paramMap()["CGS"]=AZ_cgs;
00339 paramMap()["TFQMR"]=AZ_tfqmr;
00340 paramMap()["BICGSTAB"]=AZ_bicgstab;
00341 paramMap()["Direct"]=AZ_lu;
00342 paramMap()["Precond"]=AZ_precond;
00343 paramMap()["None"]=AZ_none;
00344 paramMap()["Jacobi"]=AZ_Jacobi;
00345 paramMap()["Neumann Series"]=AZ_Neumann;
00346 paramMap()["Symmetric Gauss-Seidel"]=AZ_sym_GS;
00347 paramMap()["Least-Squares Polynomial"]=AZ_ls;
00348 paramMap()["Recursive Iterate"]=AZ_recursive_iterate;
00349 paramMap()["Domain Decomposition"]=AZ_dom_decomp;
00350 paramMap()["Subdomain Solver"]=AZ_subdomain_solve;
00351 paramMap()["Approximate Sparse LU"]=AZ_lu;
00352 paramMap()["Saad ILUT"]=AZ_ilut;
00353 paramMap()["ILU"]=AZ_ilu;
00354 paramMap()["RILU"]=AZ_rilu;
00355 paramMap()["Block ILU"]=AZ_bilu;
00356 paramMap()["Incomplete Cholesky"]=AZ_icc;
00357 paramMap()["Residual Scaling"]=AZ_conv;
00358 paramMap()["Initial"]=AZ_r0;
00359 paramMap()["RHS"]=AZ_rhs;
00360 paramMap()["Matrix"]=AZ_Anorm;
00361 paramMap()["Solution"]=AZ_sol;
00362 paramMap()["No Scaling"]=AZ_noscaled;
00363 paramMap()["Verbosity"]=AZ_output;
00364 paramMap()["All"]=AZ_all;
00365 paramMap()["Silent"]=AZ_none;
00366 paramMap()["Warnings"]=AZ_warnings;
00367 paramMap()["Final Residual"]=AZ_last;
00368 paramMap()["Graph Fill"]=AZ_graph_fill;
00369 paramMap()["Max Iterations"]=AZ_max_iter;
00370 paramMap()["Polynomial Order"]=AZ_poly_ord;
00371 paramMap()["Overlap"]=AZ_overlap;
00372 paramMap()["Overlap Type"]=AZ_type_overlap;
00373 paramMap()["Standard"]=AZ_standard;
00374 paramMap()["Symmetric"]=AZ_symmetric;
00375 paramMap()["Restart Size"]=AZ_kspace;
00376 paramMap()["Reorder ILU"]=AZ_reorder;
00377 paramMap()["Keep Factorization"]=AZ_keep_info;
00378 paramMap()["GMRES Orthogonalization"]=AZ_orthog;
00379 paramMap()["Classical Gram-Schmidt"]=AZ_classic;
00380 paramMap()["Modified Gram-Schmidt"]=AZ_modified;
00381 paramMap()["Auxiliary Vector"]=AZ_aux_vec;
00382 paramMap()["Residual"]=AZ_resid;
00383 paramMap()["Random"]=AZ_rand;
00384 paramMap()["Tolerance"]=AZ_tol;
00385 paramMap()["Drop Tolerance"]=AZ_drop;
00386 paramMap()["Fill Ratio"]=AZ_ilut_fill;
00387 paramMap()["Damping"]=AZ_omega;
00388
00389 first = false;
00390 }
00391 }
00392
00393