Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosBiCGStabSolMgr.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_BICGSTAB_SOLMGR_HPP
43 #define BELOS_BICGSTAB_SOLMGR_HPP
44 
49 #include "BelosConfigDefs.hpp"
50 #include "BelosTypes.hpp"
51 
52 #include "BelosLinearProblem.hpp"
53 #include "BelosSolverManager.hpp"
54 
55 #include "BelosBiCGStabIter.hpp"
58 #include "BelosStatusTestCombo.hpp"
60 #include "BelosOutputManager.hpp"
61 #include "Teuchos_BLAS.hpp"
62 #include "Teuchos_LAPACK.hpp"
63 #ifdef BELOS_TEUCHOS_TIME_MONITOR
64 #include "Teuchos_TimeMonitor.hpp"
65 #endif
66 
83 namespace Belos {
84 
86 
87 
95  BiCGStabSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
96  {}};
97 
104  class BiCGStabSolMgrOrthoFailure : public BelosError {public:
105  BiCGStabSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
106  {}};
107 
108  template<class ScalarType, class MV, class OP>
109  class BiCGStabSolMgr : public SolverManager<ScalarType,MV,OP> {
110 
111  private:
115  typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
117 
118  public:
119 
121 
122 
128  BiCGStabSolMgr();
129 
141 
143  virtual ~BiCGStabSolMgr() {};
144 
148  }
150 
152 
153 
154  const LinearProblem<ScalarType,MV,OP>& getProblem() const override {
155  return *problem_;
156  }
157 
161 
165 
172  return Teuchos::tuple(timerSolve_);
173  }
174 
175 
186  MagnitudeType achievedTol() const override {
187  return achievedTol_;
188  }
189 
191  int getNumIters() const override {
192  return numIters_;
193  }
194 
198  bool isLOADetected() const override { return false; }
199 
201 
203 
204 
206  void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) override { problem_ = problem; }
207 
209  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params ) override;
210 
212 
214 
215 
219  void reset( const ResetType type ) override { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
221 
223 
224 
242  ReturnType solve() override;
243 
245 
248 
250  std::string description() const override;
251 
253  private:
254 
255  // Linear problem.
257 
258  // Output manager.
260  Teuchos::RCP<std::ostream> outputStream_;
261 
262  // Status test.
267 
268  // Current parameter list.
270 
276  mutable Teuchos::RCP<const Teuchos::ParameterList> validParams_;
277 
278  // Default solver values.
279  static constexpr int maxIters_default_ = 1000;
280  static constexpr bool showMaxResNormOnly_default_ = false;
281  static constexpr int verbosity_default_ = Belos::Errors;
282  static constexpr int outputStyle_default_ = Belos::General;
283  static constexpr int outputFreq_default_ = -1;
284  static constexpr int defQuorum_default_ = 1;
285  static constexpr const char * resScale_default_ = "Norm of Initial Residual";
286  static constexpr const char * label_default_ = "Belos";
287  static constexpr std::ostream * outputStream_default_ = &std::cout;
288 
289  // Current solver values.
290  MagnitudeType convtol_,achievedTol_;
291  int maxIters_, numIters_;
292  int verbosity_, outputStyle_, outputFreq_, defQuorum_;
293  bool showMaxResNormOnly_;
294  std::string resScale_;
295 
296  // Timers.
297  std::string label_;
298  Teuchos::RCP<Teuchos::Time> timerSolve_;
299 
300  // Internal state variables.
301  bool isSet_;
302  };
303 
304 // Empty Constructor
305 template<class ScalarType, class MV, class OP>
307  outputStream_(Teuchos::rcp(outputStream_default_,false)),
308  convtol_(DefaultSolverParameters::convTol),
309  maxIters_(maxIters_default_),
310  numIters_(0),
311  verbosity_(verbosity_default_),
312  outputStyle_(outputStyle_default_),
313  outputFreq_(outputFreq_default_),
314  defQuorum_(defQuorum_default_),
315  showMaxResNormOnly_(showMaxResNormOnly_default_),
316  resScale_(resScale_default_),
317  label_(label_default_),
318  isSet_(false)
319 {}
320 
321 // Basic Constructor
322 template<class ScalarType, class MV, class OP>
326  problem_(problem),
327  outputStream_(Teuchos::rcp(outputStream_default_,false)),
328  convtol_(DefaultSolverParameters::convTol),
329  maxIters_(maxIters_default_),
330  numIters_(0),
331  verbosity_(verbosity_default_),
332  outputStyle_(outputStyle_default_),
333  outputFreq_(outputFreq_default_),
334  defQuorum_(defQuorum_default_),
335  showMaxResNormOnly_(showMaxResNormOnly_default_),
336  resScale_(resScale_default_),
337  label_(label_default_),
338  isSet_(false)
339 {
341  problem_.is_null (), std::invalid_argument,
342  "Belos::BiCGStabSolMgr two-argument constructor: "
343  "'problem' is null. You must supply a non-null Belos::LinearProblem "
344  "instance when calling this constructor.");
345 
346  if (! pl.is_null ()) {
347  // Set the parameters using the list that was passed in.
348  setParameters (pl);
349  }
350 }
351 
352 template<class ScalarType, class MV, class OP>
354 {
356  using Teuchos::parameterList;
357  using Teuchos::RCP;
358 
359  RCP<const ParameterList> defaultParams = getValidParameters();
360 
361  // Create the internal parameter list if one doesn't already exist.
362  if (params_.is_null()) {
363  params_ = parameterList (*defaultParams);
364  } else {
365  params->validateParameters (*defaultParams);
366  }
367 
368  // Check for maximum number of iterations
369  if (params->isParameter("Maximum Iterations")) {
370  maxIters_ = params->get("Maximum Iterations",maxIters_default_);
371 
372  // Update parameter in our list and in status test.
373  params_->set("Maximum Iterations", maxIters_);
374  if (maxIterTest_!=Teuchos::null)
375  maxIterTest_->setMaxIters( maxIters_ );
376  }
377 
378  // Check to see if the timer label changed.
379  if (params->isParameter("Timer Label")) {
380  std::string tempLabel = params->get("Timer Label", label_default_);
381 
382  // Update parameter in our list and solver timer
383  if (tempLabel != label_) {
384  label_ = tempLabel;
385  params_->set("Timer Label", label_);
386  std::string solveLabel = label_ + ": BiCGStabSolMgr total solve time";
387 #ifdef BELOS_TEUCHOS_TIME_MONITOR
388  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
389 #endif
390  }
391  }
392 
393  // Check for a change in verbosity level
394  if (params->isParameter("Verbosity")) {
395  if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
396  verbosity_ = params->get("Verbosity", verbosity_default_);
397  } else {
398  verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
399  }
400 
401  // Update parameter in our list.
402  params_->set("Verbosity", verbosity_);
403  if (printer_ != Teuchos::null)
404  printer_->setVerbosity(verbosity_);
405  }
406 
407  // Check for a change in output style
408  if (params->isParameter("Output Style")) {
409  if (Teuchos::isParameterType<int>(*params,"Output Style")) {
410  outputStyle_ = params->get("Output Style", outputStyle_default_);
411  } else {
412  outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
413  }
414 
415  // Reconstruct the convergence test if the explicit residual test is not being used.
416  params_->set("Output Style", outputStyle_);
417  outputTest_ = Teuchos::null;
418  }
419 
420  // output stream
421  if (params->isParameter("Output Stream")) {
422  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
423 
424  // Update parameter in our list.
425  params_->set("Output Stream", outputStream_);
426  if (printer_ != Teuchos::null)
427  printer_->setOStream( outputStream_ );
428  }
429 
430  // frequency level
431  if (verbosity_ & Belos::StatusTestDetails) {
432  if (params->isParameter("Output Frequency")) {
433  outputFreq_ = params->get("Output Frequency", outputFreq_default_);
434  }
435 
436  // Update parameter in out list and output status test.
437  params_->set("Output Frequency", outputFreq_);
438  if (outputTest_ != Teuchos::null)
439  outputTest_->setOutputFrequency( outputFreq_ );
440  }
441 
442  // Create output manager if we need to.
443  if (printer_ == Teuchos::null) {
444  printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
445  }
446 
447  // Convergence
448  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
449  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestResNorm_t;
450 
451  // Check for convergence tolerance
452  if (params->isParameter("Convergence Tolerance")) {
453  if (params->isType<MagnitudeType> ("Convergence Tolerance")) {
454  convtol_ = params->get ("Convergence Tolerance",
455  static_cast<MagnitudeType> (DefaultSolverParameters::convTol));
456  }
457  else {
458  convtol_ = params->get ("Convergence Tolerance", DefaultSolverParameters::convTol);
459  }
460 
461  // Update parameter in our list and residual tests.
462  params_->set("Convergence Tolerance", convtol_);
463  if (convTest_ != Teuchos::null)
464  convTest_->setTolerance( convtol_ );
465  }
466 
467  if (params->isParameter("Show Maximum Residual Norm Only")) {
468  showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");
469 
470  // Update parameter in our list and residual tests
471  params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
472  if (convTest_ != Teuchos::null)
473  convTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
474  }
475 
476  // Check for a change in scaling, if so we need to build new residual tests.
477  bool newResTest = false;
478  {
479  // "Residual Scaling" is the old parameter name; "Implicit
480  // Residual Scaling" is the new name. We support both options for
481  // backwards compatibility.
482  std::string tempResScale = resScale_;
483  bool implicitResidualScalingName = false;
484  if (params->isParameter ("Residual Scaling")) {
485  tempResScale = params->get<std::string> ("Residual Scaling");
486  }
487  else if (params->isParameter ("Implicit Residual Scaling")) {
488  tempResScale = params->get<std::string> ("Implicit Residual Scaling");
489  implicitResidualScalingName = true;
490  }
491 
492  // Only update the scaling if it's different.
493  if (resScale_ != tempResScale) {
494  Belos::ScaleType resScaleType = convertStringToScaleType( tempResScale );
495  resScale_ = tempResScale;
496 
497  // Update parameter in our list and residual tests, using the
498  // given parameter name.
499  if (implicitResidualScalingName) {
500  params_->set ("Implicit Residual Scaling", resScale_);
501  }
502  else {
503  params_->set ("Residual Scaling", resScale_);
504  }
505 
506  if (! convTest_.is_null()) {
507  try {
508  convTest_->defineScaleForm( resScaleType, Belos::TwoNorm );
509  }
510  catch (std::exception& e) {
511  // Make sure the convergence test gets constructed again.
512  newResTest = true;
513  }
514  }
515  }
516  }
517 
518  // Get the deflation quorum, or number of converged systems before deflation is allowed
519  if (params->isParameter("Deflation Quorum")) {
520  defQuorum_ = params->get("Deflation Quorum", defQuorum_);
521  params_->set("Deflation Quorum", defQuorum_);
522  if (convTest_ != Teuchos::null)
523  convTest_->setQuorum( defQuorum_ );
524  }
525 
526  // Create status tests if we need to.
527 
528  // Basic test checks maximum iterations and native residual.
529  if (maxIterTest_ == Teuchos::null)
530  maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
531 
532  // Implicit residual test, using the native residual to determine if convergence was achieved.
533  if (convTest_ == Teuchos::null || newResTest) {
534  convTest_ = Teuchos::rcp( new StatusTestResNorm_t( convtol_, defQuorum_, showMaxResNormOnly_ ) );
535  convTest_->defineScaleForm( convertStringToScaleType( resScale_ ), Belos::TwoNorm );
536  }
537 
538  if (sTest_ == Teuchos::null || newResTest)
539  sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
540 
541  if (outputTest_ == Teuchos::null || newResTest) {
542 
543  // Create the status test output class.
544  // This class manages and formats the output from the status test.
545  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
546  outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
547 
548  // Set the solver string for the output test
549  std::string solverDesc = " Pseudo Block BiCGStab ";
550  outputTest_->setSolverDesc( solverDesc );
551 
552  }
553 
554  // Create the timer if we need to.
555  if (timerSolve_ == Teuchos::null) {
556  std::string solveLabel = label_ + ": BiCGStabSolMgr total solve time";
557 #ifdef BELOS_TEUCHOS_TIME_MONITOR
558  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
559 #endif
560  }
561 
562  // Inform the solver manager that the current parameters were set.
563  isSet_ = true;
564 }
565 
566 
567 template<class ScalarType, class MV, class OP>
570 {
572  using Teuchos::parameterList;
573  using Teuchos::RCP;
574 
575  if (validParams_.is_null()) {
576  // Set all the valid parameters and their default values.
577  RCP<ParameterList> pl = parameterList ();
578 
579  // The static_cast is to resolve an issue with older clang versions which
580  // would cause the constexpr to link fail. With c++17 the problem is resolved.
581  pl->set("Convergence Tolerance", static_cast<MagnitudeType>(DefaultSolverParameters::convTol),
582  "The relative residual tolerance that needs to be achieved by the\n"
583  "iterative solver in order for the linera system to be declared converged.");
584  pl->set("Maximum Iterations", static_cast<int>(maxIters_default_),
585  "The maximum number of block iterations allowed for each\n"
586  "set of RHS solved.");
587  pl->set("Verbosity", static_cast<int>(verbosity_default_),
588  "What type(s) of solver information should be outputted\n"
589  "to the output stream.");
590  pl->set("Output Style", static_cast<int>(outputStyle_default_),
591  "What style is used for the solver information outputted\n"
592  "to the output stream.");
593  pl->set("Output Frequency", static_cast<int>(outputFreq_default_),
594  "How often convergence information should be outputted\n"
595  "to the output stream.");
596  pl->set("Deflation Quorum", static_cast<int>(defQuorum_default_),
597  "The number of linear systems that need to converge before\n"
598  "they are deflated. This number should be <= block size.");
599  pl->set("Output Stream", Teuchos::rcp(outputStream_default_,false),
600  "A reference-counted pointer to the output stream where all\n"
601  "solver output is sent.");
602  pl->set("Show Maximum Residual Norm Only", static_cast<bool>(showMaxResNormOnly_default_),
603  "When convergence information is printed, only show the maximum\n"
604  "relative residual norm when the block size is greater than one.");
605  pl->set("Implicit Residual Scaling", static_cast<const char *>(resScale_default_),
606  "The type of scaling used in the residual convergence test.");
607  // We leave the old name as a valid parameter for backwards
608  // compatibility (so that validateParametersAndSetDefaults()
609  // doesn't raise an exception if it encounters "Residual
610  // Scaling"). The new name was added for compatibility with other
611  // solvers, none of which use "Residual Scaling".
612  pl->set("Residual Scaling", static_cast<const char *>(resScale_default_),
613  "The type of scaling used in the residual convergence test. This "
614  "name is deprecated; the new name is \"Implicit Residual Scaling\".");
615  pl->set("Timer Label", static_cast<const char *>(label_default_),
616  "The string to use as a prefix for the timer labels.");
617  validParams_ = pl;
618  }
619  return validParams_;
620 }
621 
622 
623 template<class ScalarType, class MV, class OP>
625 {
626  // Set the current parameters if they were not set before.
627  // NOTE: This may occur if the user generated the solver manager with the default constructor and
628  // then didn't set any parameters using setParameters().
629  if (! isSet_) {
630  setParameters (params_);
631  }
632 
634 
636  (! problem_->isProblemSet (), BiCGStabSolMgrLinearProblemFailure,
637  "Belos::BiCGStabSolMgr::solve: Linear problem is not ready. "
638  "You must call setProblem() on the LinearProblem before you may solve it.");
640  (problem_->isLeftPrec (), std::logic_error, "Belos::BiCGStabSolMgr::solve: "
641  "The left-preconditioned case has not yet been implemented. Please use "
642  "right preconditioning for now. If you need to use left preconditioning, "
643  "please contact the Belos developers. Left preconditioning is more "
644  "interesting in BiCGStab because whether it works depends on the initial "
645  "guess (e.g., an initial guess of all zeros might NOT work).");
646 
647  // Create indices for the linear systems to be solved.
648  int startPtr = 0;
649  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
650  int numCurrRHS = numRHS2Solve;
651 
652  std::vector<int> currIdx( numRHS2Solve ), currIdx2( numRHS2Solve );
653  for (int i=0; i<numRHS2Solve; ++i) {
654  currIdx[i] = startPtr+i;
655  currIdx2[i]=i;
656  }
657 
658  // Inform the linear problem of the current linear system to solve.
659  problem_->setLSIndex( currIdx );
660 
662  // Parameter list (iteration)
664 
665  // Reset the status test.
666  outputTest_->reset();
667 
668  // Assume convergence is achieved, then let any failed convergence set this to false.
669  bool isConverged = true;
670 
672  // Pseudo-Block BiCGStab solver
673 
675  = Teuchos::rcp( new BiCGStabIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,plist) );
676 
677  // Enter solve() iterations
678  {
679 #ifdef BELOS_TEUCHOS_TIME_MONITOR
680  Teuchos::TimeMonitor slvtimer(*timerSolve_);
681 #endif
682 
683  //bool first_time=true;
684  while ( numRHS2Solve > 0 ) {
685  // Reset the active / converged vectors from this block
686  std::vector<int> convRHSIdx;
687  std::vector<int> currRHSIdx( currIdx );
688  currRHSIdx.resize(numCurrRHS);
689 
690  // Reset the number of iterations.
691  block_cg_iter->resetNumIters();
692 
693  // Reset the number of calls that the status test output knows about.
694  outputTest_->resetNumCalls();
695 
696  // Get the current residual for this block of linear systems.
697  Teuchos::RCP<MV> R_0 = MVT::CloneViewNonConst( *(Teuchos::rcp_const_cast<MV>(problem_->getInitResVec())), currIdx );
698 
699  // Get a new state struct and initialize the solver.
701  newState.R = R_0;
702  block_cg_iter->initializeBiCGStab(newState);
703 
704  while(1) {
705 
706  // tell block_gmres_iter to iterate
707  try {
708 
709  block_cg_iter->iterate();
710 
712  //
713  // check convergence first
714  //
716  if ( convTest_->getStatus() == Passed ) {
717 
718  // Figure out which linear systems converged.
719  std::vector<int> convIdx = Teuchos::rcp_dynamic_cast<StatusTestGenResNorm<ScalarType,MV,OP> >(convTest_)->convIndices();
720 
721  // If the number of converged linear systems is equal to the
722  // number of current linear systems, then we are done with this block.
723  if (convIdx.size() == currRHSIdx.size())
724  break; // break from while(1){block_cg_iter->iterate()}
725 
726  // Inform the linear problem that we are finished with this current linear system.
727  problem_->setCurrLS();
728 
729  // Reset currRHSIdx to have the right-hand sides that are left to converge for this block.
730  int have = 0;
731  std::vector<int> unconvIdx(currRHSIdx.size());
732  for (unsigned int i=0; i<currRHSIdx.size(); ++i) {
733  bool found = false;
734  for (unsigned int j=0; j<convIdx.size(); ++j) {
735  if (currRHSIdx[i] == convIdx[j]) {
736  found = true;
737  break;
738  }
739  }
740  if (!found) {
741  currIdx2[have] = currIdx2[i];
742  currRHSIdx[have++] = currRHSIdx[i];
743  }
744  }
745  currRHSIdx.resize(have);
746  currIdx2.resize(have);
747 
748  // Set the remaining indices after deflation.
749  problem_->setLSIndex( currRHSIdx );
750 
751  // Get the current residual vector.
752  std::vector<MagnitudeType> norms;
753  R_0 = MVT::CloneCopy( *(block_cg_iter->getNativeResiduals(&norms)),currIdx2 );
754  for (int i=0; i<have; ++i) { currIdx2[i] = i; }
755 
756  // Set the new state and initialize the solver.
758  defstate.R = R_0;
759  block_cg_iter->initializeBiCGStab(defstate);
760  }
761 
763  //
764  // check for maximum iterations
765  //
767  else if ( maxIterTest_->getStatus() == Passed ) {
768  // we don't have convergence
769  isConverged = false;
770  break; // break from while(1){block_cg_iter->iterate()}
771  }
772 
774  //
775  // we returned from iterate(), but none of our status tests Passed.
776  // something is wrong, and it is probably our fault.
777  //
779 
780  else {
781  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
782  "Belos::BiCGStabSolMgr::solve(): Invalid return from BiCGStabIter::iterate().");
783  }
784  }
785  catch (const std::exception &e) {
786  printer_->stream(Errors) << "Error! Caught std::exception in BiCGStabIter::iterate() at iteration "
787  << block_cg_iter->getNumIters() << std::endl
788  << e.what() << std::endl;
789  throw;
790  }
791  }
792 
793  // Inform the linear problem that we are finished with this block linear system.
794  problem_->setCurrLS();
795 
796  // Update indices for the linear systems to be solved.
797  startPtr += numCurrRHS;
798  numRHS2Solve -= numCurrRHS;
799 
800  if ( numRHS2Solve > 0 ) {
801 
802  numCurrRHS = numRHS2Solve;
803  currIdx.resize( numCurrRHS );
804  currIdx2.resize( numCurrRHS );
805  for (int i=0; i<numCurrRHS; ++i)
806  { currIdx[i] = startPtr+i; currIdx2[i] = i; }
807 
808  // Set the next indices.
809  problem_->setLSIndex( currIdx );
810  }
811  else {
812  currIdx.resize( numRHS2Solve );
813  }
814 
815  //first_time=false;
816  }// while ( numRHS2Solve > 0 )
817 
818  }
819 
820  // print final summary
821  sTest_->print( printer_->stream(FinalSummary) );
822 
823  // print timing information
824 #ifdef BELOS_TEUCHOS_TIME_MONITOR
825  // Calling summarize() can be expensive, so don't call unless the
826  // user wants to print out timing details. summarize() will do all
827  // the work even if it's passed a "black hole" output stream.
828  if (verbosity_ & TimingDetails)
829  Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
830 #endif
831 
832  // get iteration information for this solve
833  numIters_ = maxIterTest_->getNumIters();
834 
835 
836  // Save the convergence test value ("achieved tolerance") for this
837  // solve.
838  const std::vector<MagnitudeType>* pTestValues = convTest_->getTestValue();
839  achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
840 
841 
842  if (!isConverged ) {
843  return Unconverged; // return from BiCGStabSolMgr::solve()
844  }
845  return Converged; // return from BiCGStabSolMgr::solve()
846 }
847 
848 // This method requires the solver manager to return a std::string that describes itself.
849 template<class ScalarType, class MV, class OP>
851 {
852  std::ostringstream oss;
853  oss << "Belos::BiCGStabSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
854  oss << "{";
855  oss << "}";
856  return oss.str();
857 }
858 
859 
860 
861 } // end Belos namespace
862 
863 #endif /* BELOS_BICGSTAB_SOLMGR_HPP */
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Return a reference to the linear problem being solved by this solver manager.
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.
BiCGStabSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate orthono...
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
This class implements the pseudo-block BiCGStab iteration, where the basic BiCGStab algorithm is perf...
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.
Teuchos::RCP< const MV > R
The current residual.
ScaleType
The type of scaling to use on the residual norm value.
Definition: BelosTypes.hpp:120
T & get(const std::string &name, T def_value)
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)
A factory class for generating StatusTestOutput objects.
Structure to contain pointers to BiCGStabIteration state variables.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
An implementation of StatusTestResNorm using a family of residual norms.
void reset(const ResetType type) override
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
static const double convTol
Default convergence tolerance.
Definition: BelosTypes.hpp:293
Belos::StatusTest class for specifying a maximum number of iterations.
static std::string name()
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
A factory class for generating StatusTestOutput objects.
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
BiCGStabSolMgrLinearProblemFailure(const std::string &what_arg)
bool isParameter(const std::string &name) const
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 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)
BiCGStabSolMgr()
Empty constructor for BiCGStabSolMgr. This constructor takes no arguments and sets the default values...
virtual ~BiCGStabSolMgr()
Destructor.
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
BiCGStabSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i...
ReturnType solve() override
This method performs possibly repeated calls to the underlying linear solver&#39;s iterate() routine unti...
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:155
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Get a parameter list containing the valid parameters for this object.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
The Belos::BiCGStabSolMgr provides a powerful and fully-featured solver manager over the pseudo-block...
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.
std::string description() const override
Method to return description of the block BiCGStab solver manager.
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.
Class which defines basic traits for the operator type.
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:60
int getNumIters() const override
Get the iteration count for the most recent call to solve().
Default parameters common to most Belos solvers.
Definition: BelosTypes.hpp:283
bool isLOADetected() const override
Return whether a loss of accuracy was detected by this solver during the most current solve...
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
Belos header file which uses auto-configuration information to include necessary C++ headers...
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
BiCGStabSolMgrOrthoFailure(const std::string &what_arg)
Belos concrete class for performing the pseudo-block BiCGStab iteration.
bool is_null() const

Generated on Thu Apr 25 2024 09:27:23 for Belos by doxygen 1.8.5