Belos Version of the Day
Loading...
Searching...
No Matches
BelosBlockGmresSolMgr.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Belos: Block Linear Solvers Package
4//
5// Copyright 2004-2016 NTESS and the Belos contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef BELOS_BLOCK_GMRES_SOLMGR_HPP
11#define BELOS_BLOCK_GMRES_SOLMGR_HPP
12
17#include "BelosConfigDefs.hpp"
18#include "BelosTypes.hpp"
19
22
34#ifdef BELOS_TEUCHOS_TIME_MONITOR
35#include "Teuchos_TimeMonitor.hpp"
36#endif
37
54namespace Belos {
55
57
58
68
78
95template<class ScalarType, class MV, class OP>
96class BlockGmresSolMgr : public SolverManager<ScalarType,MV,OP> {
97
98private:
101 typedef Teuchos::ScalarTraits<ScalarType> SCT;
102 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
103 typedef Teuchos::ScalarTraits<MagnitudeType> MT;
104
105public:
106
108
109
116
138 const Teuchos::RCP<Teuchos::ParameterList> &pl );
139
141 virtual ~BlockGmresSolMgr() {};
142
144 Teuchos::RCP<SolverManager<ScalarType, MV, OP> > clone () const override {
145 return Teuchos::rcp(new BlockGmresSolMgr<ScalarType,MV,OP>);
146 }
148
150
151
155 return *problem_;
156 }
157
160 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const override;
161
164 Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const override { return params_; }
165
171 Teuchos::Array<Teuchos::RCP<Teuchos::Time> > getTimers() const {
172 return Teuchos::tuple(timerSolve_);
173 }
174
185 MagnitudeType achievedTol() const override {
186 return achievedTol_;
187 }
188
190 int getNumIters() const override {
191 return numIters_;
192 }
193
197 bool isLOADetected() const override { return loaDetected_; }
198
200
202
203
205 void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) override { problem_ = problem; isSTSet_ = false; needsIterRebuild_ = true; }
206
208 void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params ) override;
209
211 void setDebugStatusTest( const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &debugStatusTest ) override;
212
214
216
217
221 void reset( const ResetType type ) override { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
223
225
226
244 ReturnType solve() override;
245
247
250
257 void
258 describe (Teuchos::FancyOStream& out,
259 const Teuchos::EVerbosityLevel verbLevel =
260 Teuchos::Describable::verbLevel_default) const override;
261
263 std::string description () const override;
264
266
267private:
268
269 // Method for checking current status test against defined linear problem.
270 bool checkStatusTest();
271
272 // Linear problem.
273 Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
274
275 // Output manager.
276 Teuchos::RCP<OutputManager<ScalarType> > printer_;
277 Teuchos::RCP<std::ostream> outputStream_;
278
279 // Status test.
280 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > debugStatusTest_;
281 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
282 Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
283 Teuchos::RCP<StatusTest<ScalarType,MV,OP> > convTest_;
284 Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > expConvTest_, impConvTest_;
285 Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
286
287 // Orthogonalization manager.
288 Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
289
290 // Current parameter list.
291 Teuchos::RCP<Teuchos::ParameterList> params_;
292
293 // Default solver values.
294 static constexpr int maxRestarts_default_ = 20;
295 static constexpr int maxIters_default_ = 1000;
296 static constexpr bool adaptiveBlockSize_default_ = true;
297 static constexpr bool showMaxResNormOnly_default_ = false;
298 static constexpr bool flexibleGmres_default_ = false;
299 static constexpr bool keepHessenberg_default_ = false;
300 static constexpr bool expResTest_default_ = false;
301 static constexpr int blockSize_default_ = 1;
302 static constexpr int numBlocks_default_ = 300;
303 static constexpr int verbosity_default_ = Belos::Errors;
304 static constexpr int outputStyle_default_ = Belos::General;
305 static constexpr int outputFreq_default_ = -1;
306 static constexpr const char * impResScale_default_ = "Norm of Preconditioned Initial Residual";
307 static constexpr const char * expResScale_default_ = "Norm of Initial Residual";
308 static constexpr const char * label_default_ = "Belos";
309 static constexpr const char * orthoType_default_ = "ICGS";
310
311 // Current solver values.
312 MagnitudeType convtol_, orthoKappa_, achievedTol_;
313 int maxRestarts_, maxIters_, numIters_;
314 int blockSize_, numBlocks_, verbosity_, outputStyle_, outputFreq_;
315 bool adaptiveBlockSize_, showMaxResNormOnly_, isFlexible_, keepHessenberg_, expResTest_;
316 std::string orthoType_;
317 std::string impResScale_, expResScale_;
318
319 // Timers.
320 std::string label_;
321 Teuchos::RCP<Teuchos::Time> timerSolve_;
322
323 // Cached iterator (reused across solve() calls to avoid reallocating Krylov subspace).
324 Teuchos::RCP<GmresIteration<ScalarType,MV,OP> > block_gmres_iter_;
325
326 // Internal state variables.
327 bool isSet_, isSTSet_;
328 bool needsIterRebuild_;
329 bool loaDetected_;
330};
331
332
333// Empty Constructor
334template<class ScalarType, class MV, class OP>
336 outputStream_(Teuchos::rcpFromRef(std::cout)),
337 convtol_(DefaultSolverParameters::convTol),
338 orthoKappa_(DefaultSolverParameters::orthoKappa),
339 achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
340 maxRestarts_(maxRestarts_default_),
341 maxIters_(maxIters_default_),
342 numIters_(0),
343 blockSize_(blockSize_default_),
344 numBlocks_(numBlocks_default_),
345 verbosity_(verbosity_default_),
346 outputStyle_(outputStyle_default_),
347 outputFreq_(outputFreq_default_),
348 adaptiveBlockSize_(adaptiveBlockSize_default_),
349 showMaxResNormOnly_(showMaxResNormOnly_default_),
350 isFlexible_(flexibleGmres_default_),
351 keepHessenberg_(keepHessenberg_default_),
352 expResTest_(expResTest_default_),
353 orthoType_(orthoType_default_),
354 impResScale_(impResScale_default_),
355 expResScale_(expResScale_default_),
356 label_(label_default_),
357 isSet_(false),
358 isSTSet_(false),
359 needsIterRebuild_(true),
360 loaDetected_(false)
361{}
362
363
364// Basic Constructor
365template<class ScalarType, class MV, class OP>
368 const Teuchos::RCP<Teuchos::ParameterList> &pl) :
369 problem_(problem),
370 outputStream_(Teuchos::rcpFromRef(std::cout)),
371 convtol_(DefaultSolverParameters::convTol),
372 orthoKappa_(DefaultSolverParameters::orthoKappa),
373 achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
374 maxRestarts_(maxRestarts_default_),
375 maxIters_(maxIters_default_),
376 numIters_(0),
377 blockSize_(blockSize_default_),
378 numBlocks_(numBlocks_default_),
379 verbosity_(verbosity_default_),
380 outputStyle_(outputStyle_default_),
381 outputFreq_(outputFreq_default_),
382 adaptiveBlockSize_(adaptiveBlockSize_default_),
383 showMaxResNormOnly_(showMaxResNormOnly_default_),
384 isFlexible_(flexibleGmres_default_),
385 keepHessenberg_(keepHessenberg_default_),
386 expResTest_(expResTest_default_),
387 orthoType_(orthoType_default_),
388 impResScale_(impResScale_default_),
389 expResScale_(expResScale_default_),
390 label_(label_default_),
391 isSet_(false),
392 isSTSet_(false),
393 needsIterRebuild_(true),
394 loaDetected_(false)
395{
396
397 TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null, std::invalid_argument, "Problem not given to solver manager.");
398
399 // If the parameter list pointer is null, then set the current parameters to the default parameter list.
400 if ( !is_null(pl) ) {
401 setParameters( pl );
402 }
403
404}
405
406
407template<class ScalarType, class MV, class OP>
408Teuchos::RCP<const Teuchos::ParameterList>
410{
411 static Teuchos::RCP<const Teuchos::ParameterList> validPL;
412 if (is_null(validPL)) {
413 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
414
415 // The static_cast is to resolve an issue with older clang versions which
416 // would cause the constexpr to link fail. With c++17 the problem is resolved.
417 pl->set("Convergence Tolerance", static_cast<MagnitudeType>(DefaultSolverParameters::convTol),
418 "The relative residual tolerance that needs to be achieved by the\n"
419 "iterative solver in order for the linear system to be declared converged." );
420 pl->set("Maximum Restarts", static_cast<int>(maxRestarts_default_),
421 "The maximum number of restarts allowed for each\n"
422 "set of RHS solved.");
423 pl->set("Maximum Iterations", static_cast<int>(maxIters_default_),
424 "The maximum number of block iterations allowed for each\n"
425 "set of RHS solved.");
426 pl->set("Num Blocks", static_cast<int>(numBlocks_default_),
427 "The maximum number of blocks allowed in the Krylov subspace\n"
428 "for each set of RHS solved.");
429 pl->set("Block Size", static_cast<int>(blockSize_default_),
430 "The number of vectors in each block. This number times the\n"
431 "number of blocks is the total Krylov subspace dimension.");
432 pl->set("Adaptive Block Size", static_cast<bool>(adaptiveBlockSize_default_),
433 "Whether the solver manager should adapt the block size\n"
434 "based on the number of RHS to solve.");
435 pl->set("Verbosity", static_cast<int>(verbosity_default_),
436 "What type(s) of solver information should be outputted\n"
437 "to the output stream.");
438 pl->set("Output Style", static_cast<int>(outputStyle_default_),
439 "What style is used for the solver information outputted\n"
440 "to the output stream.");
441 pl->set("Output Frequency", static_cast<int>(outputFreq_default_),
442 "How often convergence information should be outputted\n"
443 "to the output stream.");
444 pl->set("Output Stream", Teuchos::rcpFromRef(std::cout),
445 "A reference-counted pointer to the output stream where all\n"
446 "solver output is sent.");
447 pl->set("Show Maximum Residual Norm Only", static_cast<bool>(showMaxResNormOnly_default_),
448 "When convergence information is printed, only show the maximum\n"
449 "relative residual norm when the block size is greater than one.");
450 pl->set("Flexible Gmres", static_cast<bool>(flexibleGmres_default_),
451 "Whether the solver manager should use the flexible variant\n"
452 "of GMRES.");
453 pl->set("Keep Hessenberg", static_cast<bool>(keepHessenberg_default_),
454 "Whether the raw upper Hessenberg matrix should be stored separately\n"
455 "from the QR-factored least squares system. Useful for harmonic Ritz\n"
456 "pair computation from the GMRES iteration state.");
457 pl->set("Explicit Residual Test", static_cast<bool>(expResTest_default_),
458 "Whether the explicitly computed residual should be used in the convergence test.");
459 pl->set("Implicit Residual Scaling", static_cast<const char *>(impResScale_default_),
460 "The type of scaling used in the implicit residual convergence test.");
461 pl->set("Explicit Residual Scaling", static_cast<const char *>(expResScale_default_),
462 "The type of scaling used in the explicit residual convergence test.");
463 pl->set("Timer Label", static_cast<const char *>(label_default_),
464 "The string to use as a prefix for the timer labels.");
465 pl->set("Orthogonalization", static_cast<const char *>(orthoType_default_),
466 "The type of orthogonalization to use: DGKS, ICGS, or IMGS.");
467 pl->set("Orthogonalization Constant",static_cast<MagnitudeType>(DefaultSolverParameters::orthoKappa),
468 "The constant used by DGKS orthogonalization to determine\n"
469 "whether another step of classical Gram-Schmidt is necessary.");
470 validPL = pl;
471 }
472 return validPL;
473}
474
475
476template<class ScalarType, class MV, class OP>
477void BlockGmresSolMgr<ScalarType,MV,OP>::setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params )
478{
479
480 // Create the internal parameter list if ones doesn't already exist.
481 if (params_ == Teuchos::null) {
482 params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
483 }
484 else {
485 params->validateParameters(*getValidParameters());
486 }
487
488 // Check for maximum number of restarts
489 if (params->isParameter("Maximum Restarts")) {
490 maxRestarts_ = params->get("Maximum Restarts",maxRestarts_default_);
491
492 // Update parameter in our list.
493 params_->set("Maximum Restarts", maxRestarts_);
494 }
495
496 // Check for maximum number of iterations
497 if (params->isParameter("Maximum Iterations")) {
498 maxIters_ = params->get("Maximum Iterations",maxIters_default_);
499
500 // Update parameter in our list and in status test.
501 params_->set("Maximum Iterations", maxIters_);
502 if (maxIterTest_!=Teuchos::null)
503 maxIterTest_->setMaxIters( maxIters_ );
504 }
505
506 // Check for blocksize
507 if (params->isParameter("Block Size")) {
508 blockSize_ = params->get("Block Size",blockSize_default_);
509 TEUCHOS_TEST_FOR_EXCEPTION(blockSize_ <= 0, std::invalid_argument,
510 "Belos::BlockGmresSolMgr: \"Block Size\" must be strictly positive.");
511
512 // Update parameter in our list.
513 params_->set("Block Size", blockSize_);
514 }
515
516 // Check if the blocksize should be adaptive
517 if (params->isParameter("Adaptive Block Size")) {
518 adaptiveBlockSize_ = params->get("Adaptive Block Size",adaptiveBlockSize_default_);
519
520 // Update parameter in our list.
521 params_->set("Adaptive Block Size", adaptiveBlockSize_);
522 }
523
524 // Check for the maximum number of blocks.
525 if (params->isParameter("Num Blocks")) {
526 numBlocks_ = params->get("Num Blocks",numBlocks_default_);
527 TEUCHOS_TEST_FOR_EXCEPTION(numBlocks_ <= 0, std::invalid_argument,
528 "Belos::BlockGmresSolMgr: \"Num Blocks\" must be strictly positive.");
529
530 // Update parameter in our list.
531 params_->set("Num Blocks", numBlocks_);
532 }
533
534 // Check to see if the timer label changed.
535 if (params->isParameter("Timer Label")) {
536 std::string tempLabel = params->get("Timer Label", label_default_);
537
538 // Update parameter in our list, solver timer, and orthogonalization label
539 if (tempLabel != label_) {
540 label_ = tempLabel;
541 params_->set("Timer Label", label_);
542 std::string solveLabel = label_ + ": BlockGmresSolMgr total solve time";
543#ifdef BELOS_TEUCHOS_TIME_MONITOR
544 timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
545#endif
546 if (ortho_ != Teuchos::null) {
547 ortho_->setLabel( label_ );
548 }
549 }
550 }
551
552 // Determine whether the raw Hessenberg should be stored separately from R.
553 if (params->isParameter("Keep Hessenberg")) {
554 keepHessenberg_ = Teuchos::getParameter<bool>(*params,"Keep Hessenberg");
555 params_->set("Keep Hessenberg", keepHessenberg_);
556 }
557
558 // Determine whether this solver should be "flexible".
559 if (params->isParameter("Flexible Gmres")) {
560 isFlexible_ = Teuchos::getParameter<bool>(*params,"Flexible Gmres");
561 params_->set("Flexible Gmres", isFlexible_);
562 if (isFlexible_ && expResTest_) {
563 // Use an implicit convergence test if the Gmres solver is flexible
564 isSTSet_ = false;
565 }
566 }
567
568 // Check for a change in verbosity level
569 if (params->isParameter("Verbosity")) {
570 if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
571 verbosity_ = params->get("Verbosity", verbosity_default_);
572 } else {
573 verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
574 }
575
576 // Update parameter in our list.
577 params_->set("Verbosity", verbosity_);
578 if (printer_ != Teuchos::null)
579 printer_->setVerbosity(verbosity_);
580 }
581
582 // Check for a change in output style
583 if (params->isParameter("Output Style")) {
584 if (Teuchos::isParameterType<int>(*params,"Output Style")) {
585 outputStyle_ = params->get("Output Style", outputStyle_default_);
586 } else {
587 outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
588 }
589
590 // Reconstruct the convergence test if the explicit residual test is not being used.
591 params_->set("Output Style", outputStyle_);
592 if (outputTest_ != Teuchos::null) {
593 isSTSet_ = false;
594 }
595 }
596
597 // output stream
598 if (params->isParameter("Output Stream")) {
599 outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
600
601 // Update parameter in our list.
602 params_->set("Output Stream", outputStream_);
603 if (printer_ != Teuchos::null)
604 printer_->setOStream( outputStream_ );
605 }
606
607 // frequency level
608 if (verbosity_ & Belos::StatusTestDetails) {
609 if (params->isParameter("Output Frequency")) {
610 outputFreq_ = params->get("Output Frequency", outputFreq_default_);
611 }
612
613 // Update parameter in out list and output status test.
614 params_->set("Output Frequency", outputFreq_);
615 if (outputTest_ != Teuchos::null)
616 outputTest_->setOutputFrequency( outputFreq_ );
617 }
618
619 // Create output manager if we need to.
620 if (printer_ == Teuchos::null) {
621 printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
622 }
623
624 // Check if the orthogonalization changed.
625 bool changedOrthoType = false;
626 if (params->isParameter("Orthogonalization")) {
627 std::string tempOrthoType = params->get("Orthogonalization",orthoType_default_);
628 if (tempOrthoType != orthoType_) {
629 orthoType_ = tempOrthoType;
630 changedOrthoType = true;
631 }
632 }
633 params_->set("Orthogonalization", orthoType_);
634
635 // Check which orthogonalization constant to use.
636 if (params->isParameter("Orthogonalization Constant")) {
637 if (params->isType<MagnitudeType> ("Orthogonalization Constant")) {
638 orthoKappa_ = params->get ("Orthogonalization Constant",
639 static_cast<MagnitudeType> (DefaultSolverParameters::orthoKappa));
640 }
641 else {
642 orthoKappa_ = params->get ("Orthogonalization Constant",
644 }
645
646 // Update parameter in our list.
647 params_->set("Orthogonalization Constant",orthoKappa_);
648 if (orthoType_=="DGKS") {
649 if (orthoKappa_ > 0 && ortho_ != Teuchos::null && !changedOrthoType) {
650 Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
651 }
652 }
653 }
654
655 // Create orthogonalization manager if we need to.
656 if (ortho_ == Teuchos::null || changedOrthoType) {
658 Teuchos::RCP<Teuchos::ParameterList> paramsOrtho;
659 if (orthoType_=="DGKS" && orthoKappa_ > 0) {
660 paramsOrtho = Teuchos::rcp(new Teuchos::ParameterList());
661 paramsOrtho->set ("depTol", orthoKappa_ );
662 }
663
664 ortho_ = factory.makeMatOrthoManager (orthoType_, Teuchos::null, printer_, label_, paramsOrtho);
665 }
666
667 // Check for convergence tolerance
668 if (params->isParameter("Convergence Tolerance")) {
669 if (params->isType<MagnitudeType> ("Convergence Tolerance")) {
670 convtol_ = params->get ("Convergence Tolerance",
671 static_cast<MagnitudeType> (DefaultSolverParameters::convTol));
672 }
673 else {
674 convtol_ = params->get ("Convergence Tolerance", DefaultSolverParameters::convTol);
675 }
676
677 // Update parameter in our list and residual tests.
678 params_->set("Convergence Tolerance", convtol_);
679 if (impConvTest_ != Teuchos::null)
680 impConvTest_->setTolerance( convtol_ );
681 if (expConvTest_ != Teuchos::null)
682 expConvTest_->setTolerance( convtol_ );
683 }
684
685 // Check for a change in scaling, if so we need to build new residual tests.
686 if (params->isParameter("Implicit Residual Scaling")) {
687 std::string tempImpResScale = Teuchos::getParameter<std::string>( *params, "Implicit Residual Scaling" );
688
689 // Only update the scaling if it's different.
690 if (impResScale_ != tempImpResScale) {
692 impResScale_ = tempImpResScale;
693
694 // Update parameter in our list and residual tests
695 params_->set("Implicit Residual Scaling", impResScale_);
696 if (impConvTest_ != Teuchos::null) {
697 try {
698 impConvTest_->defineScaleForm( impResScaleType, Belos::TwoNorm );
699 }
700 catch (std::exception& e) {
701 // Make sure the convergence test gets constructed again.
702 isSTSet_ = false;
703 }
704 }
705 }
706 }
707
708 if (params->isParameter("Explicit Residual Scaling")) {
709 std::string tempExpResScale = Teuchos::getParameter<std::string>( *params, "Explicit Residual Scaling" );
710
711 // Only update the scaling if it's different.
712 if (expResScale_ != tempExpResScale) {
714 expResScale_ = tempExpResScale;
715
716 // Update parameter in our list and residual tests
717 params_->set("Explicit Residual Scaling", expResScale_);
718 if (expConvTest_ != Teuchos::null) {
719 try {
720 expConvTest_->defineScaleForm( expResScaleType, Belos::TwoNorm );
721 }
722 catch (std::exception& e) {
723 // Make sure the convergence test gets constructed again.
724 isSTSet_ = false;
725 }
726 }
727 }
728 }
729
730 if (params->isParameter("Explicit Residual Test")) {
731 expResTest_ = Teuchos::getParameter<bool>( *params,"Explicit Residual Test" );
732
733 // Reconstruct the convergence test if the explicit residual test is not being used.
734 params_->set("Explicit Residual Test", expResTest_);
735 if (expConvTest_ == Teuchos::null) {
736 isSTSet_ = false;
737 }
738 }
739
740 if (params->isParameter("Show Maximum Residual Norm Only")) {
741 showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");
742
743 // Update parameter in our list and residual tests
744 params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
745 if (impConvTest_ != Teuchos::null)
746 impConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
747 if (expConvTest_ != Teuchos::null)
748 expConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
749 }
750
751
752 // Create the timer if we need to.
753 if (timerSolve_ == Teuchos::null) {
754 std::string solveLabel = label_ + ": BlockGmresSolMgr total solve time";
755#ifdef BELOS_TEUCHOS_TIME_MONITOR
756 timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
757#endif
758 }
759
760 // Inform the solver manager that the current parameters were set.
761 // The cached iterator must be rebuilt since printer_, outputTest_, and ortho_ may have changed.
762 needsIterRebuild_ = true;
763 isSet_ = true;
764}
765
766// Check the status test versus the defined linear problem
767template<class ScalarType, class MV, class OP>
769
770 typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
773
774 // Basic test checks maximum iterations and native residual.
775 maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
776
777 // Perform sanity checking for flexible Gmres here.
778 // NOTE: If the user requests that the solver manager use flexible GMRES, but there is no right preconditioner, don't use flexible GMRES.
779 // Throw an error is the user provided a left preconditioner, as that is inconsistent with flexible GMRES.
780 if (isFlexible_ && Teuchos::is_null(problem_->getRightPrec())) {
781 isFlexible_ = false;
782 params_->set("Flexible Gmres", isFlexible_);
783
784 // If the user specified the preconditioner as a left preconditioner, throw an error.
785 TEUCHOS_TEST_FOR_EXCEPTION( !Teuchos::is_null(problem_->getLeftPrec()),BlockGmresSolMgrLinearProblemFailure,
786 "Belos::BlockGmresSolMgr::solve(): Linear problem has a left preconditioner, not a right preconditioner, which is required for flexible GMRES.");
787 }
788
789 // If there is a left preconditioner, we create a combined status test that checks the implicit
790 // and then explicit residual norm to see if we have convergence.
791 if (!Teuchos::is_null(problem_->getLeftPrec()) && !isFlexible_) {
792 expResTest_ = true;
793 }
794
795 if (expResTest_) {
796
797 // Implicit residual test, using the native residual to determine if convergence was achieved.
798 Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
799 Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
800 tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
801 tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
802 impConvTest_ = tmpImpConvTest;
803
804 // Explicit residual test once the native residual is below the tolerance
805 Teuchos::RCP<StatusTestGenResNorm_t> tmpExpConvTest =
806 Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
807 tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit, Belos::TwoNorm );
808 tmpExpConvTest->defineScaleForm( convertStringToScaleType(expResScale_), Belos::TwoNorm );
809 tmpExpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
810 expConvTest_ = tmpExpConvTest;
811
812 // The convergence test is a combination of the "cheap" implicit test and explicit test.
813 convTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
814 }
815 else {
816
817 if (isFlexible_) {
818 // Implicit residual test, using the native residual to determine if convergence was achieved.
819 Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
820 Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
821 tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
822 tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
823 impConvTest_ = tmpImpConvTest;
824 }
825 else {
826 // Implicit residual test, using the native residual to determine if convergence was achieved.
827 // Use test that checks for loss of accuracy.
828 Teuchos::RCP<StatusTestImpResNorm_t> tmpImpConvTest =
829 Teuchos::rcp( new StatusTestImpResNorm_t( convtol_ ) );
830 tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
831 tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
832 impConvTest_ = tmpImpConvTest;
833 }
834
835 // Set the explicit and total convergence test to this implicit test that checks for accuracy loss.
836 expConvTest_ = impConvTest_;
837 convTest_ = impConvTest_;
838 }
839
840 // Create the status test.
841 sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
842
843 // Add debug status test, if one is provided by the user
844 if (nonnull(debugStatusTest_) ) {
845 // Add debug convergence test
846 Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(sTest_)->addStatusTest( debugStatusTest_ );
847 }
848
849 // Create the status test output class.
850 // This class manages and formats the output from the status test.
852 outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
853
854 // Set the solver string for the output test
855 std::string solverDesc = " Block Gmres ";
856 if (isFlexible_)
857 solverDesc = "Flexible" + solverDesc;
858 outputTest_->setSolverDesc( solverDesc );
859
860 // The status test is now set.
861 isSTSet_ = true;
862
863 return false;
864}
865
866template<class ScalarType, class MV, class OP>
873
874
875// solve()
876template<class ScalarType, class MV, class OP>
878
879 // Set the current parameters if they were not set before.
880 // NOTE: This may occur if the user generated the solver manager with the default constructor and
881 // then didn't set any parameters using setParameters().
882 if (!isSet_) {
883 setParameters(Teuchos::parameterList(*getValidParameters()));
884 }
885
887 "Belos::BlockGmresSolMgr::solve(): Linear problem is not a valid object.");
888
890 "Belos::BlockGmresSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
891
892 if (!isSTSet_ || (!expResTest_ && !Teuchos::is_null(problem_->getLeftPrec())) ) {
894 "Belos::BlockGmresSolMgr::solve(): Linear problem and requested status tests are incompatible.");
895 }
896
897 // Create indices for the linear systems to be solved.
898 int startPtr = 0;
899 int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
900 int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
901
902 std::vector<int> currIdx;
903 // If an adaptive block size is allowed then only the linear systems that need to be solved are solved.
904 // Otherwise, the index set is generated that informs the linear problem that some linear systems are augmented.
905 if ( adaptiveBlockSize_ ) {
906 blockSize_ = numCurrRHS;
907 currIdx.resize( numCurrRHS );
908 for (int i=0; i<numCurrRHS; ++i)
909 { currIdx[i] = startPtr+i; }
910
911 }
912 else {
913 currIdx.resize( blockSize_ );
914 for (int i=0; i<numCurrRHS; ++i)
915 { currIdx[i] = startPtr+i; }
916 for (int i=numCurrRHS; i<blockSize_; ++i)
917 { currIdx[i] = -1; }
918 }
919
920 // Inform the linear problem of the current linear system to solve.
921 problem_->setLSIndex( currIdx );
922
924 // Parameter list
925 Teuchos::ParameterList plist;
926 plist.set("Block Size",blockSize_);
927 plist.set("Keep Hessenberg",keepHessenberg_);
928
929 ptrdiff_t dim = MVT::GetGlobalLength( *(problem_->getRHS()) );
930 if (blockSize_*static_cast<ptrdiff_t>(numBlocks_) > dim) {
931 int tmpNumBlocks = 0;
932 if (blockSize_ == 1)
933 tmpNumBlocks = dim / blockSize_; // Allow for a good breakdown.
934 else
935 tmpNumBlocks = ( dim - blockSize_) / blockSize_; // Allow for restarting.
936 printer_->stream(Warnings) <<
937 "Belos::BlockGmresSolMgr::solve(): Warning! Requested Krylov subspace dimension is larger than operator dimension!"
938 << std::endl << " The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << tmpNumBlocks << std::endl;
939 plist.set("Num Blocks",tmpNumBlocks);
940 }
941 else
942 plist.set("Num Blocks",numBlocks_);
943
944 // Reset the status test.
945 outputTest_->reset();
946 loaDetected_ = false;
947
948 // Assume convergence is achieved, then let any failed convergence set this to false.
949 bool isConverged = true;
950
952 // BlockGmres solver
953
954 if (needsIterRebuild_) {
955 if (isFlexible_)
956 block_gmres_iter_ = Teuchos::rcp( new BlockFGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
957 else
958 block_gmres_iter_ = Teuchos::rcp( new BlockGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
959 needsIterRebuild_ = false;
960 }
961 Teuchos::RCP<GmresIteration<ScalarType,MV,OP> > &block_gmres_iter = block_gmres_iter_;
962
963 // Enter solve() iterations
964 {
965#ifdef BELOS_TEUCHOS_TIME_MONITOR
966 Teuchos::TimeMonitor slvtimer(*timerSolve_);
967#endif
968
969 while ( numRHS2Solve > 0 ) {
970
971 // Set the current number of blocks and blocksize with the Gmres iteration.
972 if (blockSize_*numBlocks_ > dim) {
973 int tmpNumBlocks = 0;
974 if (blockSize_ == 1)
975 tmpNumBlocks = dim / blockSize_; // Allow for a good breakdown.
976 else
977 tmpNumBlocks = ( dim - blockSize_) / blockSize_; // Allow for restarting.
978 block_gmres_iter->setSize( blockSize_, tmpNumBlocks );
979 }
980 else
981 block_gmres_iter->setSize( blockSize_, numBlocks_ );
982
983 // Reset the number of iterations.
984 block_gmres_iter->resetNumIters();
985
986 // Reset the number of calls that the status test output knows about.
987 outputTest_->resetNumCalls();
988
989 // Create the first block in the current Krylov basis.
990 Teuchos::RCP<MV> V_0;
991 if (isFlexible_) {
992 // Load the correct residual if the system is augmented
993 if (currIdx[blockSize_-1] == -1) {
994 V_0 = MVT::Clone( *(problem_->getInitResVec()), blockSize_ );
995 problem_->computeCurrResVec( &*V_0 );
996 }
997 else {
998 V_0 = MVT::CloneCopy( *(problem_->getInitResVec()), currIdx );
999 }
1000 }
1001 else {
1002 // Load the correct residual if the system is augmented
1003 if (currIdx[blockSize_-1] == -1) {
1004 V_0 = MVT::Clone( *(problem_->getInitPrecResVec()), blockSize_ );
1005 problem_->computeCurrPrecResVec( &*V_0 );
1006 }
1007 else {
1008 V_0 = MVT::CloneCopy( *(problem_->getInitPrecResVec()), currIdx );
1009 }
1010 }
1011
1012 // Get a matrix to hold the orthonormalization coefficients.
1013 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > z_0 =
1014 Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( blockSize_, blockSize_ ) );
1015
1016 // Orthonormalize the new V_0
1017 int rank = ortho_->normalize( *V_0, z_0 );
1019 "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors.");
1020
1021 // Set the new state and initialize the solver.
1023 newstate.V = V_0;
1024 newstate.z = z_0;
1025 newstate.curDim = 0;
1026 block_gmres_iter->initializeGmres(newstate);
1027 int numRestarts = 0;
1028
1029 while(1) {
1030 // tell block_gmres_iter to iterate
1031 try {
1032 block_gmres_iter->iterate();
1033
1035 //
1036 // check convergence first
1037 //
1039 if ( convTest_->getStatus() == Passed ) {
1040 if ( expConvTest_->getLOADetected() ) {
1041 // we don't have convergence
1042 loaDetected_ = true;
1043 printer_->stream(Warnings) <<
1044 "Belos::BlockGmresSolMgr::solve(): Warning! Solver has experienced a loss of accuracy!" << std::endl;
1045 isConverged = false;
1046 }
1047 break; // break from while(1){block_gmres_iter->iterate()}
1048 }
1050 //
1051 // check for maximum iterations
1052 //
1054 else if ( maxIterTest_->getStatus() == Passed ) {
1055 // we don't have convergence
1056 isConverged = false;
1057 break; // break from while(1){block_gmres_iter->iterate()}
1058 }
1060 //
1061 // check for restarting, i.e. the subspace is full
1062 //
1064 else if ( block_gmres_iter->getCurSubspaceDim() == block_gmres_iter->getMaxSubspaceDim() ) {
1065
1066 if ( numRestarts >= maxRestarts_ ) {
1067 isConverged = false;
1068 break; // break from while(1){block_gmres_iter->iterate()}
1069 }
1070 numRestarts++;
1071
1072 printer_->stream(Debug) << " Performing restart number " << numRestarts << " of " << maxRestarts_ << std::endl << std::endl;
1073
1074 // Update the linear problem.
1075 Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1076 if (isFlexible_) {
1077 // Update the solution manually, since the preconditioning doesn't need to be undone.
1078 Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1079 MVT::MvAddMv( 1.0, *curX, 1.0, *update, *curX );
1080 }
1081 else
1082 problem_->updateSolution( update, true );
1083
1084 // Get the state.
1086
1087 // Compute the restart std::vector.
1088 // Get a view of the current Krylov basis.
1089 V_0 = MVT::Clone( *(oldState.V), blockSize_ );
1090 if (isFlexible_)
1091 problem_->computeCurrResVec( &*V_0 );
1092 else
1093 problem_->computeCurrPrecResVec( &*V_0 );
1094
1095 // Get a view of the first block of the Krylov basis.
1096 z_0 = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( blockSize_, blockSize_ ) );
1097
1098 // Orthonormalize the new V_0
1099 rank = ortho_->normalize( *V_0, z_0 );
1101 "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors after restart.");
1102
1103 // Set the new state and initialize the solver.
1104 newstate.V = V_0;
1105 newstate.z = z_0;
1106 newstate.curDim = 0;
1107 block_gmres_iter->initializeGmres(newstate);
1108
1109 } // end of restarting
1110
1112 //
1113 // we returned from iterate(), but none of our status tests Passed.
1114 // something is wrong, and it is probably our fault.
1115 //
1117
1118 else {
1119 TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
1120 "Belos::BlockGmresSolMgr::solve(): Invalid return from BlockGmresIter::iterate().");
1121 }
1122 }
1123 catch (const GmresIterationOrthoFailure &e) {
1124 // If the block size is not one, it's not considered a lucky breakdown.
1125 if (blockSize_ != 1) {
1126 printer_->stream(Errors) << "Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
1127 << block_gmres_iter->getNumIters() << std::endl
1128 << e.what() << std::endl;
1129 if (convTest_->getStatus() != Passed)
1130 isConverged = false;
1131 break;
1132 }
1133 else {
1134 // If the block size is one, try to recover the most recent least-squares solution
1135 block_gmres_iter->updateLSQR( block_gmres_iter->getCurSubspaceDim() );
1136
1137 // Check to see if the most recent least-squares solution yielded convergence.
1138 sTest_->checkStatus( &*block_gmres_iter );
1139 if (convTest_->getStatus() != Passed)
1140 isConverged = false;
1141 break;
1142 }
1143 }
1144 catch (const StatusTestNaNError& e) {
1145 // A NaN was detected in the solver. Set the solution to zero and return unconverged.
1146 achievedTol_ = MT::one();
1147 Teuchos::RCP<MV> X = problem_->getLHS();
1148 MVT::MvInit( *X, SCT::zero() );
1149 printer_->stream(Warnings) << "Belos::BlockGmresSolMgr::solve(): Warning! NaN has been detected!"
1150 << std::endl;
1151 return Unconverged;
1152 }
1153 catch (const std::exception &e) {
1154 printer_->stream(Errors) << "Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
1155 << block_gmres_iter->getNumIters() << std::endl
1156 << e.what() << std::endl;
1157 throw;
1158 }
1159 }
1160
1161 // Compute the current solution.
1162 // Update the linear problem.
1163 if (isFlexible_) {
1164 // Update the solution manually, since the preconditioning doesn't need to be undone.
1165 Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1166 Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1167 // Update the solution only if there is a valid update from the iteration
1168 if (update != Teuchos::null)
1169 MVT::MvAddMv( 1.0, *curX, 1.0, *update, *curX );
1170 }
1171 else {
1172 // Attempt to get the current solution from the residual status test, if it has one.
1173 if ( !Teuchos::is_null(expConvTest_->getSolution()) ) {
1174 Teuchos::RCP<MV> newX = expConvTest_->getSolution();
1175 Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1176 MVT::MvAddMv( 0.0, *newX, 1.0, *newX, *curX );
1177 }
1178 else {
1179 Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1180 problem_->updateSolution( update, true );
1181 }
1182 }
1183
1184 // Inform the linear problem that we are finished with this block linear system.
1185 problem_->setCurrLS();
1186
1187 // Update indices for the linear systems to be solved.
1190 if ( numRHS2Solve > 0 ) {
1191 numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1192
1193 if ( adaptiveBlockSize_ ) {
1194 blockSize_ = numCurrRHS;
1195 currIdx.resize( numCurrRHS );
1196 for (int i=0; i<numCurrRHS; ++i)
1197 { currIdx[i] = startPtr+i; }
1198 }
1199 else {
1200 currIdx.resize( blockSize_ );
1201 for (int i=0; i<numCurrRHS; ++i)
1202 { currIdx[i] = startPtr+i; }
1203 for (int i=numCurrRHS; i<blockSize_; ++i)
1204 { currIdx[i] = -1; }
1205 }
1206 // Set the next indices.
1207 problem_->setLSIndex( currIdx );
1208 }
1209 else {
1210 currIdx.resize( numRHS2Solve );
1211 }
1212
1213 }// while ( numRHS2Solve > 0 )
1214
1215 }
1216
1217 // print final summary
1218 sTest_->print( printer_->stream(FinalSummary) );
1219
1220 // print timing information
1221#ifdef BELOS_TEUCHOS_TIME_MONITOR
1222 // Calling summarize() can be expensive, so don't call unless the
1223 // user wants to print out timing details. summarize() will do all
1224 // the work even if it's passed a "black hole" output stream.
1225 if (verbosity_ & TimingDetails)
1226 Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
1227#endif
1228
1229 // get iteration information for this solve
1230 numIters_ = maxIterTest_->getNumIters();
1231
1232 // Save the convergence test value ("achieved tolerance") for this
1233 // solve. This requires a bit more work than for BlockCGSolMgr,
1234 // since for this solver, convTest_ may either be a single residual
1235 // norm test, or a combination of two residual norm tests. In the
1236 // latter case, the master convergence test convTest_ is a SEQ combo
1237 // of the implicit resp. explicit tests. If the implicit test never
1238 // passes, then the explicit test won't ever be executed. This
1239 // manifests as expConvTest_->getTestValue()->size() < 1. We deal
1240 // with this case by using the values returned by
1241 // impConvTest_->getTestValue().
1242 {
1243 // We'll fetch the vector of residual norms one way or the other.
1244 const std::vector<MagnitudeType>* pTestValues = NULL;
1245 if (expResTest_) {
1246 pTestValues = expConvTest_->getTestValue();
1247 if (pTestValues == NULL || pTestValues->size() < 1) {
1248 pTestValues = impConvTest_->getTestValue();
1249 }
1250 }
1251 else {
1252 // Only the implicit residual norm test is being used.
1253 pTestValues = impConvTest_->getTestValue();
1254 }
1255 TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
1256 "Belos::BlockGmresSolMgr::solve(): The implicit convergence test's "
1257 "getTestValue() method returned NULL. Please report this bug to the "
1258 "Belos developers.");
1259 TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
1260 "Belos::BlockGmresSolMgr::solve(): The implicit convergence test's "
1261 "getTestValue() method returned a vector of length zero. Please report "
1262 "this bug to the Belos developers.");
1263
1264 // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
1265 // achieved tolerances for all vectors in the current solve(), or
1266 // just for the vectors from the last deflation?
1267 achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1268 }
1269
1270 if (!isConverged || loaDetected_) {
1271 return Unconverged; // return from BlockGmresSolMgr::solve()
1272 }
1273 return Converged; // return from BlockGmresSolMgr::solve()
1274}
1275
1276
1277template<class ScalarType, class MV, class OP>
1279{
1280 std::ostringstream out;
1281 out << "\"Belos::BlockGmresSolMgr\": {";
1282 if (this->getObjectLabel () != "") {
1283 out << "Label: " << this->getObjectLabel () << ", ";
1284 }
1285 out << "Flexible: " << (isFlexible_ ? "true" : "false")
1286 << ", Num Blocks: " << numBlocks_
1287 << ", Maximum Iterations: " << maxIters_
1288 << ", Maximum Restarts: " << maxRestarts_
1289 << ", Convergence Tolerance: " << convtol_
1290 << "}";
1291 return out.str ();
1292}
1293
1294
1295template<class ScalarType, class MV, class OP>
1296void
1298describe (Teuchos::FancyOStream &out,
1299 const Teuchos::EVerbosityLevel verbLevel) const
1300{
1301 using Teuchos::TypeNameTraits;
1302 using Teuchos::VERB_DEFAULT;
1303 using Teuchos::VERB_NONE;
1304 using Teuchos::VERB_LOW;
1305 // using Teuchos::VERB_MEDIUM;
1306 // using Teuchos::VERB_HIGH;
1307 // using Teuchos::VERB_EXTREME;
1308 using std::endl;
1309
1310 // Set default verbosity if applicable.
1311 const Teuchos::EVerbosityLevel vl =
1313
1314 if (vl != VERB_NONE) {
1315 Teuchos::OSTab tab0 (out);
1316
1317 out << "\"Belos::BlockGmresSolMgr\":" << endl;
1318 Teuchos::OSTab tab1 (out);
1319 out << "Template parameters:" << endl;
1320 {
1321 Teuchos::OSTab tab2 (out);
1322 out << "ScalarType: " << TypeNameTraits<ScalarType>::name () << endl
1323 << "MV: " << TypeNameTraits<MV>::name () << endl
1324 << "OP: " << TypeNameTraits<OP>::name () << endl;
1325 }
1326 if (this->getObjectLabel () != "") {
1327 out << "Label: " << this->getObjectLabel () << endl;
1328 }
1329 out << "Flexible: " << (isFlexible_ ? "true" : "false") << endl
1330 << "Num Blocks: " << numBlocks_ << endl
1331 << "Maximum Iterations: " << maxIters_ << endl
1332 << "Maximum Restarts: " << maxRestarts_ << endl
1333 << "Convergence Tolerance: " << convtol_ << endl;
1334 }
1335}
1336
1337
1338} // end Belos namespace
1339
1340#endif /* BELOS_BLOCK_GMRES_SOLMGR_HPP */
Belos concrete class for performing the block, flexible GMRES iteration.
Belos concrete class for performing the block GMRES iteration.
Belos header file which uses auto-configuration information to include necessary C++ headers.
Pure virtual base class which augments the basic interface for a Gmres linear solver iteration.
Class which describes the linear problem to be solved by the iterative solver.
Class which manages the output and verbosity of the Belos solvers.
Pure virtual base class which describes the basic interface for a solver manager.
Belos::StatusTest for logically combining several status tests.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Belos::StatusTest for specifying an implicit residual norm stopping criteria that checks for loss of ...
Belos::StatusTest class for specifying a maximum number of iterations.
Virtual base class for StatusTest that printing status tests.
A factory class for generating StatusTestOutput objects.
Collection of types and exceptions used within the Belos solvers.
Parent class to all Belos exceptions.
Interface to Block GMRES and Flexible GMRES.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
int getNumIters() const override
Get the iteration count for the most recent call to solve().
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
Print the object with the given verbosity level to a FancyOStream.
void setDebugStatusTest(const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &debugStatusTest) override
Set a debug status test that will be checked at the same time as the top-level status test.
ReturnType solve() override
This method performs possibly repeated calls to the underlying linear solver's iterate() routine unti...
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
std::string description() const override
Return a one-line description of this object.
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
Set the parameters the solver manager should use to solve the linear problem.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Get a parameter list containing the valid parameters for this object.
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
BlockGmresSolMgr()
Empty constructor for BlockGmresSolMgr. This constructor takes no arguments and sets the default valu...
virtual ~BlockGmresSolMgr()
Destructor.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
void reset(const ResetType type) override
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
bool isLOADetected() const override
Return whether a loss of accuracy was detected by this solver during the most current solve.
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Get current linear problem being solved for in this object.
BlockGmresSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i....
BlockGmresSolMgrLinearProblemFailure(const std::string &what_arg)
BlockGmresSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate ortho...
BlockGmresSolMgrOrthoFailure(const std::string &what_arg)
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
Alternative run-time polymorphic interface for operators.
Operator()
Default constructor (does nothing).
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
@ StatusTestDetails
@ FinalSummary
@ TimingDetails
ReturnType
Whether the Belos solve converged for all linear systems.
@ Unconverged
ScaleType
The type of scaling to use on the residual norm value.
ResetType
How to reset the solver.
Default parameters common to most Belos solvers.
static const double convTol
Default convergence tolerance.
static const double orthoKappa
DGKS orthogonalization constant.

Generated for Belos by doxygen 1.9.8