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