Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosPseudoBlockStochasticCGSolMgr.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_PSEUDO_BLOCK_STOCHASTIC_CG_SOLMGR_HPP
43 #define BELOS_PSEUDO_BLOCK_STOCHASTIC_CG_SOLMGR_HPP
44 
49 #include "BelosConfigDefs.hpp"
50 #include "BelosTypes.hpp"
51 
52 #include "BelosLinearProblem.hpp"
53 #include "BelosSolverManager.hpp"
54 
58 #include "BelosStatusTestCombo.hpp"
60 #include "BelosOutputManager.hpp"
61 #include "Teuchos_BLAS.hpp"
62 #ifdef BELOS_TEUCHOS_TIME_MONITOR
63 #include "Teuchos_TimeMonitor.hpp"
64 #endif
65 
75 namespace Belos {
76 
78 
79 
87  PseudoBlockStochasticCGSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
88  {}};
89 
97  PseudoBlockStochasticCGSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
98  {}};
99 
100  template<class ScalarType, class MV, class OP>
101  class PseudoBlockStochasticCGSolMgr : public SolverManager<ScalarType,MV,OP> {
102 
103  private:
107  typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
109 
110  public:
111 
113 
114 
121 
133 
136 
140  }
142 
144 
145 
146  const LinearProblem<ScalarType,MV,OP>& getProblem() const override {
147  return *problem_;
148  }
149 
153 
157 
164  return Teuchos::tuple(timerSolve_);
165  }
166 
168  int getNumIters() const override {
169  return numIters_;
170  }
171 
175  bool isLOADetected() const override { return false; }
176 
178 
180 
181 
183  void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) override { problem_ = problem; }
184 
186  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params ) override;
187 
189 
191 
192 
196  void reset( const ResetType type ) override { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
198 
200 
201 
219  ReturnType solve() override;
220 
222 
225 
228 
230  std::string description() const override;
231 
233 
234  private:
235 
236  // Linear problem.
238 
239  // Output manager.
241  Teuchos::RCP<std::ostream> outputStream_;
242 
243  // Status test.
248 
249  // Current parameter list.
251 
257  mutable Teuchos::RCP<const Teuchos::ParameterList> validParams_;
258 
259  // Default solver values.
260  static constexpr int maxIters_default_ = 1000;
261  static constexpr bool assertPositiveDefiniteness_default_ = true;
262  static constexpr bool showMaxResNormOnly_default_ = false;
263  static constexpr int verbosity_default_ = Belos::Errors;
264  static constexpr int outputStyle_default_ = Belos::General;
265  static constexpr int outputFreq_default_ = -1;
266  static constexpr int defQuorum_default_ = 1;
267  static constexpr const char * resScale_default_ = "Norm of Initial Residual";
268  static constexpr const char * label_default_ = "Belos";
269  static constexpr std::ostream * outputStream_default_ = &std::cout;
270 
271  // Current solver values.
272  MagnitudeType convtol_;
273  int maxIters_, numIters_;
274  int verbosity_, outputStyle_, outputFreq_, defQuorum_;
275  bool assertPositiveDefiniteness_, showMaxResNormOnly_;
276  std::string resScale_;
277 
278  // Timers.
279  std::string label_;
280  Teuchos::RCP<Teuchos::Time> timerSolve_;
281 
282  // Internal state variables.
283  bool isSet_;
284 
285  // Stashed copy of the stochastic vector
286  Teuchos::RCP<MV> Y_;
287 
288  };
289 
290 
291 // Empty Constructor
292 template<class ScalarType, class MV, class OP>
294  outputStream_(Teuchos::rcp(outputStream_default_,false)),
295  convtol_(DefaultSolverParameters::convTol),
296  maxIters_(maxIters_default_),
297  numIters_(0),
298  verbosity_(verbosity_default_),
299  outputStyle_(outputStyle_default_),
300  outputFreq_(outputFreq_default_),
301  defQuorum_(defQuorum_default_),
302  assertPositiveDefiniteness_(assertPositiveDefiniteness_default_),
303  showMaxResNormOnly_(showMaxResNormOnly_default_),
304  resScale_(resScale_default_),
305  label_(label_default_),
306  isSet_(false)
307 {}
308 
309 // Basic Constructor
310 template<class ScalarType, class MV, class OP>
314  problem_(problem),
315  outputStream_(Teuchos::rcp(outputStream_default_,false)),
316  convtol_(DefaultSolverParameters::convTol),
317  maxIters_(maxIters_default_),
318  numIters_(0),
319  verbosity_(verbosity_default_),
320  outputStyle_(outputStyle_default_),
321  outputFreq_(outputFreq_default_),
322  defQuorum_(defQuorum_default_),
323  assertPositiveDefiniteness_(assertPositiveDefiniteness_default_),
324  showMaxResNormOnly_(showMaxResNormOnly_default_),
325  resScale_(resScale_default_),
326  label_(label_default_),
327  isSet_(false)
328 {
330  problem_.is_null (), std::invalid_argument,
331  "Belos::PseudoBlockStochasticCGSolMgr two-argument constructor: "
332  "'problem' is null. You must supply a non-null Belos::LinearProblem "
333  "instance when calling this constructor.");
334 
335  if (! pl.is_null ()) {
336  // Set the parameters using the list that was passed in.
337  setParameters (pl);
338  }
339 }
340 
341 template<class ScalarType, class MV, class OP>
343 {
345  using Teuchos::parameterList;
346  using Teuchos::RCP;
347 
348  RCP<const ParameterList> defaultParams = getValidParameters();
349 
350  // Create the internal parameter list if one doesn't already exist.
351  if (params_.is_null()) {
352  params_ = parameterList (*defaultParams);
353  } else {
354  params->validateParameters (*defaultParams);
355  }
356 
357  // Check for maximum number of iterations
358  if (params->isParameter("Maximum Iterations")) {
359  maxIters_ = params->get("Maximum Iterations",maxIters_default_);
360 
361  // Update parameter in our list and in status test.
362  params_->set("Maximum Iterations", maxIters_);
363  if (maxIterTest_!=Teuchos::null)
364  maxIterTest_->setMaxIters( maxIters_ );
365  }
366 
367  // Check if positive definiteness assertions are to be performed
368  if (params->isParameter("Assert Positive Definiteness")) {
369  assertPositiveDefiniteness_ = params->get("Assert Positive Definiteness",assertPositiveDefiniteness_default_);
370 
371  // Update parameter in our list.
372  params_->set("Assert Positive Definiteness", assertPositiveDefiniteness_);
373  }
374 
375  // Check to see if the timer label changed.
376  if (params->isParameter("Timer Label")) {
377  std::string tempLabel = params->get("Timer Label", label_default_);
378 
379  // Update parameter in our list and solver timer
380  if (tempLabel != label_) {
381  label_ = tempLabel;
382  params_->set("Timer Label", label_);
383  std::string solveLabel = label_ + ": PseudoBlockStochasticCGSolMgr total solve time";
384 #ifdef BELOS_TEUCHOS_TIME_MONITOR
385  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
386 #endif
387  }
388  }
389 
390  // Check for a change in verbosity level
391  if (params->isParameter("Verbosity")) {
392  if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
393  verbosity_ = params->get("Verbosity", verbosity_default_);
394  } else {
395  verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
396  }
397 
398  // Update parameter in our list.
399  params_->set("Verbosity", verbosity_);
400  if (printer_ != Teuchos::null)
401  printer_->setVerbosity(verbosity_);
402  }
403 
404  // Check for a change in output style
405  if (params->isParameter("Output Style")) {
406  if (Teuchos::isParameterType<int>(*params,"Output Style")) {
407  outputStyle_ = params->get("Output Style", outputStyle_default_);
408  } else {
409  outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
410  }
411 
412  // Reconstruct the convergence test if the explicit residual test is not being used.
413  params_->set("Output Style", outputStyle_);
414  outputTest_ = Teuchos::null;
415  }
416 
417  // output stream
418  if (params->isParameter("Output Stream")) {
419  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
420 
421  // Update parameter in our list.
422  params_->set("Output Stream", outputStream_);
423  if (printer_ != Teuchos::null)
424  printer_->setOStream( outputStream_ );
425  }
426 
427  // frequency level
428  if (verbosity_ & Belos::StatusTestDetails) {
429  if (params->isParameter("Output Frequency")) {
430  outputFreq_ = params->get("Output Frequency", outputFreq_default_);
431  }
432 
433  // Update parameter in out list and output status test.
434  params_->set("Output Frequency", outputFreq_);
435  if (outputTest_ != Teuchos::null)
436  outputTest_->setOutputFrequency( outputFreq_ );
437  }
438 
439  // Create output manager if we need to.
440  if (printer_ == Teuchos::null) {
441  printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
442  }
443 
444  // Convergence
445  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
446  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestResNorm_t;
447 
448  // Check for convergence tolerance
449  if (params->isParameter("Convergence Tolerance")) {
450  if (params->isType<MagnitudeType> ("Convergence Tolerance")) {
451  convtol_ = params->get ("Convergence Tolerance",
452  static_cast<MagnitudeType> (DefaultSolverParameters::convTol));
453  }
454  else {
455  convtol_ = params->get ("Convergence Tolerance", DefaultSolverParameters::convTol);
456  }
457 
458  // Update parameter in our list and residual tests.
459  params_->set("Convergence Tolerance", convtol_);
460  if (convTest_ != Teuchos::null)
461  convTest_->setTolerance( convtol_ );
462  }
463 
464  if (params->isParameter("Show Maximum Residual Norm Only")) {
465  showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");
466 
467  // Update parameter in our list and residual tests
468  params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
469  if (convTest_ != Teuchos::null)
470  convTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
471  }
472 
473  // Check for a change in scaling, if so we need to build new residual tests.
474  bool newResTest = false;
475  {
476  // "Residual Scaling" is the old parameter name; "Implicit
477  // Residual Scaling" is the new name. We support both options for
478  // backwards compatibility.
479  std::string tempResScale = resScale_;
480  bool implicitResidualScalingName = false;
481  if (params->isParameter ("Residual Scaling")) {
482  tempResScale = params->get<std::string> ("Residual Scaling");
483  }
484  else if (params->isParameter ("Implicit Residual Scaling")) {
485  tempResScale = params->get<std::string> ("Implicit Residual Scaling");
486  implicitResidualScalingName = true;
487  }
488 
489  // Only update the scaling if it's different.
490  if (resScale_ != tempResScale) {
491  Belos::ScaleType resScaleType = convertStringToScaleType( tempResScale );
492  resScale_ = tempResScale;
493 
494  // Update parameter in our list and residual tests, using the
495  // given parameter name.
496  if (implicitResidualScalingName) {
497  params_->set ("Implicit Residual Scaling", resScale_);
498  }
499  else {
500  params_->set ("Residual Scaling", resScale_);
501  }
502 
503  if (! convTest_.is_null()) {
504  try {
505  convTest_->defineScaleForm( resScaleType, Belos::TwoNorm );
506  }
507  catch (std::exception& e) {
508  // Make sure the convergence test gets constructed again.
509  newResTest = true;
510  }
511  }
512  }
513  }
514 
515  // Get the deflation quorum, or number of converged systems before deflation is allowed
516  if (params->isParameter("Deflation Quorum")) {
517  defQuorum_ = params->get("Deflation Quorum", defQuorum_);
518  params_->set("Deflation Quorum", defQuorum_);
519  if (convTest_ != Teuchos::null)
520  convTest_->setQuorum( defQuorum_ );
521  }
522 
523  // Create status tests if we need to.
524 
525  // Basic test checks maximum iterations and native residual.
526  if (maxIterTest_ == Teuchos::null)
527  maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
528 
529  // Implicit residual test, using the native residual to determine if convergence was achieved.
530  if (convTest_ == Teuchos::null || newResTest) {
531  convTest_ = Teuchos::rcp( new StatusTestResNorm_t( convtol_, defQuorum_, showMaxResNormOnly_ ) );
532  convTest_->defineScaleForm( convertStringToScaleType( resScale_ ), Belos::TwoNorm );
533  }
534 
535  if (sTest_ == Teuchos::null || newResTest)
536  sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
537 
538  if (outputTest_ == Teuchos::null || newResTest) {
539 
540  // Create the status test output class.
541  // This class manages and formats the output from the status test.
542  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
543  outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
544 
545  // Set the solver string for the output test
546  std::string solverDesc = " Pseudo Block CG ";
547  outputTest_->setSolverDesc( solverDesc );
548 
549  }
550 
551  // Create the timer if we need to.
552  if (timerSolve_ == Teuchos::null) {
553  std::string solveLabel = label_ + ": PseudoBlockStochasticCGSolMgr total solve time";
554 #ifdef BELOS_TEUCHOS_TIME_MONITOR
555  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
556 #endif
557  }
558 
559  // Inform the solver manager that the current parameters were set.
560  isSet_ = true;
561 }
562 
563 
564 template<class ScalarType, class MV, class OP>
567 {
569  using Teuchos::parameterList;
570  using Teuchos::RCP;
571 
572  if (validParams_.is_null()) {
573  // Set all the valid parameters and their default values.
574 
575  // The static_cast is to resolve an issue with older clang versions which
576  // would cause the constexpr to link fail. With c++17 the problem is resolved.
577  RCP<ParameterList> pl = parameterList ();
578  pl->set("Convergence Tolerance", static_cast<MagnitudeType>(DefaultSolverParameters::convTol),
579  "The relative residual tolerance that needs to be achieved by the\n"
580  "iterative solver in order for the linera system to be declared converged.");
581  pl->set("Maximum Iterations", static_cast<int>(maxIters_default_),
582  "The maximum number of block iterations allowed for each\n"
583  "set of RHS solved.");
584  pl->set("Assert Positive Definiteness", static_cast<bool>(assertPositiveDefiniteness_default_),
585  "Whether or not to assert that the linear operator\n"
586  "and the preconditioner are indeed positive definite.");
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", 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", 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 // solve()
624 template<class ScalarType, class MV, class OP>
626 
627  // Set the current parameters if they were not set before.
628  // NOTE: This may occur if the user generated the solver manager with the default constructor and
629  // then didn't set any parameters using setParameters().
630  if (!isSet_) { setParameters( params_ ); }
631 
633 
635  "Belos::PseudoBlockStochasticCGSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
636 
637  // Create indices for the linear systems to be solved.
638  int startPtr = 0;
639  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
640  int numCurrRHS = numRHS2Solve;
641 
642  std::vector<int> currIdx( numRHS2Solve ), currIdx2( numRHS2Solve );
643  for (int i=0; i<numRHS2Solve; ++i) {
644  currIdx[i] = startPtr+i;
645  currIdx2[i]=i;
646  }
647 
648  // Inform the linear problem of the current linear system to solve.
649  problem_->setLSIndex( currIdx );
650 
652  // Parameter list
654 
655  plist.set("Assert Positive Definiteness",assertPositiveDefiniteness_);
656 
657  // Reset the status test.
658  outputTest_->reset();
659 
660  // Assume convergence is achieved, then let any failed convergence set this to false.
661  bool isConverged = true;
662 
664  // Pseudo-Block CG solver
665 
667  = Teuchos::rcp( new PseudoBlockStochasticCGIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,plist) );
668 
669  // Enter solve() iterations
670  {
671 #ifdef BELOS_TEUCHOS_TIME_MONITOR
672  Teuchos::TimeMonitor slvtimer(*timerSolve_);
673 #endif
674 
675  while ( numRHS2Solve > 0 ) {
676 
677  // Reset the active / converged vectors from this block
678  std::vector<int> convRHSIdx;
679  std::vector<int> currRHSIdx( currIdx );
680  currRHSIdx.resize(numCurrRHS);
681 
682  // Reset the number of iterations.
683  block_cg_iter->resetNumIters();
684 
685  // Reset the number of calls that the status test output knows about.
686  outputTest_->resetNumCalls();
687 
688  // Get the current residual for this block of linear systems.
689  Teuchos::RCP<MV> R_0 = MVT::CloneViewNonConst( *(Teuchos::rcp_const_cast<MV>(problem_->getInitResVec())), currIdx );
690 
691  // Get a new state struct and initialize the solver.
693  newState.R = R_0;
694  block_cg_iter->initializeCG(newState);
695 
696  while(1) {
697 
698  // tell block_gmres_iter to iterate
699  try {
700  block_cg_iter->iterate();
701 
703  //
704  // check convergence first
705  //
707  if ( convTest_->getStatus() == Passed ) {
708 
709  // Figure out which linear systems converged.
710  std::vector<int> convIdx = Teuchos::rcp_dynamic_cast<StatusTestGenResNorm<ScalarType,MV,OP> >(convTest_)->convIndices();
711 
712  // If the number of converged linear systems is equal to the
713  // number of current linear systems, then we are done with this block.
714  if (convIdx.size() == currRHSIdx.size())
715  break; // break from while(1){block_cg_iter->iterate()}
716 
717  // Inform the linear problem that we are finished with this current linear system.
718  problem_->setCurrLS();
719 
720  // Reset currRHSIdx to have the right-hand sides that are left to converge for this block.
721  int have = 0;
722  std::vector<int> unconvIdx(currRHSIdx.size());
723  for (unsigned int i=0; i<currRHSIdx.size(); ++i) {
724  bool found = false;
725  for (unsigned int j=0; j<convIdx.size(); ++j) {
726  if (currRHSIdx[i] == convIdx[j]) {
727  found = true;
728  break;
729  }
730  }
731  if (!found) {
732  currIdx2[have] = currIdx2[i];
733  currRHSIdx[have++] = currRHSIdx[i];
734  }
735  }
736  currRHSIdx.resize(have);
737  currIdx2.resize(have);
738 
739  // Set the remaining indices after deflation.
740  problem_->setLSIndex( currRHSIdx );
741 
742  // Get the current residual vector.
743  std::vector<MagnitudeType> norms;
744  R_0 = MVT::CloneCopy( *(block_cg_iter->getNativeResiduals(&norms)),currIdx2 );
745  for (int i=0; i<have; ++i) { currIdx2[i] = i; }
746 
747  // Set the new state and initialize the solver.
749  defstate.R = R_0;
750  block_cg_iter->initializeCG(defstate);
751  }
752 
754  //
755  // check for maximum iterations
756  //
758  else if ( maxIterTest_->getStatus() == Passed ) {
759  // we don't have convergence
760  isConverged = false;
761  break; // break from while(1){block_cg_iter->iterate()}
762  }
763 
765  //
766  // we returned from iterate(), but none of our status tests Passed.
767  // something is wrong, and it is probably our fault.
768  //
770 
771  else {
772  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
773  "Belos::PseudoBlockStochasticCGSolMgr::solve(): Invalid return from PseudoBlockStochasticCGIter::iterate().");
774  }
775  }
776  catch (const std::exception &e) {
777  printer_->stream(Errors) << "Error! Caught std::exception in PseudoBlockStochasticCGIter::iterate() at iteration "
778  << block_cg_iter->getNumIters() << std::endl
779  << e.what() << std::endl;
780  throw;
781  }
782  }
783 
784  // Inform the linear problem that we are finished with this block linear system.
785  problem_->setCurrLS();
786 
787  // Update indices for the linear systems to be solved.
788  startPtr += numCurrRHS;
789  numRHS2Solve -= numCurrRHS;
790 
791  if ( numRHS2Solve > 0 ) {
792 
793  numCurrRHS = numRHS2Solve;
794  currIdx.resize( numCurrRHS );
795  currIdx2.resize( numCurrRHS );
796  for (int i=0; i<numCurrRHS; ++i)
797  { currIdx[i] = startPtr+i; currIdx2[i] = i; }
798 
799  // Set the next indices.
800  problem_->setLSIndex( currIdx );
801  }
802  else {
803  currIdx.resize( numRHS2Solve );
804  }
805 
806  }// while ( numRHS2Solve > 0 )
807 
808  }
809 
810  // get the final stochastic vector
811  Y_=block_cg_iter->getStochasticVector();
812 
813 
814  // print final summary
815  sTest_->print( printer_->stream(FinalSummary) );
816 
817  // print timing information
818 #ifdef BELOS_TEUCHOS_TIME_MONITOR
819  // Calling summarize() can be expensive, so don't call unless the
820  // user wants to print out timing details. summarize() will do all
821  // the work even if it's passed a "black hole" output stream.
822  if (verbosity_ & TimingDetails)
823  Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
824 #endif
825 
826  // get iteration information for this solve
827  numIters_ = maxIterTest_->getNumIters();
828 
829  if (!isConverged ) {
830  return Unconverged; // return from PseudoBlockStochasticCGSolMgr::solve()
831  }
832  return Converged; // return from PseudoBlockStochasticCGSolMgr::solve()
833 }
834 
835 // This method requires the solver manager to return a std::string that describes itself.
836 template<class ScalarType, class MV, class OP>
838 {
839  std::ostringstream oss;
840  oss << "Belos::PseudoBlockStochasticCGSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
841  oss << "{";
842  oss << "}";
843  return oss.str();
844 }
845 
846 } // end Belos namespace
847 
848 #endif /* BELOS_PSEUDO_BLOCK_CG_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.
bool isLOADetected() const override
Return whether a loss of accuracy was detected by this solver during the most current 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.
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)
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)
PseudoBlockStochasticCGSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to g...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
A factory class for generating StatusTestOutput objects.
int getNumIters() const override
Get the iteration count for the most recent call to solve().
std::string description() const override
Method to return description of the block CG solver manager.
An implementation of StatusTestResNorm using a family of residual norms.
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()
A factory class for generating StatusTestOutput objects.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Return a reference to the linear problem being solved by this solver manager.
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.
Belos concrete class for performing the stochastic pseudo-block CG iteration.
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)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Get a parameter list containing the valid parameters for this object.
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
Set the parameters the solver manager should use to solve the linear problem.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
Teuchos::RCP< const MV > R
The current residual.
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 ...
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
PseudoBlockStochasticCGSolMgr()
Empty constructor for BlockStochasticCGSolMgr. This constructor takes no arguments and sets the defau...
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.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
ReturnType solve() override
This method performs possibly repeated calls to the underlying linear solver&#39;s iterate() routine unti...
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
This class implements the stochastic pseudo-block CG iteration, where the basic stochastic CG algorit...
The Belos::PseudoBlockStochasticCGSolMgr provides a powerful and fully-featured solver manager over t...
bool isType(const std::string &name) const
Teuchos::RCP< MV > getStochasticVector()
Get a copy of the final stochastic vector.
A class for extending the status testing capabilities of Belos via logical combinations.
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 defines basic traits for the operator type.
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:60
Default parameters common to most Belos solvers.
Definition: BelosTypes.hpp:283
Belos header file which uses auto-configuration information to include necessary C++ headers...
PseudoBlockStochasticCGSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i...
Structure to contain pointers to CGIteration state variables.
bool is_null() const

Generated on Thu Apr 18 2024 09:27:54 for Belos by doxygen 1.8.5