Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosTFQMRSolMgr.hpp
Go to the documentation of this file.
1 //@HEADER
2 // ************************************************************************
3 //
4 // Belos: Block Linear Solvers Package
5 // Copyright 2004 Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 //@HEADER
41 
42 #ifndef BELOS_TFQMR_SOLMGR_HPP
43 #define BELOS_TFQMR_SOLMGR_HPP
44 
49 #include "BelosConfigDefs.hpp"
50 #include "BelosTypes.hpp"
51 
52 #include "BelosLinearProblem.hpp"
53 #include "BelosSolverManager.hpp"
54 
55 #include "BelosTFQMRIter.hpp"
58 #include "BelosStatusTestCombo.hpp"
60 #include "BelosOutputManager.hpp"
61 #ifdef BELOS_TEUCHOS_TIME_MONITOR
62 #include "Teuchos_TimeMonitor.hpp"
63 #endif
64 
78 namespace Belos {
79 
81 
82 
90  TFQMRSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
91  {}};
92 
99  class TFQMRSolMgrOrthoFailure : public BelosError {public:
100  TFQMRSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
101  {}};
102 
103  template<class ScalarType, class MV, class OP>
104  class TFQMRSolMgr : public SolverManager<ScalarType,MV,OP> {
105 
106  private:
112 
113  public:
114 
116 
117 
123  TFQMRSolMgr();
124 
143 
145  virtual ~TFQMRSolMgr() {};
146 
150  }
152 
154 
155 
156  const LinearProblem<ScalarType,MV,OP>& getProblem() const override {
157  return *problem_;
158  }
159 
163 
167 
174  return Teuchos::tuple(timerSolve_);
175  }
176 
182  MagnitudeType achievedTol() const override {
183  return achievedTol_;
184  }
185 
187  int getNumIters() const override {
188  return numIters_;
189  }
190 
198  bool isLOADetected() const override { return false; }
200 
202 
203 
205  void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) override { problem_ = problem; }
206 
208  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params ) override;
209 
211 
213 
218  void reset( const ResetType type ) override { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
220 
222 
223 
241  ReturnType solve() override;
242 
244 
246 
248  std::string description() const override;
250 
251  private:
252 
253  // Method for checking current status test against defined linear problem.
254  bool checkStatusTest();
255 
256  // Linear problem.
258 
259  // Output manager.
262 
263  // Status test.
269 
270  // Current parameter list.
272 
273  // Default solver values.
274  static constexpr int maxIters_default_ = 1000;
275  static constexpr bool expResTest_default_ = false;
276  static constexpr int verbosity_default_ = Belos::Errors;
277  static constexpr int outputStyle_default_ = Belos::General;
278  static constexpr int outputFreq_default_ = -1;
279  static constexpr const char * impResScale_default_ = "Norm of Preconditioned Initial Residual";
280  static constexpr const char * expResScale_default_ = "Norm of Initial Residual";
281  static constexpr const char * label_default_ = "Belos";
282  static constexpr std::ostream * outputStream_default_ = &std::cout;
283 
284  // Current solver values.
291 
292  // Timers.
293  std::string label_;
295 
296  // Internal state variables.
298  };
299 
300 
301 // Empty Constructor
302 template<class ScalarType, class MV, class OP>
304  outputStream_(Teuchos::rcp(outputStream_default_,false)),
305  convtol_(DefaultSolverParameters::convTol),
306  impTolScale_(DefaultSolverParameters::impTolScale),
307  achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
308  maxIters_(maxIters_default_),
309  numIters_(0),
310  verbosity_(verbosity_default_),
311  outputStyle_(outputStyle_default_),
312  outputFreq_(outputFreq_default_),
313  blockSize_(1),
314  expResTest_(expResTest_default_),
315  impResScale_(impResScale_default_),
316  expResScale_(expResScale_default_),
317  label_(label_default_),
318  isSet_(false),
319  isSTSet_(false)
320 {}
321 
322 
323 // Basic Constructor
324 template<class ScalarType, class MV, class OP>
328  problem_(problem),
329  outputStream_(Teuchos::rcp(outputStream_default_,false)),
330  convtol_(DefaultSolverParameters::convTol),
331  impTolScale_(DefaultSolverParameters::impTolScale),
332  achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
333  maxIters_(maxIters_default_),
334  numIters_(0),
335  verbosity_(verbosity_default_),
336  outputStyle_(outputStyle_default_),
337  outputFreq_(outputFreq_default_),
338  blockSize_(1),
339  expResTest_(expResTest_default_),
340  impResScale_(impResScale_default_),
341  expResScale_(expResScale_default_),
342  label_(label_default_),
343  isSet_(false),
344  isSTSet_(false)
345 {
346  TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null, std::invalid_argument, "Problem not given to solver manager.");
347 
348  // If the parameter list pointer is null, then set the current parameters to the default parameter list.
349  if ( !is_null(pl) ) {
350  setParameters( pl );
351  }
352 }
353 
354 template<class ScalarType, class MV, class OP>
356 {
357  // Create the internal parameter list if ones doesn't already exist.
358  if (params_ == Teuchos::null) {
359  params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
360  }
361  else {
362  params->validateParameters(*getValidParameters());
363  }
364 
365  // Check for maximum number of iterations
366  if (params->isParameter("Maximum Iterations")) {
367  maxIters_ = params->get("Maximum Iterations",maxIters_default_);
368 
369  // Update parameter in our list and in status test.
370  params_->set("Maximum Iterations", maxIters_);
371  if (maxIterTest_!=Teuchos::null)
372  maxIterTest_->setMaxIters( maxIters_ );
373  }
374 
375  // Check for blocksize
376  if (params->isParameter("Block Size")) {
377  blockSize_ = params->get("Block Size",1);
378  TEUCHOS_TEST_FOR_EXCEPTION(blockSize_ != 1, std::invalid_argument,
379  "Belos::TFQMRSolMgr: \"Block Size\" must be 1.");
380 
381  // Update parameter in our list.
382  params_->set("Block Size", blockSize_);
383  }
384 
385  // Check to see if the timer label changed.
386  if (params->isParameter("Timer Label")) {
387  std::string tempLabel = params->get("Timer Label", label_default_);
388 
389  // Update parameter in our list and solver timer
390  if (tempLabel != label_) {
391  label_ = tempLabel;
392  params_->set("Timer Label", label_);
393  std::string solveLabel = label_ + ": TFQMRSolMgr total solve time";
394 #ifdef BELOS_TEUCHOS_TIME_MONITOR
395  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
396 #endif
397  }
398  }
399 
400  // Check for a change in verbosity level
401  if (params->isParameter("Verbosity")) {
402  if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
403  verbosity_ = params->get("Verbosity", verbosity_default_);
404  } else {
405  verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
406  }
407 
408  // Update parameter in our list.
409  params_->set("Verbosity", verbosity_);
410  if (printer_ != Teuchos::null)
411  printer_->setVerbosity(verbosity_);
412  }
413 
414  // Check for a change in output style
415  if (params->isParameter("Output Style")) {
416  if (Teuchos::isParameterType<int>(*params,"Output Style")) {
417  outputStyle_ = params->get("Output Style", outputStyle_default_);
418  } else {
419  outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
420  }
421 
422  // Reconstruct the convergence test if the explicit residual test is not being used.
423  params_->set("Output Style", outputStyle_);
424  isSTSet_ = false;
425  }
426 
427  // output stream
428  if (params->isParameter("Output Stream")) {
429  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
430 
431  // Update parameter in our list.
432  params_->set("Output Stream", outputStream_);
433  if (printer_ != Teuchos::null)
434  printer_->setOStream( outputStream_ );
435  }
436 
437  // frequency level
438  if (verbosity_ & Belos::StatusTestDetails) {
439  if (params->isParameter("Output Frequency")) {
440  outputFreq_ = params->get("Output Frequency", outputFreq_default_);
441  }
442 
443  // Update parameter in out list and output status test.
444  params_->set("Output Frequency", outputFreq_);
445  if (outputTest_ != Teuchos::null)
446  outputTest_->setOutputFrequency( outputFreq_ );
447  }
448 
449  // Create output manager if we need to.
450  if (printer_ == Teuchos::null) {
451  printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
452  }
453 
454  // Check for convergence tolerance
455  if (params->isParameter("Convergence Tolerance")) {
456  if (params->isType<MagnitudeType> ("Convergence Tolerance")) {
457  convtol_ = params->get ("Convergence Tolerance",
458  static_cast<MagnitudeType> (DefaultSolverParameters::convTol));
459  }
460  else {
461  convtol_ = params->get ("Convergence Tolerance", DefaultSolverParameters::convTol);
462  }
463 
464  // Update parameter in our list.
465  params_->set("Convergence Tolerance", convtol_);
466  isSTSet_ = false;
467  }
468 
469  // Check for implicit residual scaling
470  if (params->isParameter("Implicit Tolerance Scale Factor")) {
471  if (params->isType<MagnitudeType> ("Implicit Tolerance Scale Factor")) {
472  impTolScale_ = params->get ("Implicit Tolerance Scale Factor",
473  static_cast<MagnitudeType> (DefaultSolverParameters::impTolScale));
474 
475  }
476  else {
477  impTolScale_ = params->get ("Implicit Tolerance Scale Factor",
479  }
480 
481  // Update parameter in our list.
482  params_->set("Implicit Tolerance Scale Factor", impTolScale_);
483  isSTSet_ = false;
484  }
485 
486  // Check for a change in scaling, if so we need to build new residual tests.
487  if (params->isParameter("Implicit Residual Scaling")) {
488  std::string tempImpResScale = Teuchos::getParameter<std::string>( *params, "Implicit Residual Scaling" );
489 
490  // Only update the scaling if it's different.
491  if (impResScale_ != tempImpResScale) {
492  impResScale_ = tempImpResScale;
493 
494  // Update parameter in our list and residual tests
495  params_->set("Implicit Residual Scaling", impResScale_);
496 
497  // Make sure the convergence test gets constructed again.
498  isSTSet_ = false;
499  }
500  }
501 
502  if (params->isParameter("Explicit Residual Scaling")) {
503  std::string tempExpResScale = Teuchos::getParameter<std::string>( *params, "Explicit Residual Scaling" );
504 
505  // Only update the scaling if it's different.
506  if (expResScale_ != tempExpResScale) {
507  expResScale_ = tempExpResScale;
508 
509  // Update parameter in our list and residual tests
510  params_->set("Explicit Residual Scaling", expResScale_);
511 
512  // Make sure the convergence test gets constructed again.
513  isSTSet_ = false;
514  }
515  }
516 
517  if (params->isParameter("Explicit Residual Test")) {
518  expResTest_ = Teuchos::getParameter<bool>( *params,"Explicit Residual Test" );
519 
520  // Reconstruct the convergence test if the explicit residual test is not being used.
521  params_->set("Explicit Residual Test", expResTest_);
522  if (expConvTest_ == Teuchos::null) {
523  isSTSet_ = false;
524  }
525  }
526 
527  // Create the timer if we need to.
528  if (timerSolve_ == Teuchos::null) {
529  std::string solveLabel = label_ + ": TFQMRSolMgr total solve time";
530 #ifdef BELOS_TEUCHOS_TIME_MONITOR
531  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
532 #endif
533  }
534 
535  // Inform the solver manager that the current parameters were set.
536  isSet_ = true;
537 }
538 
539 
540 // Check the status test versus the defined linear problem
541 template<class ScalarType, class MV, class OP>
543 
544  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
545  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestGenResNorm_t;
546 
547  // Basic test checks maximum iterations and native residual.
548  maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
549 
550  if (expResTest_) {
551 
552  // Implicit residual test, using the native residual to determine if convergence was achieved.
553  Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
554  Teuchos::rcp( new StatusTestGenResNorm_t( impTolScale_*convtol_ ) );
555  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
556  impConvTest_ = tmpImpConvTest;
557 
558  // Explicit residual test once the native residual is below the tolerance
559  Teuchos::RCP<StatusTestGenResNorm_t> tmpExpConvTest =
560  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
561  tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit, Belos::TwoNorm );
562  tmpExpConvTest->defineScaleForm( convertStringToScaleType(expResScale_), Belos::TwoNorm );
563  expConvTest_ = tmpExpConvTest;
564 
565  // The convergence test is a combination of the "cheap" implicit test and explicit test.
566  convTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
567  }
568  else {
569 
570  // Implicit residual test, using the native residual to determine if convergence was achieved.
571  Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
572  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
573  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
574  impConvTest_ = tmpImpConvTest;
575 
576  // Set the explicit and total convergence test to this implicit test that checks for accuracy loss.
577  expConvTest_ = impConvTest_;
578  convTest_ = impConvTest_;
579  }
580  sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
581 
582  // Create the status test output class.
583  // This class manages and formats the output from the status test.
584  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
585  outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
586 
587  // Set the solver string for the output test
588  std::string solverDesc = " TFQMR ";
589  outputTest_->setSolverDesc( solverDesc );
590 
591 
592  // The status test is now set.
593  isSTSet_ = true;
594 
595  return false;
596 }
597 
598 
599 template<class ScalarType, class MV, class OP>
602 {
604 
605  // Set all the valid parameters and their default values.
606  if(is_null(validPL)) {
607  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
608 
609  // The static_cast is to resolve an issue with older clang versions which
610  // would cause the constexpr to link fail. With c++17 the problem is resolved.
611  pl->set("Convergence Tolerance", static_cast<MagnitudeType>(DefaultSolverParameters::convTol),
612  "The relative residual tolerance that needs to be achieved by the\n"
613  "iterative solver in order for the linear system to be declared converged.");
614  pl->set("Implicit Tolerance Scale Factor", static_cast<MagnitudeType>(DefaultSolverParameters::impTolScale),
615  "The scale factor used by the implicit residual test when explicit residual\n"
616  "testing is used. May enable faster convergence when TFQMR bound is too loose.");
617  pl->set("Maximum Iterations", static_cast<int>(maxIters_default_),
618  "The maximum number of block iterations allowed for each\n"
619  "set of RHS solved.");
620  pl->set("Verbosity", static_cast<int>(verbosity_default_),
621  "What type(s) of solver information should be outputted\n"
622  "to the output stream.");
623  pl->set("Output Style", static_cast<int>(outputStyle_default_),
624  "What style is used for the solver information outputted\n"
625  "to the output stream.");
626  pl->set("Output Frequency", static_cast<int>(outputFreq_default_),
627  "How often convergence information should be outputted\n"
628  "to the output stream.");
629  pl->set("Output Stream", Teuchos::rcp(outputStream_default_,false),
630  "A reference-counted pointer to the output stream where all\n"
631  "solver output is sent.");
632  pl->set("Explicit Residual Test", static_cast<bool>(expResTest_default_),
633  "Whether the explicitly computed residual should be used in the convergence test.");
634  pl->set("Implicit Residual Scaling", static_cast<const char *>(impResScale_default_),
635  "The type of scaling used in the implicit residual convergence test.");
636  pl->set("Explicit Residual Scaling", static_cast<const char *>(expResScale_default_),
637  "The type of scaling used in the explicit residual convergence test.");
638  pl->set("Timer Label", static_cast<const char *>(label_default_),
639  "The string to use as a prefix for the timer labels.");
640  validPL = pl;
641  }
642  return validPL;
643 }
644 
645 
646 // solve()
647 template<class ScalarType, class MV, class OP>
649 
650  // Set the current parameters if they were not set before.
651  // NOTE: This may occur if the user generated the solver manager with the default constructor and
652  // then didn't set any parameters using setParameters().
653  if (!isSet_) {
654  setParameters(Teuchos::parameterList(*getValidParameters()));
655  }
656 
658  "Belos::TFQMRSolMgr::solve(): Linear problem is not a valid object.");
659 
661  "Belos::TFQMRSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
662 
663  if (!isSTSet_) {
665  "Belos::TFQMRSolMgr::solve(): Linear problem and requested status tests are incompatible.");
666  }
667 
668  // Create indices for the linear systems to be solved.
669  int startPtr = 0;
670  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
671  int numCurrRHS = blockSize_;
672 
673  std::vector<int> currIdx, currIdx2;
674 
675  // The index set is generated that informs the linear problem that some linear systems are augmented.
676  currIdx.resize( blockSize_ );
677  currIdx2.resize( blockSize_ );
678  for (int i=0; i<numCurrRHS; ++i)
679  { currIdx[i] = startPtr+i; currIdx2[i]=i; }
680 
681  // Inform the linear problem of the current linear system to solve.
682  problem_->setLSIndex( currIdx );
683 
685  // Parameter list
687  plist.set("Block Size",blockSize_);
688 
689  // Reset the status test.
690  outputTest_->reset();
691 
692  // Assume convergence is achieved, then let any failed convergence set this to false.
693  bool isConverged = true;
694 
696  // TFQMR solver
697 
699  Teuchos::rcp( new TFQMRIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,plist) );
700 
701  // Enter solve() iterations
702  {
703 #ifdef BELOS_TEUCHOS_TIME_MONITOR
704  Teuchos::TimeMonitor slvtimer(*timerSolve_);
705 #endif
706 
707  while ( numRHS2Solve > 0 ) {
708  //
709  // Reset the active / converged vectors from this block
710  std::vector<int> convRHSIdx;
711  std::vector<int> currRHSIdx( currIdx );
712  currRHSIdx.resize(numCurrRHS);
713 
714  // Reset the number of iterations.
715  tfqmr_iter->resetNumIters();
716 
717  // Reset the number of calls that the status test output knows about.
718  outputTest_->resetNumCalls();
719 
720  // Get the current residual for this block of linear systems.
721  Teuchos::RCP<MV> R_0 = MVT::CloneViewNonConst( *(Teuchos::rcp_const_cast<MV>(problem_->getInitPrecResVec())), currIdx );
722 
723  // Set the new state and initialize the solver.
725  newstate.R = R_0;
726  tfqmr_iter->initializeTFQMR(newstate);
727 
728  while(1) {
729 
730  // tell tfqmr_iter to iterate
731  try {
732  tfqmr_iter->iterate();
733 
735  //
736  // check convergence first
737  //
739  if ( convTest_->getStatus() == Passed ) {
740  // We have convergence of the linear system.
741  break; // break from while(1){tfqmr_iter->iterate()}
742  }
744  //
745  // check for maximum iterations
746  //
748  else if ( maxIterTest_->getStatus() == Passed ) {
749  // we don't have convergence
750  isConverged = false;
751  break; // break from while(1){tfqmr_iter->iterate()}
752  }
753 
755  //
756  // we returned from iterate(), but none of our status tests Passed.
757  // something is wrong, and it is probably our fault.
758  //
760 
761  else {
762  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
763  "Belos::TFQMRSolMgr::solve(): Invalid return from TFQMRIter::iterate().");
764  }
765  }
766  catch (const std::exception &e) {
767  printer_->stream(Errors) << "Error! Caught std::exception in TFQMRIter::iterate() at iteration "
768  << tfqmr_iter->getNumIters() << std::endl
769  << e.what() << std::endl;
770  throw;
771  }
772  }
773 
774  // Update the current solution with the update computed by the iteration object.
775  problem_->updateSolution( tfqmr_iter->getCurrentUpdate(), true );
776 
777  // Inform the linear problem that we are finished with this block linear system.
778  problem_->setCurrLS();
779 
780  // Update indices for the linear systems to be solved.
781  startPtr += numCurrRHS;
782  numRHS2Solve -= numCurrRHS;
783  if ( numRHS2Solve > 0 ) {
784  numCurrRHS = blockSize_;
785 
786  currIdx.resize( blockSize_ );
787  currIdx2.resize( blockSize_ );
788  for (int i=0; i<numCurrRHS; ++i)
789  { currIdx[i] = startPtr+i; currIdx2[i] = i; }
790  // Set the next indices.
791  problem_->setLSIndex( currIdx );
792 
793  // Set the new blocksize for the solver.
794  tfqmr_iter->setBlockSize( blockSize_ );
795  }
796  else {
797  currIdx.resize( numRHS2Solve );
798  }
799 
800  }// while ( numRHS2Solve > 0 )
801 
802  }
803 
804  // print final summary
805  sTest_->print( printer_->stream(FinalSummary) );
806 
807  // print timing information
808 #ifdef BELOS_TEUCHOS_TIME_MONITOR
809  // Calling summarize() can be expensive, so don't call unless the
810  // user wants to print out timing details. summarize() will do all
811  // the work even if it's passed a "black hole" output stream.
812  if (verbosity_ & TimingDetails)
813  Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
814 #endif
815 
816  // get iteration information for this solve
817  numIters_ = maxIterTest_->getNumIters();
818 
819  // Save the convergence test value ("achieved tolerance") for this
820  // solve. For this solver, convTest_ may either be a single
821  // (implicit) residual norm test, or a combination of two residual
822  // norm tests. In the latter case, the master convergence test
823  // convTest_ is a SEQ combo of the implicit resp. explicit tests.
824  // If the implicit test never passes, then the explicit test won't
825  // ever be executed. This manifests as
826  // expConvTest_->getTestValue()->size() < 1. We deal with this case
827  // by using the values returned by impConvTest_->getTestValue().
828  {
829  // We'll fetch the vector of residual norms one way or the other.
830  const std::vector<MagnitudeType>* pTestValues = NULL;
831  if (expResTest_) {
832  pTestValues = expConvTest_->getTestValue();
833  if (pTestValues == NULL || pTestValues->size() < 1) {
834  pTestValues = impConvTest_->getTestValue();
835  }
836  }
837  else {
838  // Only the implicit residual norm test is being used.
839  pTestValues = impConvTest_->getTestValue();
840  }
841  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
842  "Belos::TFQMRSolMgr::solve(): The implicit convergence test's "
843  "getTestValue() method returned NULL. Please report this bug to the "
844  "Belos developers.");
845  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
846  "Belos::TMQMRSolMgr::solve(): The implicit convergence test's "
847  "getTestValue() method returned a vector of length zero. Please report "
848  "this bug to the Belos developers.");
849 
850  // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
851  // achieved tolerances for all vectors in the current solve(), or
852  // just for the vectors from the last deflation?
853  achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
854  }
855 
856  if (!isConverged) {
857  return Unconverged; // return from TFQMRSolMgr::solve()
858  }
859  return Converged; // return from TFQMRSolMgr::solve()
860 }
861 
862 // This method requires the solver manager to return a std::string that describes itself.
863 template<class ScalarType, class MV, class OP>
865 {
866  std::ostringstream oss;
867  oss << "Belos::TFQMRSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
868  oss << "{}";
869  return oss.str();
870 }
871 
872 } // end Belos namespace
873 
874 #endif /* BELOS_TFQMR_SOLMGR_HPP */
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Definition: BelosTypes.cpp:106
Collection of types and exceptions used within the Belos solvers.
This class implements the preconditioned transpose-free QMR algorithm for solving non-Hermitian linea...
int getNumIters() const override
Get the iteration count for the most recent call to solve().
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
Class which manages the output and verbosity of the Belos solvers.
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
Set the parameters the solver manager should use to solve the linear problem.
bool is_null(const boost::shared_ptr< T > &p)
static constexpr int outputStyle_default_
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > outputTest_
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > convTest_
static constexpr bool expResTest_default_
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
bool isLOADetected() const override
Whether loss of accuracy was detected during the last solve() invocation.
T & get(ParameterList &l, const std::string &name)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
static RCP< Time > getNewCounter(const std::string &name)
bool is_null(const std::shared_ptr< T > &p)
static constexpr int outputFreq_default_
MultiVecTraits< ScalarType, MV > MVT
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
A factory class for generating StatusTestOutput objects.
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Return a reference to the linear problem being solved by this solver manager.
Teuchos::RCP< Teuchos::Time > timerSolve_
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
static constexpr const char * label_default_
An implementation of StatusTestResNorm using a family of residual norms.
static constexpr int maxIters_default_
Belos concrete class for generating iterations with the preconditioned tranpose-free QMR (TFQMR) meth...
Teuchos::RCP< const MV > R
The current residual basis.
static const double convTol
Default convergence tolerance.
Definition: BelosTypes.hpp:293
Belos::StatusTest class for specifying a maximum number of iterations.
Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > expConvTest_
Teuchos::RCP< OutputManager< ScalarType > > printer_
std::string description() const override
Method to return description of the TFQMR solver manager.
static std::string name()
Teuchos::ScalarTraits< ScalarType > SCT
A factory class for generating StatusTestOutput objects.
Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > impConvTest_
TFQMRSolMgrOrthoFailure(const std::string &what_arg)
ReturnType solve() override
This method performs possibly repeated calls to the underlying linear solver&#39;s iterate() routine unti...
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
static constexpr int verbosity_default_
bool isParameter(const std::string &name) const
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
A Belos::StatusTest class for specifying a maximum number of iterations.
ResetType
How to reset the solver.
Definition: BelosTypes.hpp:206
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Pure virtual base class which describes the basic interface for a solver manager. ...
static const double impTolScale
&quot;Implicit Tolerance Scale Factor&quot;
Definition: BelosTypes.hpp:305
static void summarize(Ptr< const Comm< int > > comm, std::ostream &out=std::cout, const bool alwaysWriteLocal=false, const bool writeGlobalStats=true, const bool writeZeroTimers=true, const ECounterSetOp setOp=Intersection, const std::string &filter="", const bool ignoreZeroTimers=false)
Teuchos::RCP< StatusTestMaxIters< ScalarType, MV, OP > > maxIterTest_
Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > problem_
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
A linear system to solve, and its associated information.
Structure to contain pointers to TFQMRIter state variables.
Class which describes the linear problem to be solved by the iterative solver.
TFQMRSolMgr()
Empty constructor for TFQMRSolMgr. This constructor takes no arguments and sets the default values fo...
static constexpr std::ostream * outputStream_default_
Teuchos::RCP< std::ostream > outputStream_
virtual ~TFQMRSolMgr()
Destructor.
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:155
void reset(const ResetType type) override
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
The Belos::TFQMRSolMgr provides a powerful and fully-featured solver manager over the TFQMR linear so...
MagnitudeType impTolScale_
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
static constexpr const char * expResScale_default_
Teuchos::RCP< Teuchos::ParameterList > params_
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > create(const Teuchos::RCP< OutputManager< ScalarType > > &printer, Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test, int mod, int printStates)
Create the StatusTestOutput object specified by the outputStyle.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Get a parameter list containing the valid parameters for this object.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
bool isType(const std::string &name) const
A class for extending the status testing capabilities of Belos via logical combinations.
MagnitudeType achievedTol_
TFQMRSolMgrLinearProblemFailure(const std::string &what_arg)
static constexpr const char * impResScale_default_
Class which defines basic traits for the operator type.
TFQMRSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate orthonorma...
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:60
Default parameters common to most Belos solvers.
Definition: BelosTypes.hpp:283
TFQMRSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i.e.
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > sTest_
Belos header file which uses auto-configuration information to include necessary C++ headers...
OperatorTraits< ScalarType, MV, OP > OPT
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
Teuchos::ScalarTraits< MagnitudeType > MT