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