00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #ifndef BELOS_GCRODR_SOLMGR_HPP
00030 #define BELOS_GCRODR_SOLMGR_HPP
00031
00036 #include "BelosConfigDefs.hpp"
00037 #include "BelosTypes.hpp"
00038
00039 #include "BelosLinearProblem.hpp"
00040 #include "BelosSolverManager.hpp"
00041
00042 #include "BelosBlockGmresIter.hpp"
00043 #include "BelosGCRODRIter.hpp"
00044 #include "BelosBlockFGmresIter.hpp"
00045 #include "BelosDGKSOrthoManager.hpp"
00046 #include "BelosICGSOrthoManager.hpp"
00047 #include "BelosIMGSOrthoManager.hpp"
00048 #include "BelosStatusTestMaxIters.hpp"
00049 #include "BelosStatusTestGenResNorm.hpp"
00050 #include "BelosStatusTestCombo.hpp"
00051 #include "BelosStatusTestOutput.hpp"
00052 #include "BelosOutputManager.hpp"
00053 #include "Teuchos_BLAS.hpp"
00054 #include "Teuchos_LAPACK.hpp"
00055 #include "Teuchos_TimeMonitor.hpp"
00056
00072 namespace Belos {
00073
00075
00076
00083 class GCRODRSolMgrLinearProblemFailure : public BelosError {public:
00084 GCRODRSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
00085 {}};
00086
00093 class GCRODRSolMgrOrthoFailure : public BelosError {public:
00094 GCRODRSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
00095 {}};
00096
00103 class GCRODRSolMgrLAPACKFailure : public BelosError {public:
00104 GCRODRSolMgrLAPACKFailure(const std::string& what_arg) : BelosError(what_arg)
00105 {}};
00106
00113 class GCRODRSolMgrRecyclingFailure : public BelosError {public:
00114 GCRODRSolMgrRecyclingFailure(const std::string& what_arg) : BelosError(what_arg)
00115 {}};
00116
00118
00119
00120 template<class ScalarType, class MV, class OP>
00121 class GCRODRSolMgr : public SolverManager<ScalarType,MV,OP> {
00122
00123 private:
00124 typedef MultiVecTraits<ScalarType,MV> MVT;
00125 typedef OperatorTraits<ScalarType,MV,OP> OPT;
00126 typedef Teuchos::ScalarTraits<ScalarType> SCT;
00127 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
00128 typedef Teuchos::ScalarTraits<MagnitudeType> MT;
00129
00130 public:
00131
00133
00134
00140 GCRODRSolMgr();
00141
00154 GCRODRSolMgr( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
00155 const Teuchos::RCP<Teuchos::ParameterList> &pl );
00156
00158 virtual ~GCRODRSolMgr() {};
00160
00162
00163
00166 const LinearProblem<ScalarType,MV,OP>& getProblem() const {
00167 return *problem_;
00168 }
00169
00172 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
00173
00176 Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const { return params_; }
00177
00183 Teuchos::Array<Teuchos::RCP<Teuchos::Time> > getTimers() const {
00184 return tuple(timerSolve_);
00185 }
00186
00189 bool isLOADetected() const { return false; }
00190
00192
00194
00195
00196 void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) { problem_ = problem; }
00197
00198 void setParameters( const Teuchos::RCP<Teuchos::ParameterList> ¶ms );
00199
00201
00203
00204
00222 ReturnType solve();
00223
00225
00228
00230 std::string description() const;
00231
00233
00234 private:
00235
00236
00237
00238
00239
00240 int getHarmonicVecs1(int m,
00241 const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
00242 Teuchos::SerialDenseMatrix<int,ScalarType>& PP);
00243
00244
00245
00246
00247
00248
00249 int getHarmonicVecs2(int keff, int m,
00250 const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
00251 const Teuchos::RCP<const MV>& VV,
00252 Teuchos::SerialDenseMatrix<int,ScalarType>& PP);
00253
00254
00255 void sort(std::vector<ScalarType>& dlist, int n, std::vector<int>& iperm);
00256
00257
00258 Belos::ScaleType convertStringToScaleType( string& scaleType ) {
00259 if (scaleType == "Norm of Initial Residual") {
00260 return Belos::NormOfInitRes;
00261 } else if (scaleType == "Norm of Preconditioned Initial Residual") {
00262 return Belos::NormOfPrecInitRes;
00263 } else if (scaleType == "Norm of RHS") {
00264 return Belos::NormOfRHS;
00265 } else if (scaleType == "None") {
00266 return Belos::None;
00267 } else
00268 TEST_FOR_EXCEPTION( true ,std::logic_error,
00269 "Belos::GCRODRSolMgr(): Invalid residual scaling type.");
00270 }
00271
00272
00273 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
00274
00275
00276 Teuchos::RCP<OutputManager<ScalarType> > printer_;
00277 Teuchos::RCP<ostream> outputStream_;
00278
00279
00280 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
00281 Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
00282 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > convTest_;
00283 Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > expConvTest_, impConvTest_;
00284 Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
00285
00286
00287 Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
00288
00289
00290 Teuchos::RCP<ParameterList> params_;
00291
00292
00293 static const MagnitudeType convtol_default_;
00294 static const MagnitudeType orthoKappa_default_;
00295 static const int maxRestarts_default_;
00296 static const int maxIters_default_;
00297 static const int numBlocks_default_;
00298 static const int recycledBlocks_default_;
00299 static const int verbosity_default_;
00300 static const int outputFreq_default_;
00301 static const std::string impResScale_default_;
00302 static const std::string expResScale_default_;
00303 static const std::string label_default_;
00304 static const std::string orthoType_default_;
00305 static const Teuchos::RCP<ostream> outputStream_default_;
00306
00307
00308 MagnitudeType convtol_, orthoKappa_;
00309 int maxRestarts_, maxIters_;
00310 int numBlocks_, recycledBlocks_, verbosity_, outputFreq_;
00311 std::string orthoType_;
00312 std::string impResScale_, expResScale_;
00313
00314
00315 Teuchos::RCP<MV> U_, C_, r_;
00316
00317
00318 std::string label_;
00319 Teuchos::RCP<Teuchos::Time> timerSolve_;
00320
00321
00322 bool isSet_;
00323 };
00324
00325
00326
00327 template<class ScalarType, class MV, class OP>
00328 const typename Teuchos::ScalarTraits<ScalarType>::magnitudeType GCRODRSolMgr<ScalarType,MV,OP>::convtol_default_ = 1e-8;
00329
00330 template<class ScalarType, class MV, class OP>
00331 const typename Teuchos::ScalarTraits<ScalarType>::magnitudeType GCRODRSolMgr<ScalarType,MV,OP>::orthoKappa_default_ = -1.0;
00332
00333 template<class ScalarType, class MV, class OP>
00334 const int GCRODRSolMgr<ScalarType,MV,OP>::maxRestarts_default_ = 20;
00335
00336 template<class ScalarType, class MV, class OP>
00337 const int GCRODRSolMgr<ScalarType,MV,OP>::maxIters_default_ = 1000;
00338
00339 template<class ScalarType, class MV, class OP>
00340 const int GCRODRSolMgr<ScalarType,MV,OP>::numBlocks_default_ = 300;
00341
00342 template<class ScalarType, class MV, class OP>
00343 const int GCRODRSolMgr<ScalarType,MV,OP>::recycledBlocks_default_ = 5;
00344
00345 template<class ScalarType, class MV, class OP>
00346 const int GCRODRSolMgr<ScalarType,MV,OP>::verbosity_default_ = Belos::Errors;
00347
00348 template<class ScalarType, class MV, class OP>
00349 const int GCRODRSolMgr<ScalarType,MV,OP>::outputFreq_default_ = -1;
00350
00351 template<class ScalarType, class MV, class OP>
00352 const std::string GCRODRSolMgr<ScalarType,MV,OP>::impResScale_default_ = "Norm of Preconditioned Initial Residual";
00353
00354 template<class ScalarType, class MV, class OP>
00355 const std::string GCRODRSolMgr<ScalarType,MV,OP>::expResScale_default_ = "Norm of Initial Residual";
00356
00357 template<class ScalarType, class MV, class OP>
00358 const std::string GCRODRSolMgr<ScalarType,MV,OP>::label_default_ = "Belos";
00359
00360 template<class ScalarType, class MV, class OP>
00361 const std::string GCRODRSolMgr<ScalarType,MV,OP>::orthoType_default_ = "DGKS";
00362
00363 template<class ScalarType, class MV, class OP>
00364 const Teuchos::RCP<ostream> GCRODRSolMgr<ScalarType,MV,OP>::outputStream_default_ = Teuchos::rcp(&std::cout,false);
00365
00366
00367
00368 template<class ScalarType, class MV, class OP>
00369 GCRODRSolMgr<ScalarType,MV,OP>::GCRODRSolMgr() :
00370 outputStream_(outputStream_default_),
00371 convtol_(convtol_default_),
00372 orthoKappa_(orthoKappa_default_),
00373 maxRestarts_(maxRestarts_default_),
00374 maxIters_(maxIters_default_),
00375 numBlocks_(numBlocks_default_),
00376 recycledBlocks_(recycledBlocks_default_),
00377 verbosity_(verbosity_default_),
00378 outputFreq_(outputFreq_default_),
00379 orthoType_(orthoType_default_),
00380 impResScale_(impResScale_default_),
00381 expResScale_(expResScale_default_),
00382 label_(label_default_),
00383 isSet_(false)
00384 {}
00385
00386
00387
00388 template<class ScalarType, class MV, class OP>
00389 GCRODRSolMgr<ScalarType,MV,OP>::GCRODRSolMgr(
00390 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
00391 const Teuchos::RCP<Teuchos::ParameterList> &pl ) :
00392 problem_(problem),
00393 outputStream_(outputStream_default_),
00394 convtol_(convtol_default_),
00395 orthoKappa_(orthoKappa_default_),
00396 maxRestarts_(maxRestarts_default_),
00397 maxIters_(maxIters_default_),
00398 numBlocks_(numBlocks_default_),
00399 recycledBlocks_(recycledBlocks_default_),
00400 verbosity_(verbosity_default_),
00401 outputFreq_(outputFreq_default_),
00402 orthoType_(orthoType_default_),
00403 impResScale_(impResScale_default_),
00404 expResScale_(expResScale_default_),
00405 label_(label_default_),
00406 isSet_(false)
00407 {
00408 TEST_FOR_EXCEPTION(problem_ == Teuchos::null, std::invalid_argument, "Problem not given to solver manager.");
00409
00410 if (!is_null(pl)) {
00411
00412 setParameters( pl );
00413 }
00414 }
00415
00416
00417 template<class ScalarType, class MV, class OP>
00418 void GCRODRSolMgr<ScalarType,MV,OP>::setParameters( const Teuchos::RCP<Teuchos::ParameterList> ¶ms )
00419 {
00420
00421 if (params_ == Teuchos::null) {
00422 params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
00423 }
00424 else {
00425 params->validateParameters(*getValidParameters());
00426 }
00427
00428
00429 if (params->isParameter("Maximum Restarts")) {
00430 maxRestarts_ = params->get("Maximum Restarts",maxRestarts_default_);
00431
00432
00433 params_->set("Maximum Restarts", maxRestarts_);
00434 }
00435
00436
00437 if (params->isParameter("Maximum Iterations")) {
00438 maxIters_ = params->get("Maximum Iterations",maxIters_default_);
00439
00440
00441 params_->set("Maximum Iterations", maxIters_);
00442 if (maxIterTest_!=Teuchos::null)
00443 maxIterTest_->setMaxIters( maxIters_ );
00444 }
00445
00446
00447 if (params->isParameter("Num Blocks")) {
00448 numBlocks_ = params->get("Num Blocks",numBlocks_default_);
00449 TEST_FOR_EXCEPTION(numBlocks_ <= 0, std::invalid_argument,
00450 "Belos::GCRODRSolMgr: \"Num Blocks\" must be strictly positive.");
00451
00452
00453 params_->set("Num Blocks", numBlocks_);
00454 }
00455
00456
00457 if (params->isParameter("Num Recycled Blocks")) {
00458 recycledBlocks_ = params->get("Num Recycled Blocks",recycledBlocks_default_);
00459 TEST_FOR_EXCEPTION(recycledBlocks_ <= 0, std::invalid_argument,
00460 "Belos::GCRODRSolMgr: \"Num Recycled Blocks\" must be strictly positive.");
00461
00462 TEST_FOR_EXCEPTION(recycledBlocks_ >= numBlocks_, std::invalid_argument,
00463 "Belos::GCRODRSolMgr: \"Num Recycled Blocks\" must be less than \"Num Blocks\".");
00464
00465
00466 params_->set("Num Recycled Blocks", recycledBlocks_);
00467 }
00468
00469
00470 if (params->isParameter("Timer Label")) {
00471 string tempLabel = params->get("Timer Label", label_default_);
00472
00473
00474 if (tempLabel != label_) {
00475 label_ = tempLabel;
00476 params_->set("Timer Label", label_);
00477 string solveLabel = label_ + ": GCRODRSolMgr total solve time";
00478 timerSolve_ = Teuchos::TimeMonitor::getNewTimer(solveLabel);
00479 }
00480 }
00481
00482
00483 if (params->isParameter("Orthogonalization")) {
00484 string tempOrthoType = params->get("Orthogonalization",orthoType_default_);
00485 TEST_FOR_EXCEPTION( tempOrthoType != "DGKS" && tempOrthoType != "ICGS" && tempOrthoType != "IMGS",
00486 std::invalid_argument,
00487 "Belos::GCRODRSolMgr: \"Orthogonalization\" must be either \"DGKS\", \"ICGS\", or \"IMGS\".");
00488 if (tempOrthoType != orthoType_) {
00489 orthoType_ = tempOrthoType;
00490
00491 if (orthoType_=="DGKS") {
00492 if (orthoKappa_ <= 0) {
00493 ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00494 }
00495 else {
00496 ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00497 Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
00498 }
00499 }
00500 else if (orthoType_=="ICGS") {
00501 ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00502 }
00503 else if (orthoType_=="IMGS") {
00504 ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00505 }
00506 }
00507 }
00508
00509
00510 if (params->isParameter("Orthogonalization Constant")) {
00511 orthoKappa_ = params->get("Orthogonalization Constant",orthoKappa_default_);
00512
00513
00514 params_->set("Orthogonalization Constant",orthoKappa_);
00515 if (orthoType_=="DGKS") {
00516 if (orthoKappa_ > 0 && ortho_ != Teuchos::null) {
00517 Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
00518 }
00519 }
00520 }
00521
00522
00523 if (params->isParameter("Verbosity")) {
00524 if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
00525 verbosity_ = params->get("Verbosity", verbosity_default_);
00526 } else {
00527 verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
00528 }
00529
00530
00531 params_->set("Verbosity", verbosity_);
00532 if (printer_ != Teuchos::null)
00533 printer_->setVerbosity(verbosity_);
00534 }
00535
00536
00537 if (params->isParameter("Output Stream")) {
00538 outputStream_ = Teuchos::getParameter<Teuchos::RCP<ostream> >(*params,"Output Stream");
00539
00540
00541 params_->set("Output Stream", outputStream_);
00542 if (printer_ != Teuchos::null)
00543 printer_->setOStream( outputStream_ );
00544 }
00545
00546
00547 if (verbosity_ & Belos::StatusTestDetails) {
00548 if (params->isParameter("Output Frequency")) {
00549 outputFreq_ = params->get("Output Frequency", outputFreq_default_);
00550 }
00551
00552
00553 params_->set("Output Frequency", outputFreq_);
00554 if (outputTest_ != Teuchos::null)
00555 outputTest_->setOutputFrequency( outputFreq_ );
00556 }
00557
00558
00559 if (printer_ == Teuchos::null) {
00560 printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
00561 }
00562
00563
00564 typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
00565 typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestResNorm_t;
00566
00567
00568 if (params->isParameter("Convergence Tolerance")) {
00569 convtol_ = params->get("Convergence Tolerance",convtol_default_);
00570
00571
00572 params_->set("Convergence Tolerance", convtol_);
00573 if (impConvTest_ != Teuchos::null)
00574 impConvTest_->setTolerance( convtol_ );
00575 if (expConvTest_ != Teuchos::null)
00576 expConvTest_->setTolerance( convtol_ );
00577 }
00578
00579
00580 if (params->isParameter("Implicit Residual Scaling")) {
00581 string tempImpResScale = Teuchos::getParameter<string>( *params, "Implicit Residual Scaling" );
00582
00583
00584 if (impResScale_ != tempImpResScale) {
00585 Belos::ScaleType impResScaleType = convertStringToScaleType( tempImpResScale );
00586 impResScale_ = tempImpResScale;
00587
00588
00589 params_->set("Implicit Residual Scaling", impResScale_);
00590 if (impConvTest_ != Teuchos::null) {
00591 try {
00592 impConvTest_->defineScaleForm( impResScaleType, Belos::TwoNorm );
00593 }
00594 catch (std::exception& e) {
00595
00596 impConvTest_ = Teuchos::null;
00597 convTest_ = Teuchos::null;
00598 }
00599 }
00600 }
00601 }
00602
00603 if (params->isParameter("Explicit Residual Scaling")) {
00604 string tempExpResScale = Teuchos::getParameter<string>( *params, "Explicit Residual Scaling" );
00605
00606
00607 if (expResScale_ != tempExpResScale) {
00608 Belos::ScaleType expResScaleType = convertStringToScaleType( tempExpResScale );
00609 expResScale_ = tempExpResScale;
00610
00611
00612 params_->set("Explicit Residual Scaling", expResScale_);
00613 if (expConvTest_ != Teuchos::null) {
00614 try {
00615 expConvTest_->defineScaleForm( expResScaleType, Belos::TwoNorm );
00616 }
00617 catch (std::exception& e) {
00618
00619 expConvTest_ = Teuchos::null;
00620 convTest_ = Teuchos::null;
00621 }
00622 }
00623 }
00624 }
00625
00626
00627
00628
00629 if (maxIterTest_ == Teuchos::null)
00630 maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
00631
00632
00633 if (impConvTest_ == Teuchos::null) {
00634 impConvTest_ = Teuchos::rcp( new StatusTestResNorm_t( convtol_ ) );
00635 impConvTest_->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
00636 }
00637
00638
00639 if (expConvTest_ == Teuchos::null) {
00640 expConvTest_ = Teuchos::rcp( new StatusTestResNorm_t( convtol_ ) );
00641 expConvTest_->defineResForm( StatusTestResNorm_t::Explicit, Belos::TwoNorm );
00642 expConvTest_->defineScaleForm( convertStringToScaleType(expResScale_), Belos::TwoNorm );
00643 }
00644
00645 if (convTest_ == Teuchos::null) {
00646 convTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
00647 }
00648
00649 sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
00650
00651 if (outputFreq_ > 0) {
00652 outputTest_ = Teuchos::rcp( new StatusTestOutput<ScalarType,MV,OP>( printer_,
00653 sTest_,
00654 outputFreq_,
00655 Passed+Failed+Undefined ) );
00656 }
00657 else {
00658 outputTest_ = Teuchos::rcp( new StatusTestOutput<ScalarType,MV,OP>( printer_,
00659 sTest_, 1 ) );
00660 }
00661
00662
00663 if (ortho_ == Teuchos::null) {
00664 if (orthoType_=="DGKS") {
00665 if (orthoKappa_ <= 0) {
00666 ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00667 }
00668 else {
00669 ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
00670 Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
00671 }
00672 }
00673 else if (orthoType_=="ICGS") {
00674 ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00675 }
00676 else if (orthoType_=="IMGS") {
00677 ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
00678 }
00679 else {
00680 TEST_FOR_EXCEPTION(orthoType_!="ICGS"&&orthoType_!="DGKS"&&orthoType_!="IMGS",std::logic_error,
00681 "Belos::GCRODRSolMgr(): Invalid orthogonalization type.");
00682 }
00683 }
00684
00685
00686 if (timerSolve_ == Teuchos::null) {
00687 string solveLabel = label_ + ": GCRODRSolMgr total solve time";
00688 timerSolve_ = Teuchos::TimeMonitor::getNewTimer(solveLabel);
00689 }
00690
00691
00692 isSet_ = true;
00693 }
00694
00695
00696 template<class ScalarType, class MV, class OP>
00697 Teuchos::RCP<const Teuchos::ParameterList>
00698 GCRODRSolMgr<ScalarType,MV,OP>::getValidParameters() const
00699 {
00700 static Teuchos::RCP<const Teuchos::ParameterList> validPL;
00701 if (is_null(validPL)) {
00702 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
00703
00704 pl->set("Convergence Tolerance", convtol_default_,
00705 "The relative residual tolerance that needs to be achieved by the\n"
00706 "iterative solver in order for the linear system to be declared converged.");
00707 pl->set("Maximum Restarts", maxRestarts_default_,
00708 "The maximum number of cycles allowed for each\n"
00709 "set of RHS solved.");
00710 pl->set("Maximum Iterations", maxIters_default_,
00711 "The maximum number of iterations allowed for each\n"
00712 "set of RHS solved.");
00713 pl->set("Num Blocks", numBlocks_default_,
00714 "The maximum number of vectors allowed in the Krylov subspace\n"
00715 "for each set of RHS solved.");
00716 pl->set("Num Recycled Blocks", recycledBlocks_default_,
00717 "The maximum number of vectors in the recycled subspace." );
00718 pl->set("Verbosity", verbosity_default_,
00719 "What type(s) of solver information should be outputted\n"
00720 "to the output stream.");
00721 pl->set("Output Frequency", outputFreq_default_,
00722 "How often convergence information should be outputted\n"
00723 "to the output stream.");
00724 pl->set("Output Stream", outputStream_default_,
00725 "A reference-counted pointer to the output stream where all\n"
00726 "solver output is sent.");
00727 pl->set("Implicit Residual Scaling", impResScale_default_,
00728 "The type of scaling used in the implicit residual convergence test.");
00729 pl->set("Explicit Residual Scaling", expResScale_default_,
00730 "The type of scaling used in the explicit residual convergence test.");
00731 pl->set("Timer Label", label_default_,
00732 "The string to use as a prefix for the timer labels.");
00733
00734 pl->set("Orthogonalization", orthoType_default_,
00735 "The type of orthogonalization to use: DGKS, ICGS, IMGS");
00736 pl->set("Orthogonalization Constant",orthoKappa_default_,
00737 "The constant used by DGKS orthogonalization to determine\n"
00738 "whether another step of classical Gram-Schmidt is necessary.");
00739 validPL = pl;
00740 }
00741 return validPL;
00742 }
00743
00744
00745
00746 template<class ScalarType, class MV, class OP>
00747 ReturnType GCRODRSolMgr<ScalarType,MV,OP>::solve() {
00748
00749
00750
00751
00752 if (!isSet_) { setParameters( params_ ); }
00753
00754 Teuchos::BLAS<int,ScalarType> blas;
00755 Teuchos::LAPACK<int,ScalarType> lapack;
00756 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
00757 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
00758
00759 TEST_FOR_EXCEPTION(problem_ == Teuchos::null,GCRODRSolMgrLinearProblemFailure,
00760 "Belos::GCRODRSolMgr::solve(): Linear problem is not a valid object.");
00761
00762 TEST_FOR_EXCEPTION(!problem_->isProblemSet(),GCRODRSolMgrLinearProblemFailure,
00763 "Belos::GCRODRSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
00764
00765
00766 int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
00767 std::vector<int> currIdx(1);
00768 currIdx[0] = 0;
00769
00770
00771 problem_->setLSIndex( currIdx );
00772
00774 int keff = 0;
00775 if (U_ != Teuchos::null) {
00776 keff = MVT::GetNumberVecs( *U_ );
00777 }
00778
00779
00780 int dim = MVT::GetVecLength( *(problem_->getRHS()) );
00781 if (numBlocks_ > dim) {
00782 numBlocks_ = dim;
00783 printer_->stream(Warnings) <<
00784 "Warning! Requested Krylov subspace dimension is larger that operator dimension!" << endl <<
00785 " The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << numBlocks_ << endl;
00786 params_->set("Num Blocks", numBlocks_);
00787 }
00788
00789
00790 bool isConverged = true;
00791
00793
00794 Teuchos::ParameterList plist;
00795
00796 plist.set("Num Blocks",numBlocks_);
00797 plist.set("Recycled Blocks",recycledBlocks_);
00798
00800
00801
00802 Teuchos::RCP<GCRODRIter<ScalarType,MV,OP> > gcrodr_iter;
00803 gcrodr_iter = Teuchos::rcp( new GCRODRIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
00804
00805 int prime_iterations = 0;
00806
00807
00808 {
00809 Teuchos::TimeMonitor slvtimer(*timerSolve_);
00810
00811 while ( numRHS2Solve > 0 ) {
00812
00813
00814 outputTest_->reset();
00815
00817
00818
00819
00820 if (keff > 0) {
00821 TEST_FOR_EXCEPTION(keff < recycledBlocks_,GCRODRSolMgrRecyclingFailure,
00822 "Belos::GCRODRSolMgr::solve(): Requested size of recycled subspace is not consistent with the current recycle subspace.");
00823
00824 printer_->stream(Debug) << " Now solving RHS index " << currIdx[0] << " using recycled subspace of dimension " << keff << std::endl << std::endl;
00825
00826 std::vector<int> index(keff);
00827 for (int i=0; i<keff; ++i) { index[i] = i; }
00828 Teuchos::RCP<MV> Ckeff = MVT::CloneView( *C_, index );
00829 Teuchos::RCP<MV> Ukeff = MVT::CloneView( *U_, index );
00830 problem_->apply( *Ukeff, *Ckeff );
00831
00832
00833
00834 Teuchos::SerialDenseMatrix<int,ScalarType> R(keff,keff);
00835 int rank = ortho_->normalize(*Ckeff, Teuchos::rcp(&R,false));
00836
00837
00838 TEST_FOR_EXCEPTION(rank != keff,GCRODRSolMgrOrthoFailure,
00839 "Belos::GCRODRSolMgr::solve(): Failed to compute orthonormal basis for initial recycled subspace.");
00840
00841
00842
00843 int info = 0;
00844 std::vector<int> ipiv(R.numRows());
00845 lapack.GETRF(R.numRows(),R.numCols(),R.values(),R.numRows(),&ipiv[0],&info);
00846 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
00847 "Belos::GCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
00848
00849
00850 int lwork = R.numRows();
00851 std::vector<ScalarType> work(lwork);
00852 lapack.GETRI(R.numRows(),R.values(),R.numRows(),&ipiv[0],&work[0],lwork,&info);
00853 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
00854 "Belos::GCRODRSolMgr::solve(): LAPACK _GETRI failed to compute an LU factorization.");
00855
00856 Teuchos::RCP<MV> tempU = MVT::Clone( *U_, keff );
00857 MVT::MvTimesMatAddMv( one, *U_, R, zero, *tempU );
00858 U_ = tempU;
00859
00860
00861 Teuchos::SerialDenseMatrix<int,ScalarType> Ctr(keff,1);
00862 problem_->computeCurrPrecResVec( &*r_ );
00863 MVT::MvTransMv( one, *C_, *r_, Ctr );
00864
00865
00866 MVT::MvTimesMatAddMv( one, *U_, Ctr, one, *problem_->getCurrLHSVec() );
00867
00868
00869 MVT::MvTimesMatAddMv( -one, *C_, Ctr, one, *r_ );
00870
00871
00872 prime_iterations = 0;
00873
00874 }
00875 else {
00876
00877
00878 printer_->stream(Debug) << " No recycled subspace available for RHS index " << currIdx[0] << std::endl << std::endl;
00879
00880 Teuchos::ParameterList primeList;
00881
00882
00883 primeList.set("Block Size",1);
00884 primeList.set("Num Blocks",numBlocks_);
00885 primeList.set("Keep Hessenberg", true);
00886 primeList.set("Initialize Hessenberg", true);
00887
00888
00889 Teuchos::RCP<BlockGmresIter<ScalarType,MV,OP> > gmres_iter;
00890 gmres_iter = Teuchos::rcp( new BlockGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,primeList) );
00891
00892
00893 if (r_ == Teuchos::null)
00894 r_ = MVT::Clone( *(problem_->getRHS()), 1 );
00895 problem_->computeCurrPrecResVec( &*r_ );
00896 Teuchos::RCP<MV> V_0 = MVT::CloneCopy( *r_ );
00897
00898
00899 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > z_0 =
00900 rcp( new Teuchos::SerialDenseVector<int,ScalarType>(1) );
00901
00902
00903 int rank = ortho_->normalize( *V_0, z_0 );
00904 TEST_FOR_EXCEPTION(rank != 1,GCRODRSolMgrOrthoFailure,
00905 "Belos::GCRODRSolMgr::solve(): Failed to compute initial block of orthonormal vectors for priming solve.");
00906
00907
00908 GmresIterationState<ScalarType,MV> newstate;
00909 newstate.V = V_0;
00910 newstate.z = z_0;
00911 newstate.curDim = 0;
00912 gmres_iter->initialize(newstate);
00913
00914
00915 bool primeConverged = false;
00916 try {
00917 gmres_iter->iterate();
00918
00919
00920 if ( convTest_->getStatus() == Passed ) {
00921
00922 primeConverged = true;
00923 }
00924 }
00925 catch (GmresIterationOrthoFailure e) {
00926
00927 gmres_iter->updateLSQR( gmres_iter->getCurSubspaceDim() );
00928
00929
00930 sTest_->checkStatus( &*gmres_iter );
00931 if (convTest_->getStatus() == Passed)
00932 primeConverged = true;
00933 }
00934 catch (std::exception e) {
00935 printer_->stream(Errors) << "Error! Caught exception in GCRODRIter::iterate() at iteration "
00936 << gmres_iter->getNumIters() << endl
00937 << e.what() << endl;
00938 throw;
00939 }
00940
00941
00942 prime_iterations = gmres_iter->getNumIters();
00943
00944
00945 Teuchos::RCP<MV> update = gmres_iter->getCurrentUpdate();
00946 problem_->updateSolution( update, true );
00947
00948
00949 GmresIterationState<ScalarType,MV> oldState = gmres_iter->getState();
00950 int p = oldState.curDim;
00951
00952
00953
00954
00955
00956
00957 if (recycledBlocks_ < p+1) {
00958 int info = 0;
00959 Teuchos::SerialDenseMatrix<int,ScalarType> PP( p, recycledBlocks_+1 );
00960 keff = getHarmonicVecs1( p, *oldState.H, PP );
00961
00962
00963 if (U_ == Teuchos::null) {
00964 U_ = MVT::Clone( *problem_->getRHS(), keff );
00965 }
00966 else {
00967 if (MVT::GetNumberVecs( *U_ ) < keff) {
00968 U_ = MVT::Clone( *problem_->getRHS(), keff );
00969 }
00970 }
00971 if (C_ == Teuchos::null) {
00972 C_ = MVT::Clone( *problem_->getRHS(), keff );
00973 }
00974 else {
00975 if (MVT::GetNumberVecs( *C_ ) < keff) {
00976 C_ = MVT::Clone( *problem_->getRHS(), keff );
00977 }
00978 }
00979
00980
00981
00982 std::vector<int> index( p );
00983 for (int i=0; i < p; ++i) { index[i] = i; }
00984 Teuchos::RCP<const MV> Vp = MVT::CloneView( *oldState.V, index );
00985 index.resize(keff);
00986 Teuchos::RCP<MV> Up = MVT::CloneView( *U_, index );
00987 const Teuchos::SerialDenseMatrix<int,ScalarType> PPview( Teuchos::View, PP, p, keff );
00988 MVT::MvTimesMatAddMv( one, *Vp, PPview, zero, *Up );
00989 Vp = null;
00990
00991
00992
00993 const Teuchos::SerialDenseMatrix<int,ScalarType> Hview( Teuchos::View, *oldState.H, p+1, p );
00994 Teuchos::SerialDenseMatrix<int,ScalarType> HP( p+1, keff );
00995 HP.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, one, Hview, PPview, zero );
00996
00997
00998 int lwork = -1;
00999 std::vector<ScalarType> tau(keff);
01000 lapack.GEQRF(HP.numRows(),HP.numCols(),HP.values(),HP.numRows(),&tau[0],
01001 &tau[0],lwork,&info);
01002 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01003 "Belos::GCRODRSolMgr::solve(): LAPACK _GEQRF failed to compute a workspace size.");
01004
01005 lwork = (int)tau[0];
01006 std::vector<ScalarType> work(lwork);
01007 lapack.GEQRF(HP.numRows(),HP.numCols(),HP.values(),HP.numRows(),&tau[0],
01008 &work[0],lwork,&info);
01009 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01010 "Belos::GCRODRSolMgr::solve(): LAPACK _GEQRF failed to compute a QR factorization.");
01011
01012
01013
01014 Teuchos::SerialDenseMatrix<int,ScalarType> R( Teuchos::Copy, HP, keff, keff );
01015 for(int i=1;i<keff;i++) { for(int j=0;j<i;j++) R(i,j) = zero; }
01016 lapack.ORGQR(HP.numRows(),HP.numCols(),HP.numCols(),HP.values(),HP.numRows(),&tau[0],
01017 &work[0],lwork,&info);
01018 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01019 "Belos::GCRODRSolMgr::solve(): LAPACK _ORGQR failed to construct the Q factor.");
01020
01021
01022 index.resize( p+1 );
01023 for (int i=0; i < p+1; ++i) { index[i] = i; }
01024 Vp = MVT::CloneView( *oldState.V, index );
01025 index.resize(keff);
01026 for (int i=0; i < keff; ++i) { index[i] = i; }
01027 Teuchos::RCP<MV> Cp = MVT::CloneView( *C_, index );
01028 MVT::MvTimesMatAddMv( one, *Vp, HP, zero, *Cp );
01029 Vp = null;
01030 Cp = null;
01031
01032
01033
01034 std::vector<int> ipiv(R.numRows());
01035 lapack.GETRF(R.numRows(),R.numCols(),R.values(),R.numRows(),&ipiv[0],&info);
01036 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01037 "Belos::GCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
01038
01039
01040 lwork = R.numRows();
01041 work.resize(lwork);
01042 lapack.GETRI(R.numRows(),R.values(),R.numRows(),&ipiv[0],&work[0],lwork,&info);
01043 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01044 "Belos::GCRODRSolMgr::solve(): LAPACK _GETRI failed to compute an LU factorization.");
01045
01046 Teuchos::RCP<MV> tempU = MVT::Clone( *U_, keff );
01047 MVT::MvTimesMatAddMv( one, *Up, R, zero, *tempU );
01048 U_ = tempU;
01049
01050 printer_->stream(Debug) << " Generated recycled subspace using RHS index " << currIdx[0] << " of dimension " << keff << std::endl << std::endl;
01051 }
01052
01053
01054 if (primeConverged) {
01055
01056 problem_->setCurrLS();
01057
01058
01059 numRHS2Solve -= 1;
01060 if ( numRHS2Solve > 0 ) {
01061 currIdx[0]++;
01062
01063
01064 problem_->setLSIndex( currIdx );
01065 }
01066 else {
01067 currIdx.resize( numRHS2Solve );
01068 }
01069
01070 continue;
01071 }
01072 }
01073
01074
01075
01076
01077 gcrodr_iter->setSize( keff, numBlocks_ );
01078
01079
01080 gcrodr_iter->resetNumIters(prime_iterations);
01081
01082
01083 outputTest_->resetNumCalls();
01084
01085
01086 if (r_ == Teuchos::null)
01087 r_ = MVT::Clone( *(problem_->getRHS()), 1 );
01088 problem_->computeCurrPrecResVec( &*r_ );
01089
01090
01091 Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > z_0 =
01092 rcp( new Teuchos::SerialDenseVector<int,ScalarType>(1) );
01093
01094
01095 int rank = ortho_->normalize( *r_, z_0 );
01096 TEST_FOR_EXCEPTION(rank != 1,GCRODRSolMgrOrthoFailure,
01097 "Belos::GCRODRSolMgr::solve(): Failed to compute initial block of orthonormal vectors.");
01098
01099
01100 GCRODRIterState<ScalarType,MV> newstate;
01101 newstate.V = r_;
01102 newstate.z = z_0;
01103 newstate.U = U_;
01104 newstate.C = C_;
01105 newstate.curDim = 0;
01106 gcrodr_iter->initialize(newstate);
01107 int numRestarts = 0;
01108 while(1) {
01109
01110
01111 try {
01112 gcrodr_iter->iterate();
01113
01115
01116
01117
01119 if ( convTest_->getStatus() == Passed ) {
01120
01121 break;
01122 }
01124
01125
01126
01128 else if ( maxIterTest_->getStatus() == Passed ) {
01129
01130 isConverged = false;
01131 break;
01132 }
01134
01135
01136
01138 else if ( gcrodr_iter->getCurSubspaceDim() == gcrodr_iter->getMaxSubspaceDim() ) {
01139
01140
01141
01142
01143 GCRODRIterState<ScalarType,MV> oldState = gcrodr_iter->getState();
01144 int p = oldState.curDim;
01145
01146
01147 Teuchos::RCP<MV> update = gcrodr_iter->getCurrentUpdate();
01148 problem_->updateSolution( update, true );
01149
01150
01151 std::vector<MagnitudeType> d(keff);
01152 MVT::MvNorm( *U_, d );
01153 for (int i=0; i<keff; ++i) {
01154 d[i] = one / d[i];
01155 }
01156 MVT::MvScale( *U_, d );
01157
01158
01159 Teuchos::SerialDenseMatrix<int,ScalarType> H2(p+keff+1, p+keff);
01160
01161
01162 for (int i=0; i<keff; ++i) {
01163 H2(i,i) = d[i];
01164 }
01165
01166
01167 Teuchos::SerialDenseMatrix<int,ScalarType> H2temp(Teuchos::View, H2, keff, p, 0, keff);
01168 H2temp.assign(*oldState.B);
01169
01170
01171 Teuchos::SerialDenseMatrix<int,ScalarType> H2temp2(Teuchos::View, H2, p+1, p, keff, keff);
01172 H2temp2.assign(*oldState.H);
01173
01174
01175 Teuchos::SerialDenseMatrix<int,ScalarType> PP( p+keff, recycledBlocks_+1 );
01176 int keff_new = getHarmonicVecs2( keff, p, H2, oldState.V, PP );
01177 printer_->stream(Debug) << " Generated new recycled subspace using RHS index " << currIdx[0] << " of dimension " << keff_new << std::endl << std::endl;
01178
01179
01180
01181
01182
01183 Teuchos::RCP<MV> tempU = MVT::Clone( *U_, keff_new );
01184 Teuchos::SerialDenseMatrix<int,ScalarType> tempPP( Teuchos::View, PP, keff, keff_new );
01185 MVT::MvTimesMatAddMv( one, *U_, tempPP, zero, *tempU );
01186
01187
01188 std::vector<int> index(p);
01189 for (int i=0; i < p; i++) { index[i] = i; }
01190 Teuchos::RCP<const MV> Vp = MVT::CloneView( *oldState.V, index );
01191 Teuchos::SerialDenseMatrix<int,ScalarType> tempPP2( Teuchos::View, PP, p, keff_new, keff );
01192 MVT::MvTimesMatAddMv( one, *Vp, tempPP2, one, *tempU );
01193
01194
01195 Teuchos::SerialDenseMatrix<int,ScalarType> tempPP3( Teuchos::View, PP, p+keff, keff_new );
01196 Teuchos::SerialDenseMatrix<int,ScalarType> HP(H2.numRows(),tempPP3.numCols());
01197 HP.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,H2,tempPP3,zero);
01198
01199
01200 int info = 0, lwork = -1;
01201 std::vector<ScalarType> tau(keff_new);
01202 lapack.GEQRF(HP.numRows(),HP.numCols(),HP.values(),HP.numRows(),&tau[0],
01203 &tau[0],lwork,&info);
01204 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01205 "Belos::GCRODRSolMgr::solve(): LAPACK _GEQRF failed to compute a workspace size.");
01206
01207 lwork = (int)tau[0];
01208 std::vector<ScalarType> work(lwork);
01209 lapack.GEQRF(HP.numRows(),HP.numCols(),HP.values(),HP.numRows(),&tau[0],
01210 &work[0],lwork,&info);
01211 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01212 "Belos::GCRODRSolMgr::solve(): LAPACK _GEQRF failed to compute a QR factorization.");
01213
01214
01215
01216 Teuchos::SerialDenseMatrix<int,ScalarType> R( Teuchos::Copy, HP, keff_new, keff_new );
01217 for(int i=1;i<keff_new;i++) { for(int j=0;j<i;j++) R(i,j) = zero; }
01218 lapack.ORGQR(HP.numRows(),HP.numCols(),HP.numCols(),HP.values(),HP.numRows(),&tau[0],
01219 &work[0],lwork,&info);
01220 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01221 "Belos::GCRODRSolMgr::solve(): LAPACK _ORGQR failed to construct the Q factor.");
01222
01223
01224
01225
01226
01227 index.resize(keff_new);
01228 for (int i=0; i < keff_new; i++) { index[i] = i; }
01229 Teuchos::RCP<MV> tempC = MVT::Clone( *C_, keff_new );
01230 Teuchos::SerialDenseMatrix<int,ScalarType> tempQ( Teuchos::View, HP, keff, keff_new );
01231 MVT::MvTimesMatAddMv( one, *C_, tempQ, zero, *tempC );
01232
01233
01234 index.resize( p+1 );
01235 for (int i=0; i < p+1; ++i) { index[i] = i; }
01236 Vp = MVT::CloneView( *oldState.V, index );
01237 Teuchos::SerialDenseMatrix<int,ScalarType> tempQ2( Teuchos::View, HP, p+1, keff_new, keff );
01238 MVT::MvTimesMatAddMv( one, *Vp, tempQ2, one, *tempC );
01239 C_ = tempC;
01240
01241
01242
01243 std::vector<int> ipiv(R.numRows());
01244 lapack.GETRF(R.numRows(),R.numCols(),R.values(),R.numRows(),&ipiv[0],&info);
01245 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01246 "Belos::GCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
01247
01248
01249 lwork = R.numRows();
01250 work.resize(lwork);
01251 lapack.GETRI(R.numRows(),R.values(),R.numRows(),&ipiv[0],&work[0],lwork,&info);
01252 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01253 "Belos::GCRODRSolMgr::solve(): LAPACK _GETRI failed to compute an LU factorization.");
01254
01255 if (keff != keff_new) {
01256 U_ = MVT::Clone( *problem_->getRHS(), keff_new );
01257 }
01258 MVT::MvTimesMatAddMv( one, *tempU, R, zero, *U_ );
01259
01260
01261 if ( numRestarts >= maxRestarts_ ) {
01262 isConverged = false;
01263 break;
01264 }
01265 numRestarts++;
01266
01267 printer_->stream(Debug) << " Performing restart number " << numRestarts << " of " << maxRestarts_ << std::endl << std::endl;
01268
01269
01270 if (r_ == Teuchos::null)
01271 r_ = MVT::Clone( *(problem_->getRHS()), 1 );
01272 problem_->computeCurrPrecResVec( &*r_ );
01273
01274
01275 Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > z_0 =
01276 rcp( new Teuchos::SerialDenseVector<int,ScalarType>(1) );
01277
01278
01279 int rank = ortho_->normalize( *r_, z_0 );
01280 TEST_FOR_EXCEPTION(rank != 1,GCRODRSolMgrOrthoFailure,
01281 "Belos::GCRODRSolMgr::solve(): Failed to compute initial block of orthonormal vectors after restart.");
01282
01283
01284 keff = keff_new;
01285 gcrodr_iter->setSize( keff, numBlocks_ );
01286
01287
01288 GCRODRIterState<ScalarType,MV> newstate;
01289 newstate.V = r_;
01290 newstate.z = z_0;
01291 newstate.U = U_;
01292 newstate.C = C_;
01293 newstate.curDim = 0;
01294 gcrodr_iter->initialize(newstate);
01295
01296 }
01297
01299
01300
01301
01302
01304
01305 else {
01306 TEST_FOR_EXCEPTION(true,std::logic_error,
01307 "Belos::GCRODRSolMgr::solve(): Invalid return from GCRODRIter::iterate().");
01308 }
01309 }
01310 catch (GCRODRIterOrthoFailure e) {
01311
01312 gcrodr_iter->updateLSQR( gcrodr_iter->getCurSubspaceDim() );
01313
01314
01315 sTest_->checkStatus( &*gcrodr_iter );
01316 if (convTest_->getStatus() != Passed)
01317 isConverged = false;
01318 break;
01319 }
01320 catch (std::exception e) {
01321 printer_->stream(Errors) << "Error! Caught exception in GCRODRIter::iterate() at iteration "
01322 << gcrodr_iter->getNumIters() << endl
01323 << e.what() << endl;
01324 throw;
01325 }
01326 }
01327
01328
01329
01330 Teuchos::RCP<MV> update = gcrodr_iter->getCurrentUpdate();
01331 problem_->updateSolution( update, true );
01332
01333
01334 problem_->setCurrLS();
01335
01336
01337 numRHS2Solve -= 1;
01338 if ( numRHS2Solve > 0 ) {
01339 currIdx[0]++;
01340
01341
01342 problem_->setLSIndex( currIdx );
01343 }
01344 else {
01345 currIdx.resize( numRHS2Solve );
01346 }
01347
01348 }
01349
01350 }
01351
01352
01353 sTest_->print( printer_->stream(FinalSummary) );
01354
01355
01356 Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
01357
01358 if (!isConverged) {
01359 return Unconverged;
01360 }
01361 return Converged;
01362 }
01363
01364
01365
01366 template<class ScalarType, class MV, class OP>
01367 int GCRODRSolMgr<ScalarType,MV,OP>::getHarmonicVecs1(int m,
01368 const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
01369 Teuchos::SerialDenseMatrix<int,ScalarType>& PP)
01370 {
01371 int i, j;
01372 bool xtraVec = false;
01373 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
01374 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
01375
01376
01377 Teuchos::LAPACK<int,ScalarType> lapack;
01378
01379
01380 std::vector<MagnitudeType> wr(m), wi(m);
01381
01382
01383 Teuchos::SerialDenseMatrix<int,ScalarType> vr(m,m);
01384
01385
01386 std::vector<MagnitudeType> w(m);
01387
01388
01389 std::vector<int> iperm(m);
01390
01391
01392 int lwork = 4*m;
01393 std::vector<ScalarType> work(lwork);
01394
01395
01396 int info = 0;
01397
01398
01399 Teuchos::SerialDenseMatrix<int, ScalarType> HHt( HH, Teuchos::TRANS );
01400 Teuchos::SerialDenseVector<int, ScalarType> e_m( m );
01401 e_m[m-1] = one;
01402 lapack.GESV(m, 1, HHt.values(), HHt.stride(), &iperm[0], e_m.values(), e_m.stride(), &info);
01403 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01404 "Belos::GCRODRSolMgr::solve(): LAPACK GESV failed to compute a solution.");
01405
01406
01407 ScalarType d = HH(m, m-1) * HH(m, m-1);
01408 Teuchos::SerialDenseMatrix<int, ScalarType> harmHH( HH );
01409 for( i=0; i<m; ++i ) {
01410 harmHH(i, m-1) += d * e_m[i];
01411 }
01412
01413
01414
01415 const int ldvl = m;
01416 ScalarType* vl = 0;
01417 lapack.GEEV('N', 'V', m, harmHH.values(), harmHH.stride(), &wr[0], &wi[0],
01418 vl, ldvl, vr.values(), vr.stride(), &work[0], lwork, &info);
01419 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01420 "Belos::GCRODRSolMgr::solve(): LAPACK GEEV failed to compute eigensolutions.");
01421
01422
01423 for( i=0; i<m; ++i ) {
01424 w[i] = Teuchos::ScalarTraits<ScalarType>::squareroot( wr[i]*wr[i] + wi[i]*wi[i] );
01425 }
01426
01427
01428 this->sort(w, m, iperm);
01429
01430
01431 if (wi[iperm[recycledBlocks_-1]] != zero) {
01432 int countImag = 0;
01433 for ( i=0; i<recycledBlocks_; ++i ) {
01434 if (wi[iperm[i]] != zero)
01435 countImag++;
01436 }
01437
01438 if (countImag % 2)
01439 xtraVec = true;
01440 }
01441
01442
01443 for( i=0; i<recycledBlocks_; ++i ) {
01444 for( j=0; j<m; j++ ) {
01445 PP(j,i) = vr(j,iperm[i]);
01446 }
01447 }
01448
01449 if (xtraVec) {
01450 if (wi[iperm[recycledBlocks_-1]] > 0) {
01451 for( j=0; j<m; ++j ) {
01452 PP(j,recycledBlocks_) = vr(j,iperm[recycledBlocks_-1]+1);
01453 }
01454 }
01455 else {
01456 for( j=0; j<m; ++j ) {
01457 PP(j,recycledBlocks_) = vr(j,iperm[recycledBlocks_-1]-1);
01458 }
01459 }
01460 }
01461
01462
01463 if (xtraVec) {
01464 return recycledBlocks_+1;
01465 }
01466 return recycledBlocks_;
01467 }
01468
01469
01470 template<class ScalarType, class MV, class OP>
01471 int GCRODRSolMgr<ScalarType,MV,OP>::getHarmonicVecs2(int keff, int m,
01472 const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
01473 const Teuchos::RCP<const MV>& VV,
01474 Teuchos::SerialDenseMatrix<int,ScalarType>& PP)
01475 {
01476 int i, j;
01477 int m2 = HH.numCols();
01478 bool xtraVec = false;
01479 ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
01480 ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
01481 std::vector<int> index;
01482
01483
01484 Teuchos::LAPACK<int,ScalarType> lapack;
01485
01486
01487 std::vector<ScalarType> wr(m2), wi(m2);
01488
01489
01490 std::vector<MagnitudeType> w(m2);
01491
01492
01493 Teuchos::SerialDenseMatrix<int,ScalarType> vr(m2,m2);
01494
01495
01496 std::vector<int> iperm(m2);
01497
01498
01499
01500
01501 Teuchos::SerialDenseMatrix<int,ScalarType> B(m2,m2);
01502 B.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,HH,HH,zero);
01503
01504
01505
01506 Teuchos::SerialDenseMatrix<int,ScalarType> A_tmp( keff+m+1, keff+m );
01507
01508
01509 Teuchos::SerialDenseMatrix<int,ScalarType> A11( Teuchos::View, A_tmp, keff, keff );
01510 MVT::MvTransMv( one, *C_, *U_, A11 );
01511
01512
01513 Teuchos::SerialDenseMatrix<int,ScalarType> A21( Teuchos::View, A_tmp, m+1, keff, keff );
01514 index.resize(m+1);
01515 for (i=0; i < m+1; i++) { index[i] = i; }
01516 Teuchos::RCP<const MV> Vp = MVT::CloneView( *VV, index );
01517 MVT::MvTransMv( one, *Vp, *U_, A21 );
01518
01519
01520 for( i=keff; i<keff+m; i++ ) {
01521 A_tmp(i,i) = one;
01522 }
01523
01524
01525 Teuchos::SerialDenseMatrix<int,ScalarType> A( m2, A_tmp.numCols() );
01526 A.multiply( Teuchos::TRANS, Teuchos::NO_TRANS, one, HH, A_tmp, zero );
01527
01528
01529
01530
01531
01532
01533
01534 char balanc='P', jobvl='N', jobvr='V', sense='N';
01535 int ld = A.numRows();
01536 int lwork = 6*ld;
01537 int ldvl = ld, ldvr = ld;
01538 int info = 0,ilo = 0,ihi = 0;
01539 ScalarType abnrm = zero, bbnrm = zero;
01540 ScalarType *vl = 0;
01541 std::vector<ScalarType> beta(ld);
01542 std::vector<ScalarType> work(lwork);
01543 std::vector<MagnitudeType> lscale(ld), rscale(ld);
01544 std::vector<MagnitudeType> rconde(ld), rcondv(ld);
01545 std::vector<int> iwork(ld+6);
01546 int *bwork = 0;
01547 lapack.GGEVX(balanc, jobvl, jobvr, sense, ld, A.values(), ld, B.values(), ld, &wr[0], &wi[0],
01548 &beta[0], vl, ldvl, vr.values(), ldvr, &ilo, &ihi, &lscale[0], &rscale[0],
01549 &abnrm, &bbnrm, &rconde[0], &rcondv[0], &work[0], lwork, &iwork[0], bwork, &info);
01550 TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,
01551 "Belos::GCRODRSolMgr::solve(): LAPACK GGEVX failed to compute eigensolutions.");
01552
01553
01554
01555 for( i=0; i<ld; i++ ) {
01556 w[i] = Teuchos::ScalarTraits<ScalarType>::squareroot( (wr[i]/beta[i])*(wr[i]/beta[i]) + (wi[i]/beta[i])*(wi[i]/beta[i]) );
01557 }
01558
01559
01560 this->sort(w,ld,iperm);
01561
01562
01563 if (wi[iperm[ld-recycledBlocks_]] != zero) {
01564 int countImag = 0;
01565 for ( i=ld-recycledBlocks_; i<ld; i++ ) {
01566 if (wi[iperm[i]] != zero)
01567 countImag++;
01568 }
01569
01570 if (countImag % 2)
01571 xtraVec = true;
01572 }
01573
01574
01575 for( i=0; i<recycledBlocks_; i++ ) {
01576 for( j=0; j<ld; j++ ) {
01577 PP(j,i) = vr(j,iperm[ld-recycledBlocks_+i]);
01578 }
01579 }
01580
01581 if (xtraVec) {
01582 if (wi[iperm[ld-recycledBlocks_]] > 0) {
01583 for( j=0; j<ld; j++ ) {
01584 PP(j,recycledBlocks_) = vr(j,iperm[ld-recycledBlocks_]+1);
01585 }
01586 }
01587 else {
01588 for( j=0; j<ld; j++ ) {
01589 PP(j,recycledBlocks_) = vr(j,iperm[ld-recycledBlocks_]-1);
01590 }
01591 }
01592 }
01593
01594
01595 if (xtraVec) {
01596 return recycledBlocks_+1;
01597 }
01598 return recycledBlocks_;
01599 }
01600
01601
01602
01603 template<class ScalarType, class MV, class OP>
01604 void GCRODRSolMgr<ScalarType,MV,OP>::sort(std::vector<ScalarType>& dlist, int n, std::vector<int>& iperm)
01605 {
01606 int l, r, j, i, flag;
01607 int RR2;
01608 double dRR, dK;
01609
01610
01611 for(j=0;j<n;j++)
01612 iperm[j] = j;
01613
01614 if (n <= 1) return;
01615
01616 l = n / 2 + 1;
01617 r = n - 1;
01618 l = l - 1;
01619 dRR = dlist[l - 1];
01620 dK = dlist[l - 1];
01621
01622 RR2 = iperm[l - 1];
01623 while (r != 0) {
01624 j = l;
01625 flag = 1;
01626
01627 while (flag == 1) {
01628 i = j;
01629 j = j + j;
01630
01631 if (j > r + 1)
01632 flag = 0;
01633 else {
01634 if (j < r + 1)
01635 if (dlist[j] > dlist[j - 1]) j = j + 1;
01636
01637 if (dlist[j - 1] > dK) {
01638 dlist[i - 1] = dlist[j - 1];
01639 iperm[i - 1] = iperm[j - 1];
01640 }
01641 else {
01642 flag = 0;
01643 }
01644 }
01645 }
01646 dlist[i - 1] = dRR;
01647 iperm[i - 1] = RR2;
01648
01649 if (l == 1) {
01650 dRR = dlist [r];
01651 RR2 = iperm[r];
01652 dK = dlist[r];
01653 dlist[r] = dlist[0];
01654 iperm[r] = iperm[0];
01655 r = r - 1;
01656 }
01657 else {
01658 l = l - 1;
01659 dRR = dlist[l - 1];
01660 RR2 = iperm[l - 1];
01661 dK = dlist[l - 1];
01662 }
01663 }
01664 dlist[0] = dRR;
01665 iperm[0] = RR2;
01666 }
01667
01668
01669
01670 template<class ScalarType, class MV, class OP>
01671 std::string GCRODRSolMgr<ScalarType,MV,OP>::description() const
01672 {
01673 std::ostringstream oss;
01674 oss << "Belos::GCRODRSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
01675 oss << "{";
01676 oss << "Ortho Type='"<<orthoType_;
01677 oss << ", Num Blocks=" <<numBlocks_<< ", Max Restarts=" << maxRestarts_;
01678 oss << "}";
01679 return oss.str();
01680 }
01681
01682 }
01683
01684 #endif