|
Stratimikos Package Browser (Single Doxygen Collection)
Version of the Day
|
00001 /* 00002 // @HEADER 00003 // *********************************************************************** 00004 // 00005 // Stratimikos: Thyra-based strategies for linear solvers 00006 // Copyright (2006) Sandia Corporation 00007 // 00008 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive 00009 // license for use of this work by or on behalf of the U.S. Government. 00010 // 00011 // Redistribution and use in source and binary forms, with or without 00012 // modification, are permitted provided that the following conditions are 00013 // met: 00014 // 00015 // 1. Redistributions of source code must retain the above copyright 00016 // notice, this list of conditions and the following disclaimer. 00017 // 00018 // 2. Redistributions in binary form must reproduce the above copyright 00019 // notice, this list of conditions and the following disclaimer in the 00020 // documentation and/or other materials provided with the distribution. 00021 // 00022 // 3. Neither the name of the Corporation nor the names of the 00023 // contributors may be used to endorse or promote products derived from 00024 // this software without specific prior written permission. 00025 // 00026 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00027 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00028 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00029 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00030 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00031 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00032 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00033 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00034 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00035 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00036 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00037 // 00038 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov) 00039 // 00040 // *********************************************************************** 00041 // @HEADER 00042 */ 00043 00044 00045 #ifndef THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP 00046 #define THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP 00047 00048 #include "Thyra_BelosLinearOpWithSolve_decl.hpp" 00049 #include "Thyra_GeneralSolveCriteriaBelosStatusTest.hpp" 00050 #include "Thyra_LinearOpWithSolveHelpers.hpp" 00051 #include "Teuchos_DebugDefaultAsserts.hpp" 00052 #include "Teuchos_Assert.hpp" 00053 #include "Teuchos_TimeMonitor.hpp" 00054 #include "Teuchos_TypeTraits.hpp" 00055 00056 namespace { 00057 // Set the Belos solver's parameter list to scale its residual norms 00058 // in the specified way. 00059 // 00060 // We break this out in a separate function because the parameters 00061 // to set depend on which parameters the Belos solver supports. Not 00062 // all Belos solvers support both the "Implicit Residual Scaling" 00063 // and "Explicit Residual Scaling" parameters, so we have to check 00064 // the solver's list of valid parameters for the existence of these. 00065 // 00066 // Scaling options: Belos lets you decide whether the solver will 00067 // scale residual norms by the (left-)preconditioned initial 00068 // residual norms (residualScalingType = "Norm of Initial 00069 // Residual"), or by the unpreconditioned initial residual norms 00070 // (residualScalingType = "Norm of RHS"). Usually you want to scale 00071 // by the unpreconditioned initial residual norms. This is because 00072 // preconditioning is just an optimization, and you really want to 00073 // make ||B - AX|| small, rather than ||M B - M (A X)||. If you're 00074 // measuring ||B - AX|| and scaling by the initial residual, you 00075 // should use the unpreconditioned initial residual to match it. 00076 // 00077 // Note, however, that the implicit residual test computes 00078 // left-preconditioned residuals, if a left preconditioner was 00079 // provided. That's OK because when Belos solvers (at least the 00080 // GMRES variants) are given a left preconditioner, they first check 00081 // the implicit residuals. If those converge, they then check the 00082 // explicit residuals. The explicit residual test does _not_ apply 00083 // the left preconditioner when computing the residual. The 00084 // implicit residual test is just an optimization so that Belos 00085 // doesn't have to compute explicit residuals B - A*X at every 00086 // iteration. This is why we use the same scaling factor for both 00087 // the implicit and explicit residuals. 00088 // 00089 // Arguments: 00090 // 00091 // solverParams [in/out] Parameters for the current solve. 00092 // 00093 // solverValidParams [in] Valid parameter list for the Belos solver. 00094 // Result of calling the solver's getValidParameters() method. 00095 // 00096 // residualScalingType [in] String describing how the solver should 00097 // scale residuals. Valid values include "Norm of RHS" and "Norm 00098 // of Initial Residual" (these are the only two options this file 00099 // currently uses, though Belos offers other options). 00100 void 00101 setResidualScalingType (const Teuchos::RCP<Teuchos::ParameterList>& solverParams, 00102 const Teuchos::RCP<const Teuchos::ParameterList>& solverValidParams, 00103 const std::string& residualScalingType) 00104 { 00105 // Many Belos solvers (especially the GMRES variants) define both 00106 // "Implicit Residual Scaling" and "Explicit Residual Scaling" 00107 // options. 00108 // 00109 // "Implicit" means "the left-preconditioned approximate 00110 // a.k.a. 'recursive' residual as computed by the Krylov method." 00111 // 00112 // "Explicit" means ||B - A*X||, the unpreconditioned, "exact" 00113 // residual. 00114 // 00115 // Belos' GMRES implementations chain these two tests in sequence. 00116 // Implicit comes first, and explicit is not evaluated unless 00117 // implicit passes. In some cases (e.g., no left preconditioner), 00118 // GMRES _only_ uses the implicit tests. This means that only 00119 // setting "Explicit Residual Scaling" won't change the solver's 00120 // behavior. Stratimikos tends to prefer using a right 00121 // preconditioner, in which case setting only the "Explicit 00122 // Residual Scaling" argument has no effect. Furthermore, if 00123 // "Explicit Residual Scaling" is set to something other than the 00124 // default (initial residual norm), without "Implicit Residual 00125 // Scaling" getting the same setting, then the implicit residual 00126 // test will be using a radically different scaling factor than 00127 // the user wanted. 00128 // 00129 // Not all Belos solvers support both options. We check the 00130 // solver's valid parameter list first before attempting to set 00131 // the option. 00132 if (solverValidParams->isParameter ("Implicit Residual Scaling")) { 00133 solverParams->set ("Implicit Residual Scaling", residualScalingType); 00134 } 00135 if (solverValidParams->isParameter ("Explicit Residual Scaling")) { 00136 solverParams->set ("Explicit Residual Scaling", residualScalingType); 00137 } 00138 } 00139 00140 } // namespace (anonymous) 00141 00142 00143 namespace Thyra { 00144 00145 00146 // Constructors/initializers/accessors 00147 00148 00149 template<class Scalar> 00150 BelosLinearOpWithSolve<Scalar>::BelosLinearOpWithSolve() 00151 :convergenceTestFrequency_(-1), 00152 isExternalPrec_(false), 00153 supportSolveUse_(SUPPORT_SOLVE_UNSPECIFIED), 00154 defaultTol_ (-1.0) 00155 {} 00156 00157 00158 template<class Scalar> 00159 void BelosLinearOpWithSolve<Scalar>::initialize( 00160 const RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > &lp, 00161 const RCP<Teuchos::ParameterList> &solverPL, 00162 const RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > &iterativeSolver, 00163 const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc, 00164 const RCP<const PreconditionerBase<Scalar> > &prec, 00165 const bool isExternalPrec_in, 00166 const RCP<const LinearOpSourceBase<Scalar> > &approxFwdOpSrc, 00167 const ESupportSolveUse &supportSolveUse_in, 00168 const int convergenceTestFrequency 00169 ) 00170 { 00171 using Teuchos::as; 00172 using Teuchos::TypeNameTraits; 00173 using Teuchos::Exceptions::InvalidParameterType; 00174 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitude_type; 00175 00176 this->setLinePrefix("BELOS/T"); 00177 // ToDo: Validate input 00178 lp_ = lp; 00179 solverPL_ = solverPL; 00180 iterativeSolver_ = iterativeSolver; 00181 fwdOpSrc_ = fwdOpSrc; 00182 prec_ = prec; 00183 isExternalPrec_ = isExternalPrec_in; 00184 approxFwdOpSrc_ = approxFwdOpSrc; 00185 supportSolveUse_ = supportSolveUse_in; 00186 convergenceTestFrequency_ = convergenceTestFrequency; 00187 // Check if "Convergence Tolerance" is in the solver parameter list. If 00188 // not, use the default from the solver. 00189 if ( !is_null(solverPL_) ) { 00190 if (solverPL_->isParameter("Convergence Tolerance")) { 00191 00192 // Stratimikos prefers tolerances as double, no matter the 00193 // Scalar type. However, we also want it to accept the 00194 // tolerance as magnitude_type, for example float if the Scalar 00195 // type is float or std::complex<float>. 00196 if (solverPL_->isType<double> ("Convergence Tolerance")) { 00197 defaultTol_ = 00198 as<magnitude_type> (solverPL_->get<double> ("Convergence Tolerance")); 00199 } 00200 else if (Teuchos::TypeTraits::is_same<double, magnitude_type>::value) { 00201 // magnitude_type == double in this case, and we've already 00202 // checked double above. 00203 TEUCHOS_TEST_FOR_EXCEPTION( 00204 true, std::invalid_argument, "BelosLinearOpWithSolve::initialize: " 00205 "The \"Convergence Tolerance\" parameter, which you provided, must " 00206 "have type double (the type of the magnitude of Scalar = double)."); 00207 } 00208 else if (solverPL_->isType<magnitude_type> ("Convergence Tolerance")) { 00209 defaultTol_ = solverPL_->get<magnitude_type> ("Convergence Tolerance"); 00210 } 00211 else { 00212 // Throwing InvalidParameterType ensures that the exception's 00213 // type is consistent both with what this method would have 00214 // thrown before for an unrecognized type, and with what the 00215 // user expects in general when the parameter doesn't have the 00216 // right type. 00217 TEUCHOS_TEST_FOR_EXCEPTION( 00218 true, InvalidParameterType, "BelosLinearOpWithSolve::initialize: " 00219 "The \"Convergence Tolerance\" parameter, which you provided, must " 00220 "have type double (preferred) or the type of the magnitude of Scalar " 00221 "= " << TypeNameTraits<Scalar>::name () << ", which is " << 00222 TypeNameTraits<magnitude_type>::name () << " in this case. You can " 00223 "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType."); 00224 } 00225 } 00226 } 00227 else { 00228 RCP<const Teuchos::ParameterList> defaultPL = 00229 iterativeSolver->getValidParameters(); 00230 00231 // Stratimikos prefers tolerances as double, no matter the 00232 // Scalar type. However, we also want it to accept the 00233 // tolerance as magnitude_type, for example float if the Scalar 00234 // type is float or std::complex<float>. 00235 if (defaultPL->isType<double> ("Convergence Tolerance")) { 00236 defaultTol_ = 00237 as<magnitude_type> (defaultPL->get<double> ("Convergence Tolerance")); 00238 } 00239 else if (Teuchos::TypeTraits::is_same<double, magnitude_type>::value) { 00240 // magnitude_type == double in this case, and we've already 00241 // checked double above. 00242 TEUCHOS_TEST_FOR_EXCEPTION( 00243 true, std::invalid_argument, "BelosLinearOpWithSolve::initialize: " 00244 "The \"Convergence Tolerance\" parameter, which you provided, must " 00245 "have type double (the type of the magnitude of Scalar = double)."); 00246 } 00247 else if (defaultPL->isType<magnitude_type> ("Convergence Tolerance")) { 00248 defaultTol_ = defaultPL->get<magnitude_type> ("Convergence Tolerance"); 00249 } 00250 else { 00251 // Throwing InvalidParameterType ensures that the exception's 00252 // type is consistent both with what this method would have 00253 // thrown before for an unrecognized type, and with what the 00254 // user expects in general when the parameter doesn't have the 00255 // right type. 00256 TEUCHOS_TEST_FOR_EXCEPTION( 00257 true, InvalidParameterType, "BelosLinearOpWithSolve::initialize: " 00258 "The \"Convergence Tolerance\" parameter, which you provided, must " 00259 "have type double (preferred) or the type of the magnitude of Scalar " 00260 "= " << TypeNameTraits<Scalar>::name () << ", which is " << 00261 TypeNameTraits<magnitude_type>::name () << " in this case. You can " 00262 "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType."); 00263 } 00264 } 00265 } 00266 00267 00268 template<class Scalar> 00269 RCP<const LinearOpSourceBase<Scalar> > 00270 BelosLinearOpWithSolve<Scalar>::extract_fwdOpSrc() 00271 { 00272 RCP<const LinearOpSourceBase<Scalar> > 00273 _fwdOpSrc = fwdOpSrc_; 00274 fwdOpSrc_ = Teuchos::null; 00275 return _fwdOpSrc; 00276 } 00277 00278 00279 template<class Scalar> 00280 RCP<const PreconditionerBase<Scalar> > 00281 BelosLinearOpWithSolve<Scalar>::extract_prec() 00282 { 00283 RCP<const PreconditionerBase<Scalar> > 00284 _prec = prec_; 00285 prec_ = Teuchos::null; 00286 return _prec; 00287 } 00288 00289 00290 template<class Scalar> 00291 bool BelosLinearOpWithSolve<Scalar>::isExternalPrec() const 00292 { 00293 return isExternalPrec_; 00294 } 00295 00296 00297 template<class Scalar> 00298 RCP<const LinearOpSourceBase<Scalar> > 00299 BelosLinearOpWithSolve<Scalar>::extract_approxFwdOpSrc() 00300 { 00301 RCP<const LinearOpSourceBase<Scalar> > 00302 _approxFwdOpSrc = approxFwdOpSrc_; 00303 approxFwdOpSrc_ = Teuchos::null; 00304 return _approxFwdOpSrc; 00305 } 00306 00307 00308 template<class Scalar> 00309 ESupportSolveUse BelosLinearOpWithSolve<Scalar>::supportSolveUse() const 00310 { 00311 return supportSolveUse_; 00312 } 00313 00314 00315 template<class Scalar> 00316 void BelosLinearOpWithSolve<Scalar>::uninitialize( 00317 RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > *lp, 00318 RCP<Teuchos::ParameterList> *solverPL, 00319 RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > *iterativeSolver, 00320 RCP<const LinearOpSourceBase<Scalar> > *fwdOpSrc, 00321 RCP<const PreconditionerBase<Scalar> > *prec, 00322 bool *isExternalPrec_in, 00323 RCP<const LinearOpSourceBase<Scalar> > *approxFwdOpSrc, 00324 ESupportSolveUse *supportSolveUse_in 00325 ) 00326 { 00327 if (lp) *lp = lp_; 00328 if (solverPL) *solverPL = solverPL_; 00329 if (iterativeSolver) *iterativeSolver = iterativeSolver_; 00330 if (fwdOpSrc) *fwdOpSrc = fwdOpSrc_; 00331 if (prec) *prec = prec_; 00332 if (isExternalPrec_in) *isExternalPrec_in = isExternalPrec_; 00333 if (approxFwdOpSrc) *approxFwdOpSrc = approxFwdOpSrc_; 00334 if (supportSolveUse_in) *supportSolveUse_in = supportSolveUse_; 00335 00336 lp_ = Teuchos::null; 00337 solverPL_ = Teuchos::null; 00338 iterativeSolver_ = Teuchos::null; 00339 fwdOpSrc_ = Teuchos::null; 00340 prec_ = Teuchos::null; 00341 isExternalPrec_ = false; 00342 approxFwdOpSrc_ = Teuchos::null; 00343 supportSolveUse_ = SUPPORT_SOLVE_UNSPECIFIED; 00344 } 00345 00346 00347 // Overridden from LinearOpBase 00348 00349 00350 template<class Scalar> 00351 RCP< const VectorSpaceBase<Scalar> > 00352 BelosLinearOpWithSolve<Scalar>::range() const 00353 { 00354 if (!is_null(lp_)) 00355 return lp_->getOperator()->range(); 00356 return Teuchos::null; 00357 } 00358 00359 00360 template<class Scalar> 00361 RCP< const VectorSpaceBase<Scalar> > 00362 BelosLinearOpWithSolve<Scalar>::domain() const 00363 { 00364 if (!is_null(lp_)) 00365 return lp_->getOperator()->domain(); 00366 return Teuchos::null; 00367 } 00368 00369 00370 template<class Scalar> 00371 RCP<const LinearOpBase<Scalar> > 00372 BelosLinearOpWithSolve<Scalar>::clone() const 00373 { 00374 return Teuchos::null; // Not supported yet but could be 00375 } 00376 00377 00378 // Overridden from Teuchos::Describable 00379 00380 00381 template<class Scalar> 00382 std::string BelosLinearOpWithSolve<Scalar>::description() const 00383 { 00384 std::ostringstream oss; 00385 oss << Teuchos::Describable::description(); 00386 if ( !is_null(lp_) && !is_null(lp_->getOperator()) ) { 00387 oss << "{"; 00388 oss << "iterativeSolver=\'"<<iterativeSolver_->description()<<"\'"; 00389 oss << ",fwdOp=\'"<<lp_->getOperator()->description()<<"\'"; 00390 if (lp_->getLeftPrec().get()) 00391 oss << ",leftPrecOp=\'"<<lp_->getLeftPrec()->description()<<"\'"; 00392 if (lp_->getRightPrec().get()) 00393 oss << ",rightPrecOp=\'"<<lp_->getRightPrec()->description()<<"\'"; 00394 oss << "}"; 00395 } 00396 // ToDo: Make Belos::SolverManager derive from Teuchos::Describable so 00397 // that we can get better information. 00398 return oss.str(); 00399 } 00400 00401 00402 template<class Scalar> 00403 void BelosLinearOpWithSolve<Scalar>::describe( 00404 Teuchos::FancyOStream &out_arg, 00405 const Teuchos::EVerbosityLevel verbLevel 00406 ) const 00407 { 00408 using Teuchos::FancyOStream; 00409 using Teuchos::OSTab; 00410 using Teuchos::describe; 00411 RCP<FancyOStream> out = rcp(&out_arg,false); 00412 OSTab tab(out); 00413 switch (verbLevel) { 00414 case Teuchos::VERB_LOW: 00415 break; 00416 case Teuchos::VERB_DEFAULT: 00417 case Teuchos::VERB_MEDIUM: 00418 *out << this->description() << std::endl; 00419 break; 00420 case Teuchos::VERB_HIGH: 00421 case Teuchos::VERB_EXTREME: 00422 { 00423 *out 00424 << Teuchos::Describable::description()<< "{" 00425 << "rangeDim=" << this->range()->dim() 00426 << ",domainDim=" << this->domain()->dim() << "}\n"; 00427 if (lp_->getOperator().get()) { 00428 OSTab tab1(out); 00429 *out 00430 << "iterativeSolver = "<<describe(*iterativeSolver_,verbLevel) 00431 << "fwdOp = " << describe(*lp_->getOperator(),verbLevel); 00432 if (lp_->getLeftPrec().get()) 00433 *out << "leftPrecOp = "<<describe(*lp_->getLeftPrec(),verbLevel); 00434 if (lp_->getRightPrec().get()) 00435 *out << "rightPrecOp = "<<describe(*lp_->getRightPrec(),verbLevel); 00436 } 00437 break; 00438 } 00439 default: 00440 TEUCHOS_TEST_FOR_EXCEPT(true); // Should never get here! 00441 } 00442 } 00443 00444 00445 // protected 00446 00447 00448 // Overridden from LinearOpBase 00449 00450 00451 template<class Scalar> 00452 bool BelosLinearOpWithSolve<Scalar>::opSupportedImpl(EOpTransp M_trans) const 00453 { 00454 return ::Thyra::opSupported(*lp_->getOperator(),M_trans); 00455 } 00456 00457 00458 template<class Scalar> 00459 void BelosLinearOpWithSolve<Scalar>::applyImpl( 00460 const EOpTransp M_trans, 00461 const MultiVectorBase<Scalar> &X, 00462 const Ptr<MultiVectorBase<Scalar> > &Y, 00463 const Scalar alpha, 00464 const Scalar beta 00465 ) const 00466 { 00467 ::Thyra::apply<Scalar>(*lp_->getOperator(), M_trans, X, Y, alpha, beta); 00468 } 00469 00470 00471 // Overridden from LinearOpWithSolveBase 00472 00473 00474 template<class Scalar> 00475 bool 00476 BelosLinearOpWithSolve<Scalar>::solveSupportsImpl(EOpTransp M_trans) const 00477 { 00478 return solveSupportsNewImpl(M_trans, Teuchos::null); 00479 } 00480 00481 00482 template<class Scalar> 00483 bool 00484 BelosLinearOpWithSolve<Scalar>::solveSupportsNewImpl(EOpTransp transp, 00485 const Ptr<const SolveCriteria<Scalar> > solveCriteria) const 00486 { 00487 // Only support forward solve right now! 00488 if (real_trans(transp)==NOTRANS) return true; 00489 return false; // ToDo: Support adjoint solves! 00490 // Otherwise, Thyra/Belos now supports every solve criteria type that exists 00491 // because of the class Thyra::GeneralSolveCriteriaBelosStatusTest! 00492 return true; 00493 /* 00494 if (real_trans(M_trans)==NOTRANS) { 00495 return ( 00496 solveMeasureType.useDefault() 00497 || 00498 solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL,SOLVE_MEASURE_NORM_RHS) 00499 || 00500 solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL,SOLVE_MEASURE_NORM_INIT_RESIDUAL) 00501 ); 00502 } 00503 */ 00504 } 00505 00506 00507 template<class Scalar> 00508 bool 00509 BelosLinearOpWithSolve<Scalar>::solveSupportsSolveMeasureTypeImpl( 00510 EOpTransp M_trans, const SolveMeasureType& solveMeasureType) const 00511 { 00512 SolveCriteria<Scalar> solveCriteria(solveMeasureType, SolveCriteria<Scalar>::unspecifiedTolerance()); 00513 return solveSupportsNewImpl(M_trans, Teuchos::constOptInArg(solveCriteria)); 00514 } 00515 00516 00517 template<class Scalar> 00518 SolveStatus<Scalar> 00519 BelosLinearOpWithSolve<Scalar>::solveImpl( 00520 const EOpTransp M_trans, 00521 const MultiVectorBase<Scalar> &B, 00522 const Ptr<MultiVectorBase<Scalar> > &X, 00523 const Ptr<const SolveCriteria<Scalar> > solveCriteria 00524 ) const 00525 { 00526 00527 THYRA_FUNC_TIME_MONITOR("Stratimikos: BelosLOWS"); 00528 00529 using Teuchos::rcp; 00530 using Teuchos::rcpFromRef; 00531 using Teuchos::rcpFromPtr; 00532 using Teuchos::FancyOStream; 00533 using Teuchos::OSTab; 00534 using Teuchos::ParameterList; 00535 using Teuchos::parameterList; 00536 using Teuchos::describe; 00537 typedef Teuchos::ScalarTraits<Scalar> ST; 00538 typedef typename ST::magnitudeType ScalarMag; 00539 Teuchos::Time totalTimer(""), timer(""); 00540 totalTimer.start(true); 00541 00542 assertSolveSupports(*this, M_trans, solveCriteria); 00543 // 2010/08/22: rabartl: Bug 4915 ToDo: Move the above into the NIV function 00544 // solve(...). 00545 00546 const RCP<FancyOStream> out = this->getOStream(); 00547 const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel(); 00548 OSTab tab = this->getOSTab(); 00549 if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW)) { 00550 *out << "\nStarting iterations with Belos:\n"; 00551 OSTab tab2(out); 00552 *out << "Using forward operator = " << describe(*fwdOpSrc_->getOp(),verbLevel); 00553 *out << "Using iterative solver = " << describe(*iterativeSolver_,verbLevel); 00554 *out << "With #Eqns="<<B.range()->dim()<<", #RHSs="<<B.domain()->dim()<<" ...\n"; 00555 } 00556 00557 // 00558 // Set RHS and LHS 00559 // 00560 00561 bool ret = lp_->setProblem( rcpFromPtr(X), rcpFromRef(B) ); 00562 TEUCHOS_TEST_FOR_EXCEPTION( 00563 ret == false, CatastrophicSolveFailure 00564 ,"Error, the Belos::LinearProblem could not be set for the current solve!" 00565 ); 00566 00567 // 00568 // Set the solution criteria 00569 // 00570 00571 // Parameter list for the current solve. 00572 const RCP<ParameterList> tmpPL = Teuchos::parameterList(); 00573 00574 // The solver's valid parameter list. 00575 RCP<const ParameterList> validPL = iterativeSolver_->getValidParameters(); 00576 00577 SolveMeasureType solveMeasureType; 00578 RCP<GeneralSolveCriteriaBelosStatusTest<Scalar> > generalSolveCriteriaBelosStatusTest; 00579 if (nonnull(solveCriteria)) { 00580 solveMeasureType = solveCriteria->solveMeasureType; 00581 const ScalarMag requestedTol = solveCriteria->requestedTol; 00582 if (solveMeasureType.useDefault()) { 00583 tmpPL->set("Convergence Tolerance", defaultTol_); 00584 } 00585 else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_RHS)) { 00586 if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) { 00587 tmpPL->set("Convergence Tolerance", requestedTol); 00588 } 00589 else { 00590 tmpPL->set("Convergence Tolerance", defaultTol_); 00591 } 00592 setResidualScalingType (tmpPL, validPL, "Norm of RHS"); 00593 } 00594 else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_INIT_RESIDUAL)) { 00595 if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) { 00596 tmpPL->set("Convergence Tolerance", requestedTol); 00597 } 00598 else { 00599 tmpPL->set("Convergence Tolerance", defaultTol_); 00600 } 00601 setResidualScalingType (tmpPL, validPL, "Norm of Initial Residual"); 00602 } 00603 else { 00604 // Set the most generic (and inefficient) solve criteria 00605 generalSolveCriteriaBelosStatusTest = createGeneralSolveCriteriaBelosStatusTest( 00606 *solveCriteria, convergenceTestFrequency_); 00607 // Set the verbosity level (one level down) 00608 generalSolveCriteriaBelosStatusTest->setOStream(out); 00609 generalSolveCriteriaBelosStatusTest->setVerbLevel(incrVerbLevel(verbLevel, -1)); 00610 // Set the default convergence tolerance to always converged to allow 00611 // the above status test to control things. 00612 tmpPL->set("Convergence Tolerance", 1.0); 00613 } 00614 // maximum iterations 00615 if (nonnull(solveCriteria->extraParameters)) { 00616 if (Teuchos::isParameterType<int>(*solveCriteria->extraParameters,"Maximum Iterations")) { 00617 tmpPL->set("Maximum Iterations", Teuchos::get<int>(*solveCriteria->extraParameters,"Maximum Iterations")); 00618 } 00619 } 00620 } 00621 else { 00622 // No solveCriteria was even passed in! 00623 tmpPL->set("Convergence Tolerance", defaultTol_); 00624 } 00625 00626 // 00627 // Solve the linear system 00628 // 00629 00630 Belos::ReturnType belosSolveStatus; 00631 { 00632 // Write detailed convergence information if requested for levels >= VERB_LOW 00633 RCP<std::ostream> 00634 outUsed = 00635 ( static_cast<int>(verbLevel) >= static_cast<int>(Teuchos::VERB_LOW) 00636 ? out 00637 : rcp(new FancyOStream(rcp(new Teuchos::oblackholestream()))) 00638 ); 00639 Teuchos::OSTab tab1(outUsed,1,"BELOS"); 00640 tmpPL->set("Output Stream", outUsed); 00641 iterativeSolver_->setParameters(tmpPL); 00642 if (nonnull(generalSolveCriteriaBelosStatusTest)) { 00643 iterativeSolver_->setUserConvStatusTest(generalSolveCriteriaBelosStatusTest); 00644 } 00645 try { 00646 belosSolveStatus = iterativeSolver_->solve(); 00647 } 00648 catch (Belos::BelosError&) { 00649 belosSolveStatus = Belos::Unconverged; 00650 } 00651 } 00652 00653 // 00654 // Report the solve status 00655 // 00656 00657 totalTimer.stop(); 00658 00659 SolveStatus<Scalar> solveStatus; 00660 00661 switch (belosSolveStatus) { 00662 case Belos::Unconverged: { 00663 solveStatus.solveStatus = SOLVE_STATUS_UNCONVERGED; 00664 // Set achievedTol even if the solver did not converge. This is 00665 // helpful for things like nonlinear solvers, which might be 00666 // able to use a partially converged result, and which would 00667 // like to know the achieved convergence tolerance for use in 00668 // computing bounds. It's also helpful for estimating whether a 00669 // small increase in the maximum iteration count might be 00670 // helpful next time. 00671 try { 00672 // Some solvers might not have implemented achievedTol(). 00673 // The default implementation throws std::runtime_error. 00674 solveStatus.achievedTol = iterativeSolver_->achievedTol(); 00675 } catch (std::runtime_error&) { 00676 // Do nothing; use the default value of achievedTol. 00677 } 00678 break; 00679 } 00680 case Belos::Converged: { 00681 solveStatus.solveStatus = SOLVE_STATUS_CONVERGED; 00682 if (nonnull(generalSolveCriteriaBelosStatusTest)) { 00683 // The user set a custom status test. This means that we 00684 // should ask the custom status test itself, rather than the 00685 // Belos solver, what the final achieved convergence tolerance 00686 // was. 00687 const ArrayView<const ScalarMag> achievedTol = 00688 generalSolveCriteriaBelosStatusTest->achievedTol(); 00689 solveStatus.achievedTol = ST::zero(); 00690 for (Ordinal i = 0; i < achievedTol.size(); ++i) { 00691 solveStatus.achievedTol = std::max(solveStatus.achievedTol, achievedTol[i]); 00692 } 00693 } 00694 else { 00695 try { 00696 // Some solvers might not have implemented achievedTol(). 00697 // The default implementation throws std::runtime_error. 00698 solveStatus.achievedTol = iterativeSolver_->achievedTol(); 00699 } catch (std::runtime_error&) { 00700 // Use the default convergence tolerance. This is a correct 00701 // upper bound, since we did actually converge. 00702 solveStatus.achievedTol = tmpPL->get("Convergence Tolerance", defaultTol_); 00703 } 00704 } 00705 break; 00706 } 00707 TEUCHOS_SWITCH_DEFAULT_DEBUG_ASSERT(); 00708 } 00709 00710 std::ostringstream ossmessage; 00711 ossmessage 00712 << "The Belos solver of type \""<<iterativeSolver_->description() 00713 <<"\" returned a solve status of \""<< toString(solveStatus.solveStatus) << "\"" 00714 << " in " << iterativeSolver_->getNumIters() << " iterations" 00715 << " with total CPU time of " << totalTimer.totalElapsedTime() << " sec" ; 00716 if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW)) 00717 *out << "\n" << ossmessage.str() << "\n"; 00718 00719 solveStatus.message = ossmessage.str(); 00720 00721 // Dump the getNumIters() and the achieved convergence tolerance 00722 // into solveStatus.extraParameters, as the "Belos/Iteration Count" 00723 // resp. "Belos/Achieved Tolerance" parameters. 00724 if (solveStatus.extraParameters.is_null()) { 00725 solveStatus.extraParameters = parameterList (); 00726 } 00727 solveStatus.extraParameters->set ("Belos/Iteration Count", 00728 iterativeSolver_->getNumIters());\ 00729 // package independent version of the same 00730 solveStatus.extraParameters->set ("Iteration Count", 00731 iterativeSolver_->getNumIters());\ 00732 // NOTE (mfh 13 Dec 2011) Though the most commonly used Belos 00733 // solvers do implement achievedTol(), some Belos solvers currently 00734 // do not. In the latter case, if the solver did not converge, the 00735 // reported achievedTol() value may just be the default "invalid" 00736 // value -1, and if the solver did converge, the reported value will 00737 // just be the convergence tolerance (a correct upper bound). 00738 solveStatus.extraParameters->set ("Belos/Achieved Tolerance", 00739 solveStatus.achievedTol); 00740 00741 // This information is in the previous line, which is printed anytime the verbosity 00742 // is not set to Teuchos::VERB_NONE, so I'm commenting this out for now. 00743 // if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_NONE)) 00744 // *out << "\nTotal solve time in Belos = "<<totalTimer.totalElapsedTime()<<" sec\n"; 00745 00746 return solveStatus; 00747 00748 } 00749 00750 00751 } // end namespace Thyra 00752 00753 00754 #endif // THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP
1.7.6.1