Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosPseudoBlockTFQMRSolMgr.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_PSEUDO_BLOCK_TFQMR_SOLMGR_HPP
11 #define BELOS_PSEUDO_BLOCK_TFQMR_SOLMGR_HPP
12 
17 #include "BelosConfigDefs.hpp"
18 #include "BelosTypes.hpp"
19 
20 #include "BelosLinearProblem.hpp"
21 #include "BelosSolverManager.hpp"
22 
26 #include "BelosStatusTestCombo.hpp"
28 #include "BelosOutputManager.hpp"
29 #ifdef BELOS_TEUCHOS_TIME_MONITOR
30 #include "Teuchos_TimeMonitor.hpp"
31 #endif
32 
49 namespace Belos {
50 
52 
53 
61  PseudoBlockTFQMRSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
62  {}};
63 
64  template<class ScalarType, class MV, class OP>
65  class PseudoBlockTFQMRSolMgr : public SolverManager<ScalarType,MV,OP> {
66 
67  private:
71  typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
73 
74  public:
75 
77 
78 
85 
104 
107 
111  }
113 
115 
116 
117  const LinearProblem<ScalarType,MV,OP>& getProblem() const override {
118  return *problem_;
119  }
120 
124 
128 
135  return Teuchos::tuple(timerSolve_);
136  }
137 
143  MagnitudeType achievedTol() const override {
144  return achievedTol_;
145  }
146 
148  int getNumIters() const override {
149  return numIters_;
150  }
151 
159  bool isLOADetected() const override { return false; }
161 
163 
164 
166  void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) override { problem_ = problem; }
167 
169  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params ) override;
170 
172 
174 
179  void reset( const ResetType type ) override { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
181 
183 
184 
202  ReturnType solve() override;
203 
205 
207 
209  std::string description() const override;
210 
212  private:
213 
214  // Method for checking current status test against defined linear problem.
215  bool checkStatusTest();
216 
217  // Linear problem.
219 
220  // Output manager.
222  Teuchos::RCP<std::ostream> outputStream_;
223 
224  // Status test.
228  Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > expConvTest_, impConvTest_;
230 
231  // Current parameter list.
233 
234  // Default solver values.
235  static constexpr int maxIters_default_ = 1000;
236  static constexpr bool expResTest_default_ = false;
237  static constexpr int verbosity_default_ = Belos::Errors;
238  static constexpr int outputStyle_default_ = Belos::General;
239  static constexpr int outputFreq_default_ = -1;
240  static constexpr int defQuorum_default_ = 1;
241  static constexpr const char * impResScale_default_ = "Norm of Preconditioned Initial Residual";
242  static constexpr const char * expResScale_default_ = "Norm of Initial Residual";
243  static constexpr const char * label_default_ = "Belos";
244 
245  // Current solver values.
246  MagnitudeType convtol_, impTolScale_, achievedTol_;
247  int maxIters_, numIters_;
248  int verbosity_, outputStyle_, outputFreq_, defQuorum_;
249  bool expResTest_;
250  std::string impResScale_, expResScale_;
251 
252  // Timers.
253  std::string label_;
254  Teuchos::RCP<Teuchos::Time> timerSolve_;
255 
256  // Internal state variables.
257  bool isSet_, isSTSet_;
258  };
259 
260 
261 // Empty Constructor
262 template<class ScalarType, class MV, class OP>
264  outputStream_(Teuchos::rcpFromRef(std::cout)),
265  convtol_(DefaultSolverParameters::convTol),
266  impTolScale_(DefaultSolverParameters::impTolScale),
267  achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
268  maxIters_(maxIters_default_),
269  numIters_(0),
270  verbosity_(verbosity_default_),
271  outputStyle_(outputStyle_default_),
272  outputFreq_(outputFreq_default_),
273  defQuorum_(defQuorum_default_),
274  expResTest_(expResTest_default_),
275  impResScale_(impResScale_default_),
276  expResScale_(expResScale_default_),
277  label_(label_default_),
278  isSet_(false),
279  isSTSet_(false)
280 {}
281 
282 
283 // Basic Constructor
284 template<class ScalarType, class MV, class OP>
288  problem_(problem),
289  outputStream_(Teuchos::rcpFromRef(std::cout)),
290  convtol_(DefaultSolverParameters::convTol),
291  impTolScale_(DefaultSolverParameters::impTolScale),
292  achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
293  maxIters_(maxIters_default_),
294  numIters_(0),
295  verbosity_(verbosity_default_),
296  outputStyle_(outputStyle_default_),
297  outputFreq_(outputFreq_default_),
298  defQuorum_(defQuorum_default_),
299  expResTest_(expResTest_default_),
300  impResScale_(impResScale_default_),
301  expResScale_(expResScale_default_),
302  label_(label_default_),
303  isSet_(false),
304  isSTSet_(false)
305 {
306  TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null, std::invalid_argument, "Problem not given to solver manager.");
307 
308  // If the parameter list pointer is null, then set the current parameters to the default parameter list.
309  if ( !is_null(pl) ) {
310  setParameters( pl );
311  }
312 }
313 
314 template<class ScalarType, class MV, class OP>
316 {
317  // Create the internal parameter list if ones doesn't already exist.
318  if (params_ == Teuchos::null) {
319  params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
320  }
321  else {
322  params->validateParameters(*getValidParameters());
323  }
324 
325  // Check for maximum number of iterations
326  if (params->isParameter("Maximum Iterations")) {
327  maxIters_ = params->get("Maximum Iterations",maxIters_default_);
328 
329  // Update parameter in our list and in status test.
330  params_->set("Maximum Iterations", maxIters_);
331  if (maxIterTest_!=Teuchos::null)
332  maxIterTest_->setMaxIters( maxIters_ );
333  }
334 
335  // Check to see if the timer label changed.
336  if (params->isParameter("Timer Label")) {
337  std::string tempLabel = params->get("Timer Label", label_default_);
338 
339  // Update parameter in our list and solver timer
340  if (tempLabel != label_) {
341  label_ = tempLabel;
342  params_->set("Timer Label", label_);
343  std::string solveLabel = label_ + ": PseudoBlockTFQMRSolMgr total solve time";
344 #ifdef BELOS_TEUCHOS_TIME_MONITOR
345  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
346 #endif
347  }
348  }
349 
350  // Check for a change in verbosity level
351  if (params->isParameter("Verbosity")) {
352  if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
353  verbosity_ = params->get("Verbosity", verbosity_default_);
354  } else {
355  verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
356  }
357 
358  // Update parameter in our list.
359  params_->set("Verbosity", verbosity_);
360  if (printer_ != Teuchos::null)
361  printer_->setVerbosity(verbosity_);
362  }
363 
364  // Check for a change in output style
365  if (params->isParameter("Output Style")) {
366  if (Teuchos::isParameterType<int>(*params,"Output Style")) {
367  outputStyle_ = params->get("Output Style", outputStyle_default_);
368  } else {
369  outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
370  }
371 
372  // Reconstruct the convergence test if the explicit residual test is not being used.
373  params_->set("Output Style", outputStyle_);
374  isSTSet_ = false;
375  }
376 
377  // output stream
378  if (params->isParameter("Output Stream")) {
379  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
380 
381  // Update parameter in our list.
382  params_->set("Output Stream", outputStream_);
383  if (printer_ != Teuchos::null)
384  printer_->setOStream( outputStream_ );
385  }
386 
387  // frequency level
388  if (verbosity_ & Belos::StatusTestDetails) {
389  if (params->isParameter("Output Frequency")) {
390  outputFreq_ = params->get("Output Frequency", outputFreq_default_);
391  }
392 
393  // Update parameter in out list and output status test.
394  params_->set("Output Frequency", outputFreq_);
395  if (outputTest_ != Teuchos::null)
396  outputTest_->setOutputFrequency( outputFreq_ );
397  }
398 
399  // Create output manager if we need to.
400  if (printer_ == Teuchos::null) {
401  printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
402  }
403 
404  // Check for convergence tolerance
405  if (params->isParameter("Convergence Tolerance")) {
406  if (params->isType<MagnitudeType> ("Convergence Tolerance")) {
407  convtol_ = params->get ("Convergence Tolerance",
408  static_cast<MagnitudeType> (DefaultSolverParameters::convTol));
409  }
410  else {
411  convtol_ = params->get ("Convergence Tolerance", DefaultSolverParameters::convTol);
412  }
413 
414  // Update parameter in our list and residual tests.
415  params_->set("Convergence Tolerance", convtol_);
416  isSTSet_ = false;
417  }
418 
419  if (params->isParameter("Implicit Tolerance Scale Factor")) {
420  if (params->isType<MagnitudeType> ("Implicit Tolerance Scale Factor")) {
421  impTolScale_ = params->get ("Implicit Tolerance Scale Factor",
422  static_cast<MagnitudeType> (DefaultSolverParameters::impTolScale));
423 
424  }
425  else {
426  impTolScale_ = params->get ("Implicit Tolerance Scale Factor",
428  }
429 
430  // Update parameter in our list.
431  params_->set("Implicit Tolerance Scale Factor", impTolScale_);
432  isSTSet_ = false;
433  }
434 
435  if (params->isParameter("Implicit Residual Scaling")) {
436  std::string tempImpResScale = Teuchos::getParameter<std::string>( *params, "Implicit Residual Scaling" );
437 
438  // Only update the scaling if it's different.
439  if (impResScale_ != tempImpResScale) {
440  impResScale_ = tempImpResScale;
441 
442  // Update parameter in our list.
443  params_->set("Implicit Residual Scaling", impResScale_);
444  isSTSet_ = false;
445  }
446  }
447 
448  if (params->isParameter("Explicit Residual Scaling")) {
449  std::string tempExpResScale = Teuchos::getParameter<std::string>( *params, "Explicit Residual Scaling" );
450 
451  // Only update the scaling if it's different.
452  if (expResScale_ != tempExpResScale) {
453  expResScale_ = tempExpResScale;
454 
455  // Update parameter in our list.
456  params_->set("Explicit Residual Scaling", expResScale_);
457  isSTSet_ = false;
458  }
459  }
460 
461  if (params->isParameter("Explicit Residual Test")) {
462  expResTest_ = Teuchos::getParameter<bool>( *params,"Explicit Residual Test" );
463 
464  // Reconstruct the convergence test if the explicit residual test is not being used.
465  params_->set("Explicit Residual Test", expResTest_);
466  if (expConvTest_ == Teuchos::null) {
467  isSTSet_ = false;
468  }
469  }
470 
471  // Get the deflation quorum, or number of converged systems before deflation is allowed
472  if (params->isParameter("Deflation Quorum")) {
473  defQuorum_ = params->get("Deflation Quorum", defQuorum_);
474  params_->set ("Deflation Quorum", defQuorum_);
475  if (! impConvTest_.is_null ()) {
476  impConvTest_->setQuorum (defQuorum_);
477  }
478  if (! expConvTest_.is_null ()) {
479  expConvTest_->setQuorum (defQuorum_);
480  }
481  }
482 
483  // Create the timer if we need to.
484  if (timerSolve_ == Teuchos::null) {
485  std::string solveLabel = label_ + ": PseudoBlockTFQMRSolMgr total solve time";
486 #ifdef BELOS_TEUCHOS_TIME_MONITOR
487  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
488 #endif
489  }
490 
491  // Inform the solver manager that the current parameters were set.
492  isSet_ = true;
493 }
494 
495 
496 // Check the status test versus the defined linear problem
497 template<class ScalarType, class MV, class OP>
499 
500  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
501  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestGenResNorm_t;
502 
503  // Basic test checks maximum iterations and native residual.
504  maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
505 
506  if (expResTest_) {
507 
508  // Implicit residual test, using the native residual to determine if convergence was achieved.
509  Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
510  Teuchos::rcp( new StatusTestGenResNorm_t( impTolScale_*convtol_, defQuorum_ ) );
511  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
512  impConvTest_ = tmpImpConvTest;
513 
514  // Explicit residual test once the native residual is below the tolerance
515  Teuchos::RCP<StatusTestGenResNorm_t> tmpExpConvTest =
516  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_, defQuorum_ ) );
517  tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit, Belos::TwoNorm );
518  tmpExpConvTest->defineScaleForm( convertStringToScaleType(expResScale_), Belos::TwoNorm );
519  expConvTest_ = tmpExpConvTest;
520 
521  // The convergence test is a combination of the "cheap" implicit test and explicit test.
522  convTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
523  }
524  else {
525 
526  // Implicit residual test, using the native residual to determine if convergence was achieved.
527  Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
528  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_, defQuorum_ ) );
529  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
530  impConvTest_ = tmpImpConvTest;
531 
532  // Set the explicit and total convergence test to this implicit test that checks for accuracy loss.
533  expConvTest_ = impConvTest_;
534  convTest_ = impConvTest_;
535  }
536  sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
537 
538  // Create the status test output class.
539  // This class manages and formats the output from the status test.
540  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
541  outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
542 
543  // Set the solver string for the output test
544  std::string solverDesc = " Pseudo Block TFQMR ";
545  outputTest_->setSolverDesc( solverDesc );
546 
547 
548  // The status test is now set.
549  isSTSet_ = true;
550 
551  return false;
552 }
553 
554 
555 template<class ScalarType, class MV, class OP>
558 {
560 
561  // Set all the valid parameters and their default values.
562  if(is_null(validPL)) {
563  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
564 
565  // The static_cast is to resolve an issue with older clang versions which
566  // would cause the constexpr to link fail. With c++17 the problem is resolved.
567  pl->set("Convergence Tolerance", static_cast<MagnitudeType>(DefaultSolverParameters::convTol),
568  "The relative residual tolerance that needs to be achieved by the\n"
569  "iterative solver in order for the linear system to be declared converged.");
570  pl->set("Implicit Tolerance Scale Factor", static_cast<MagnitudeType>(DefaultSolverParameters::impTolScale),
571  "The scale factor used by the implicit residual test when explicit residual\n"
572  "testing is used. May enable faster convergence when TFQMR bound is too loose.");
573  pl->set("Maximum Iterations", static_cast<int>(maxIters_default_),
574  "The maximum number of block iterations allowed for each\n"
575  "set of RHS solved.");
576  pl->set("Verbosity", static_cast<int>(verbosity_default_),
577  "What type(s) of solver information should be outputted\n"
578  "to the output stream.");
579  pl->set("Output Style", static_cast<int>(outputStyle_default_),
580  "What style is used for the solver information outputted\n"
581  "to the output stream.");
582  pl->set("Output Frequency", static_cast<int>(outputFreq_default_),
583  "How often convergence information should be outputted\n"
584  "to the output stream.");
585  pl->set("Deflation Quorum", static_cast<int>(defQuorum_default_),
586  "The number of linear systems that need to converge before they are deflated.");
587  pl->set("Output Stream", Teuchos::rcpFromRef(std::cout),
588  "A reference-counted pointer to the output stream where all\n"
589  "solver output is sent.");
590  pl->set("Explicit Residual Test", static_cast<bool>(expResTest_default_),
591  "Whether the explicitly computed residual should be used in the convergence test.");
592  pl->set("Implicit Residual Scaling", static_cast<const char *>(impResScale_default_),
593  "The type of scaling used in the implicit residual convergence test.");
594  pl->set("Explicit Residual Scaling", static_cast<const char *>(expResScale_default_),
595  "The type of scaling used in the explicit residual convergence test.");
596  pl->set("Timer Label", static_cast<const char *>(label_default_),
597  "The string to use as a prefix for the timer labels.");
598  validPL = pl;
599  }
600  return validPL;
601 }
602 
603 
604 // solve()
605 template<class ScalarType, class MV, class OP>
607 
608  // Set the current parameters if they were not set before.
609  // NOTE: This may occur if the user generated the solver manager with the default constructor and
610  // then didn't set any parameters using setParameters().
611  if (!isSet_) {
612  setParameters(Teuchos::parameterList(*getValidParameters()));
613  }
614 
616  "Belos::PseudoBlockTFQMRSolMgr::solve(): Linear problem is not a valid object.");
617 
619  "Belos::PseudoBlockTFQMRSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
620 
621  if (!isSTSet_) {
623  "Belos::PseudoBlockTFQMRSolMgr::solve(): Linear problem and requested status tests are incompatible.");
624  }
625 
626  // Create indices for the linear systems to be solved.
627  int startPtr = 0;
628  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
629  int numCurrRHS = numRHS2Solve;
630 
631  std::vector<int> currIdx( numRHS2Solve );
632  for (int i=0; i<numRHS2Solve; ++i) {
633  currIdx[i] = startPtr+i;
634  }
635 
636  // Inform the linear problem of the current linear system to solve.
637  problem_->setLSIndex( currIdx );
638 
640  // Parameter list
642 
643  // Reset the status test.
644  outputTest_->reset();
645 
646  // Assume convergence is achieved, then let any failed convergence set this to false.
647  bool isConverged = true;
648 
650  // TFQMR solver
651 
653  Teuchos::rcp( new PseudoBlockTFQMRIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,plist) );
654 
655  // Enter solve() iterations
656  {
657 #ifdef BELOS_TEUCHOS_TIME_MONITOR
658  Teuchos::TimeMonitor slvtimer(*timerSolve_);
659 #endif
660 
661  while ( numRHS2Solve > 0 ) {
662  //
663  // Reset the active / converged vectors from this block
664  std::vector<int> convRHSIdx;
665  std::vector<int> currRHSIdx( currIdx );
666  currRHSIdx.resize(numCurrRHS);
667 
668  // Reset the number of iterations.
669  block_tfqmr_iter->resetNumIters();
670 
671  // Reset the number of calls that the status test output knows about.
672  outputTest_->resetNumCalls();
673 
674  // Get the current residual for this block of linear systems.
675  Teuchos::RCP<MV> R_0 = MVT::CloneViewNonConst( *(Teuchos::rcp_const_cast<MV>(problem_->getInitPrecResVec())), currIdx );
676 
677  // Set the new state and initialize the solver.
679  newstate.Rtilde = R_0;
680  block_tfqmr_iter->initializeTFQMR(newstate);
681 
682  while(1) {
683 
684  // tell block_tfqmr_iter to iterate
685  try {
686  block_tfqmr_iter->iterate();
687 
689  //
690  // check convergence first
691  //
693  if ( convTest_->getStatus() == Passed ) {
694 
695  // Figure out which linear systems converged.
696  std::vector<int> convIdx = expConvTest_->convIndices();
697 
698  // If the number of converged linear systems is equal to the
699  // number of current linear systems, then we are done with this block.
700  if (convIdx.size() == currRHSIdx.size())
701  break; // break from while(1){block_tfqmr_iter->iterate()}
702 
703  // Update the current solution with the update computed by the iteration object.
704  problem_->updateSolution( block_tfqmr_iter->getCurrentUpdate(), true );
705 
706  // Inform the linear problem that we are finished with this current linear system.
707  problem_->setCurrLS();
708 
709  // Reset currRHSIdx to have the right-hand sides that are left to converge for this block.
710  int have = 0;
711  std::vector<int> unconvIdx( currRHSIdx.size() );
712  for (unsigned int i=0; i<currRHSIdx.size(); ++i) {
713  bool found = false;
714  for (unsigned int j=0; j<convIdx.size(); ++j) {
715  if (currRHSIdx[i] == convIdx[j]) {
716  found = true;
717  break;
718  }
719  }
720  if (!found) {
721  unconvIdx[have] = i;
722  currRHSIdx[have++] = currRHSIdx[i];
723  }
724  }
725  unconvIdx.resize(have);
726  currRHSIdx.resize(have);
727 
728  // Set the remaining indices after deflation.
729  problem_->setLSIndex( currRHSIdx );
730 
731  // Get the current residual vector.
732  // Set the new state and initialize the solver.
733  PseudoBlockTFQMRIterState<ScalarType,MV> currentState = block_tfqmr_iter->getState();
734 
735  // Set the new state and initialize the solver.
737 
738  // Copy over the vectors.
739  defstate.Rtilde = MVT::CloneView( *currentState.Rtilde, unconvIdx);
740  defstate.U = MVT::CloneView( *currentState.U, unconvIdx );
741  defstate.AU = MVT::CloneView( *currentState.AU, unconvIdx );
742  defstate.V = MVT::CloneView( *currentState.V, unconvIdx );
743  defstate.W = MVT::CloneView( *currentState.W, unconvIdx );
744  defstate.D = MVT::CloneView( *currentState.D, unconvIdx );
745 
746  // Copy over the scalars.
747  for (std::vector<int>::iterator uIter = unconvIdx.begin(); uIter != unconvIdx.end(); uIter++)
748  {
749  defstate.alpha.push_back( currentState.alpha[ *uIter ] );
750  defstate.eta.push_back( currentState.eta[ *uIter ] );
751  defstate.rho.push_back( currentState.rho[ *uIter ] );
752  defstate.tau.push_back( currentState.tau[ *uIter ] );
753  defstate.theta.push_back( currentState.theta[ *uIter ] );
754  }
755 
756  block_tfqmr_iter->initializeTFQMR(defstate);
757  }
759  //
760  // check for maximum iterations
761  //
763  else if ( maxIterTest_->getStatus() == Passed ) {
764  // we don't have convergence
765  isConverged = false;
766  break; // break from while(1){block_tfqmr_iter->iterate()}
767  }
768 
770  //
771  // we returned from iterate(), but none of our status tests Passed.
772  // something is wrong, and it is probably our fault.
773  //
775 
776  else {
777  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
778  "Belos::PseudoBlockTFQMRSolMgr::solve(): Invalid return from PseudoBlockTFQMRIter::iterate().");
779  }
780  }
781  catch (const StatusTestNaNError& e) {
782  // A NaN was detected in the solver. Set the solution to zero and return unconverged.
783  achievedTol_ = MT::one();
784  Teuchos::RCP<MV> X = problem_->getLHS();
785  MVT::MvInit( *X, SCT::zero() );
786  printer_->stream(Warnings) << "Belos::PseudoBlockTFQMRSolMgr::solve(): Warning! NaN has been detected!"
787  << std::endl;
788  return Unconverged;
789  }
790  catch (const std::exception &e) {
791  printer_->stream(Errors) << "Error! Caught std::exception in PseudoBlockTFQMRIter::iterate() at iteration "
792  << block_tfqmr_iter->getNumIters() << std::endl
793  << e.what() << std::endl;
794  throw;
795  }
796  }
797 
798  // Update the current solution with the update computed by the iteration object.
799  problem_->updateSolution( block_tfqmr_iter->getCurrentUpdate(), true );
800 
801  // Inform the linear problem that we are finished with this block linear system.
802  problem_->setCurrLS();
803 
804  // Update indices for the linear systems to be solved.
805  startPtr += numCurrRHS;
806  numRHS2Solve -= numCurrRHS;
807  if ( numRHS2Solve > 0 ) {
808  numCurrRHS = numRHS2Solve;
809  currIdx.resize( numCurrRHS );
810  for (int i=0; i<numCurrRHS; ++i)
811  { currIdx[i] = startPtr+i; }
812 
813  // Adapt the status test quorum if we need to.
814  if (defQuorum_ > numCurrRHS) {
815  if (impConvTest_ != Teuchos::null)
816  impConvTest_->setQuorum( numCurrRHS );
817  if (expConvTest_ != Teuchos::null)
818  expConvTest_->setQuorum( numCurrRHS );
819  }
820 
821  // Set the next indices.
822  problem_->setLSIndex( currIdx );
823  }
824  else {
825  currIdx.resize( numRHS2Solve );
826  }
827 
828  }// while ( numRHS2Solve > 0 )
829 
830  }
831 
832  // print final summary
833  sTest_->print( printer_->stream(FinalSummary) );
834 
835  // print timing information
836 #ifdef BELOS_TEUCHOS_TIME_MONITOR
837  // Calling summarize() can be expensive, so don't call unless the
838  // user wants to print out timing details. summarize() will do all
839  // the work even if it's passed a "black hole" output stream.
840  if (verbosity_ & TimingDetails)
841  Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
842 #endif
843 
844  // get iteration information for this solve
845  numIters_ = maxIterTest_->getNumIters();
846 
847  // Save the convergence test value ("achieved tolerance") for this
848  // solve. For this solver, convTest_ may either be a single
849  // (implicit) residual norm test, or a combination of two residual
850  // norm tests. In the latter case, the master convergence test
851  // convTest_ is a SEQ combo of the implicit resp. explicit tests.
852  // If the implicit test never passes, then the explicit test won't
853  // ever be executed. This manifests as
854  // expConvTest_->getTestValue()->size() < 1. We deal with this case
855  // by using the values returned by impConvTest_->getTestValue().
856  {
857  // We'll fetch the vector of residual norms one way or the other.
858  const std::vector<MagnitudeType>* pTestValues = NULL;
859  if (expResTest_) {
860  pTestValues = expConvTest_->getTestValue();
861  if (pTestValues == NULL || pTestValues->size() < 1) {
862  pTestValues = impConvTest_->getTestValue();
863  }
864  }
865  else {
866  // Only the implicit residual norm test is being used.
867  pTestValues = impConvTest_->getTestValue();
868  }
869  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
870  "Belos::PseudoBlockTFQMRSolMgr::solve(): The implicit convergence test's "
871  "getTestValue() method returned NULL. Please report this bug to the "
872  "Belos developers.");
873  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
874  "Belos::TMQMRSolMgr::solve(): The implicit convergence test's "
875  "getTestValue() method returned a vector of length zero. Please report "
876  "this bug to the Belos developers.");
877 
878  // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
879  // achieved tolerances for all vectors in the current solve(), or
880  // just for the vectors from the last deflation?
881  achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
882  }
883 
884  if (!isConverged) {
885  return Unconverged; // return from PseudoBlockTFQMRSolMgr::solve()
886  }
887  return Converged; // return from PseudoBlockTFQMRSolMgr::solve()
888 }
889 
890 // This method requires the solver manager to return a std::string that describes itself.
891 template<class ScalarType, class MV, class OP>
893 {
894  std::ostringstream oss;
895  oss << "Belos::PseudoBlockTFQMRSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
896  oss << "{}";
897  return oss.str();
898 }
899 
900 } // end Belos namespace
901 
902 #endif /* BELOS_PSEUDO_BLOCK_TFQMR_SOLMGR_HPP */
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Definition: BelosTypes.cpp:74
Collection of types and exceptions used within the Belos solvers.
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
void reset(const ResetType type) override
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
Class which manages the output and verbosity of the Belos solvers.
bool is_null(const boost::shared_ptr< T > &p)
This class implements the preconditioned transpose-free QMR algorithm for solving non-Hermitian linea...
T & get(const std::string &name, T def_value)
Structure to contain pointers to PseudoBlockTFQMRIter state variables.
static RCP< Time > getNewCounter(const std::string &name)
bool is_null(const std::shared_ptr< T > &p)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Belos concrete class for generating iterations with the preconditioned tranpose-free QMR (TFQMR) meth...
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
An implementation of StatusTestResNorm using a family of residual norms.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
static const double convTol
Default convergence tolerance.
Definition: BelosTypes.hpp:261
Belos::StatusTest class for specifying a maximum number of iterations.
static std::string name()
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Return a reference to the linear problem being solved by this solver manager.
A factory class for generating StatusTestOutput objects.
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
bool isParameter(const std::string &name) const
A Belos::StatusTest class for specifying a maximum number of iterations.
PseudoBlockTFQMRSolMgr()
Empty constructor for PseudoBlockTFQMRSolMgr. This constructor takes no arguments and sets the defaul...
PseudoBlockTFQMRSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i...
ResetType
How to reset the solver.
Definition: BelosTypes.hpp:174
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. ...
bool isLOADetected() const override
Whether loss of accuracy was detected during the last solve() invocation.
static const double impTolScale
&quot;Implicit Tolerance Scale Factor&quot;
Definition: BelosTypes.hpp:273
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Get a parameter list containing the valid parameters for this object.
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)
std::string description() const override
Method to return description of the pseudo-block TFQMR solver manager.
ReturnType solve() override
This method performs possibly repeated calls to the underlying linear solver&#39;s iterate() routine unti...
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
The Belos::PseudoBlockTFQMRSolMgr provides a powerful and fully-featured solver manager over the pseu...
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:123
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().
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Teuchos::RCP< const MV > W
The current residual basis.
bool isType(const std::string &name) const
A class for extending the status testing capabilities of Belos via logical combinations.
Class which defines basic traits for the operator type.
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:28
Default parameters common to most Belos solvers.
Definition: BelosTypes.hpp:251
Belos header file which uses auto-configuration information to include necessary C++ headers...
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
Set the parameters the solver manager should use to solve the linear problem.

Generated on Fri Nov 22 2024 09:23:06 for Belos by doxygen 1.8.5