Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosBlockCGSolMgr.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_BLOCK_CG_SOLMGR_HPP
43 #define BELOS_BLOCK_CG_SOLMGR_HPP
44 
49 #include "BelosConfigDefs.hpp"
50 #include "BelosTypes.hpp"
51 
52 #include "BelosLinearProblem.hpp"
53 #include "BelosSolverManager.hpp"
54 
55 #include "BelosCGIter.hpp"
56 #include "BelosCGSingleRedIter.hpp"
57 #include "BelosBlockCGIter.hpp"
61 #include "BelosStatusTestCombo.hpp"
63 #include "BelosOutputManager.hpp"
64 #include "Teuchos_LAPACK.hpp"
65 #ifdef BELOS_TEUCHOS_TIME_MONITOR
66 # include "Teuchos_TimeMonitor.hpp"
67 #endif
68 #if defined(HAVE_TEUCHOSCORE_CXX11)
69 # include <type_traits>
70 #endif // defined(HAVE_TEUCHOSCORE_CXX11)
71 #include <algorithm>
72 
89 namespace Belos {
90 
92 
93 
101  BlockCGSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
102  {}};
103 
104  template<class ScalarType, class MV, class OP,
105  const bool lapackSupportsScalarType =
108  public Details::SolverManagerRequiresLapack<ScalarType,MV,OP>
109  {
110  static const bool requiresLapack =
112  typedef Details::SolverManagerRequiresLapack<ScalarType, MV, OP,
114 
115  public:
117  base_type ()
118  {}
121  base_type ()
122  {}
123  virtual ~BlockCGSolMgr () {}
124  };
125 
126 
127  // Partial specialization for ScalarType types for which
128  // Teuchos::LAPACK has a valid implementation. This contains the
129  // actual working implementation of BlockCGSolMgr.
130  template<class ScalarType, class MV, class OP>
131  class BlockCGSolMgr<ScalarType, MV, OP, true> :
132  public Details::SolverManagerRequiresLapack<ScalarType, MV, OP, true>
133  {
134  // This partial specialization is already chosen for those scalar types
135  // that support lapack, so we don't need to have an additional compile-time
136  // check that the scalar type is float/double/complex.
137 // #if defined(HAVE_TEUCHOSCORE_CXX11)
138 // # if defined(HAVE_TEUCHOS_COMPLEX)
139 // static_assert (std::is_same<ScalarType, std::complex<float> >::value ||
140 // std::is_same<ScalarType, std::complex<double> >::value ||
141 // std::is_same<ScalarType, float>::value ||
142 // std::is_same<ScalarType, double>::value,
143 // "Belos::GCRODRSolMgr: ScalarType must be one of the four "
144 // "types (S,D,C,Z) supported by LAPACK.");
145 // # else
146 // static_assert (std::is_same<ScalarType, float>::value ||
147 // std::is_same<ScalarType, double>::value,
148 // "Belos::GCRODRSolMgr: ScalarType must be float or double. "
149 // "Complex arithmetic support is currently disabled. To "
150 // "enable it, set Teuchos_ENABLE_COMPLEX=ON.");
151 // # endif // defined(HAVE_TEUCHOS_COMPLEX)
152 // #endif // defined(HAVE_TEUCHOSCORE_CXX11)
153 
154  private:
160 
161  public:
162 
164 
165 
171  BlockCGSolMgr();
172 
211 
213  virtual ~BlockCGSolMgr() {};
214 
218  }
220 
222 
223 
224  const LinearProblem<ScalarType,MV,OP>& getProblem() const override {
225  return *problem_;
226  }
227 
230  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const override;
231 
235 
242  return Teuchos::tuple(timerSolve_);
243  }
244 
250  MagnitudeType achievedTol() const override {
251  return achievedTol_;
252  }
253 
255  int getNumIters() const override {
256  return numIters_;
257  }
258 
261  bool isLOADetected() const override { return false; }
262 
264 
266 
267 
269  void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) override { problem_ = problem; }
270 
272  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params ) override;
273 
275 
277 
278 
282  void reset( const ResetType type ) override { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
284 
286 
287 
305  ReturnType solve() override;
306 
308 
311 
313  std::string description() const override;
314 
316 
317  private:
318 
321 
326 
332 
335 
338 
341 
344 
347 
348  //
349  // Default solver parameters.
350  //
351  static constexpr int maxIters_default_ = 1000;
352  static constexpr bool adaptiveBlockSize_default_ = true;
353  static constexpr bool showMaxResNormOnly_default_ = false;
354  static constexpr bool useSingleReduction_default_ = false;
355  static constexpr int blockSize_default_ = 1;
356  static constexpr int verbosity_default_ = Belos::Errors;
357  static constexpr int outputStyle_default_ = Belos::General;
358  static constexpr int outputFreq_default_ = -1;
359  static constexpr const char * resNorm_default_ = "TwoNorm";
360  static constexpr bool foldConvergenceDetectionIntoAllreduce_default_ = false;
361  static constexpr const char * resScale_default_ = "Norm of Initial Residual";
362  static constexpr const char * label_default_ = "Belos";
363  static constexpr const char * orthoType_default_ = "ICGS";
364  static constexpr bool assertPositiveDefiniteness_default_ = true;
365 
366  //
367  // Current solver parameters and other values.
368  //
369 
372 
375 
382 
385 
388 
390  int blockSize_, verbosity_, outputStyle_, outputFreq_;
391  bool adaptiveBlockSize_, showMaxResNormOnly_, useSingleReduction_;
392  std::string orthoType_, resScale_;
395 
397  std::string label_;
398 
401 
403  bool isSet_;
404  };
405 
406 
407 // Empty Constructor
408 template<class ScalarType, class MV, class OP>
410  outputStream_(Teuchos::rcpFromRef(std::cout)),
411  convtol_(DefaultSolverParameters::convTol),
412  orthoKappa_(DefaultSolverParameters::orthoKappa),
413  achievedTol_(Teuchos::ScalarTraits<MagnitudeType>::zero()),
414  maxIters_(maxIters_default_),
415  numIters_(0),
416  blockSize_(blockSize_default_),
417  verbosity_(verbosity_default_),
418  outputStyle_(outputStyle_default_),
419  outputFreq_(outputFreq_default_),
420  adaptiveBlockSize_(adaptiveBlockSize_default_),
421  showMaxResNormOnly_(showMaxResNormOnly_default_),
422  useSingleReduction_(useSingleReduction_default_),
423  orthoType_(orthoType_default_),
424  resScale_(resScale_default_),
425  assertPositiveDefiniteness_(assertPositiveDefiniteness_default_),
426  foldConvergenceDetectionIntoAllreduce_(foldConvergenceDetectionIntoAllreduce_default_),
427  label_(label_default_),
428  isSet_(false)
429 {}
430 
431 
432 // Basic Constructor
433 template<class ScalarType, class MV, class OP>
437  problem_(problem),
438  outputStream_(Teuchos::rcpFromRef(std::cout)),
439  convtol_(DefaultSolverParameters::convTol),
440  orthoKappa_(DefaultSolverParameters::orthoKappa),
441  achievedTol_(Teuchos::ScalarTraits<MagnitudeType>::zero()),
442  maxIters_(maxIters_default_),
443  numIters_(0),
444  blockSize_(blockSize_default_),
445  verbosity_(verbosity_default_),
446  outputStyle_(outputStyle_default_),
447  outputFreq_(outputFreq_default_),
448  adaptiveBlockSize_(adaptiveBlockSize_default_),
449  showMaxResNormOnly_(showMaxResNormOnly_default_),
450  useSingleReduction_(useSingleReduction_default_),
451  orthoType_(orthoType_default_),
452  resScale_(resScale_default_),
453  assertPositiveDefiniteness_(assertPositiveDefiniteness_default_),
454  foldConvergenceDetectionIntoAllreduce_(foldConvergenceDetectionIntoAllreduce_default_),
455  label_(label_default_),
456  isSet_(false)
457 {
458  TEUCHOS_TEST_FOR_EXCEPTION(problem_.is_null(), std::invalid_argument,
459  "BlockCGSolMgr's constructor requires a nonnull LinearProblem instance.");
460 
461  // If the user passed in a nonnull parameter list, set parameters.
462  // Otherwise, the next solve() call will use default parameters,
463  // unless the user calls setParameters() first.
464  if (! pl.is_null()) {
465  setParameters (pl);
466  }
467 }
468 
469 template<class ScalarType, class MV, class OP>
470 void
473 {
474  // Create the internal parameter list if one doesn't already exist.
475  if (params_ == Teuchos::null) {
476  params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
477  }
478  else {
479  params->validateParameters(*getValidParameters());
480  }
481 
482  // Check for maximum number of iterations
483  if (params->isParameter("Maximum Iterations")) {
484  maxIters_ = params->get("Maximum Iterations",maxIters_default_);
485 
486  // Update parameter in our list and in status test.
487  params_->set("Maximum Iterations", maxIters_);
488  if (maxIterTest_!=Teuchos::null)
489  maxIterTest_->setMaxIters( maxIters_ );
490  }
491 
492  // Check for blocksize
493  if (params->isParameter("Block Size")) {
494  blockSize_ = params->get("Block Size",blockSize_default_);
495  TEUCHOS_TEST_FOR_EXCEPTION(blockSize_ <= 0, std::invalid_argument,
496  "Belos::BlockCGSolMgr: \"Block Size\" must be strictly positive.");
497 
498  // Update parameter in our list.
499  params_->set("Block Size", blockSize_);
500  }
501 
502  // Check if the blocksize should be adaptive
503  if (params->isParameter("Adaptive Block Size")) {
504  adaptiveBlockSize_ = params->get("Adaptive Block Size",adaptiveBlockSize_default_);
505 
506  // Update parameter in our list.
507  params_->set("Adaptive Block Size", adaptiveBlockSize_);
508  }
509 
510  // Check if the user is requesting the single-reduction version of CG (only for blocksize == 1)
511  if (params->isParameter("Use Single Reduction")) {
512  useSingleReduction_ = params->get("Use Single Reduction", useSingleReduction_default_);
513  }
514 
515  if (params->isParameter("Fold Convergence Detection Into Allreduce")) {
516  foldConvergenceDetectionIntoAllreduce_ = params->get("Fold Convergence Detection Into Allreduce",
517  foldConvergenceDetectionIntoAllreduce_default_);
518  }
519 
520  // Check to see if the timer label changed.
521  if (params->isParameter("Timer Label")) {
522  std::string tempLabel = params->get("Timer Label", label_default_);
523 
524  // Update parameter in our list and solver timer
525  if (tempLabel != label_) {
526  label_ = tempLabel;
527  params_->set("Timer Label", label_);
528  std::string solveLabel = label_ + ": BlockCGSolMgr total solve time";
529 #ifdef BELOS_TEUCHOS_TIME_MONITOR
530  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
531 #endif
532  if (ortho_ != Teuchos::null) {
533  ortho_->setLabel( label_ );
534  }
535  }
536  }
537 
538  // Check for a change in verbosity level
539  if (params->isParameter("Verbosity")) {
540  if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
541  verbosity_ = params->get("Verbosity", verbosity_default_);
542  } else {
543  verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
544  }
545 
546  // Update parameter in our list.
547  params_->set("Verbosity", verbosity_);
548  if (printer_ != Teuchos::null)
549  printer_->setVerbosity(verbosity_);
550  }
551 
552  // Check for a change in output style
553  if (params->isParameter("Output Style")) {
554  if (Teuchos::isParameterType<int>(*params,"Output Style")) {
555  outputStyle_ = params->get("Output Style", outputStyle_default_);
556  } else {
557  outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
558  }
559 
560  // Update parameter in our list.
561  params_->set("Output Style", outputStyle_);
562  outputTest_ = Teuchos::null;
563  }
564 
565  // output stream
566  if (params->isParameter("Output Stream")) {
567  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
568 
569  // Update parameter in our list.
570  params_->set("Output Stream", outputStream_);
571  if (printer_ != Teuchos::null)
572  printer_->setOStream( outputStream_ );
573  }
574 
575  // frequency level
576  if (verbosity_ & Belos::StatusTestDetails) {
577  if (params->isParameter("Output Frequency")) {
578  outputFreq_ = params->get("Output Frequency", outputFreq_default_);
579  }
580 
581  // Update parameter in out list and output status test.
582  params_->set("Output Frequency", outputFreq_);
583  if (outputTest_ != Teuchos::null)
584  outputTest_->setOutputFrequency( outputFreq_ );
585  }
586 
587  // Create output manager if we need to.
588  if (printer_ == Teuchos::null) {
589  printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
590  }
591 
592  // Check if the orthogonalization changed.
593  bool changedOrthoType = false;
594  if (params->isParameter("Orthogonalization")) {
595  std::string tempOrthoType = params->get("Orthogonalization",orthoType_default_);
596  if (tempOrthoType != orthoType_) {
597  orthoType_ = tempOrthoType;
598  changedOrthoType = true;
599  }
600  }
601  params_->set("Orthogonalization", orthoType_);
602 
603  // Check which orthogonalization constant to use.
604  if (params->isParameter("Orthogonalization Constant")) {
605  if (params->isType<MagnitudeType> ("Orthogonalization Constant")) {
606  orthoKappa_ = params->get ("Orthogonalization Constant",
607  static_cast<MagnitudeType> (DefaultSolverParameters::orthoKappa));
608  }
609  else {
610  orthoKappa_ = params->get ("Orthogonalization Constant",
612  }
613 
614  // Update parameter in our list.
615  params_->set("Orthogonalization Constant",orthoKappa_);
616  if (orthoType_=="DGKS") {
617  if (orthoKappa_ > 0 && ortho_ != Teuchos::null && !changedOrthoType) {
618  Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
619  }
620  }
621  }
622 
623  // Create orthogonalization manager if we need to.
624  if (ortho_ == Teuchos::null || changedOrthoType) {
627  if (orthoType_=="DGKS" && orthoKappa_ > 0) {
628  paramsOrtho = Teuchos::rcp(new Teuchos::ParameterList());
629  paramsOrtho->set ("depTol", orthoKappa_ );
630  }
631 
632  ortho_ = factory.makeMatOrthoManager (orthoType_, Teuchos::null, printer_, label_, paramsOrtho);
633  }
634 
635  // Convergence
636  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
637  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestResNorm_t;
638 
639  // Check for convergence tolerance
640  if (params->isParameter("Convergence Tolerance")) {
641  if (params->isType<MagnitudeType> ("Convergence Tolerance")) {
642  convtol_ = params->get ("Convergence Tolerance",
643  static_cast<MagnitudeType> (DefaultSolverParameters::convTol));
644  }
645  else {
646  convtol_ = params->get ("Convergence Tolerance", DefaultSolverParameters::convTol);
647  }
648 
649  // Update parameter in our list and residual tests.
650  params_->set("Convergence Tolerance", convtol_);
651  if (convTest_ != Teuchos::null)
652  convTest_->setTolerance( convtol_ );
653  }
654 
655  if (params->isParameter("Show Maximum Residual Norm Only")) {
656  showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");
657 
658  // Update parameter in our list and residual tests
659  params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
660  if (convTest_ != Teuchos::null)
661  convTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
662  }
663 
664  // Check for a change in scaling, if so we need to build new residual tests.
665  bool newResTest = false;
666  {
667  std::string tempResScale = resScale_;
668  if (params->isParameter ("Implicit Residual Scaling")) {
669  tempResScale = params->get<std::string> ("Implicit Residual Scaling");
670  }
671 
672  // Only update the scaling if it's different.
673  if (resScale_ != tempResScale) {
674  const Belos::ScaleType resScaleType =
675  convertStringToScaleType (tempResScale);
676  resScale_ = tempResScale;
677 
678  // Update parameter in our list and residual tests
679  params_->set ("Implicit Residual Scaling", resScale_);
680 
681  if (! convTest_.is_null ()) {
682  try {
683  NormType normType = Belos::TwoNorm;
684  if (params->isParameter("Residual Norm")) {
685  if (params->isType<std::string> ("Residual Norm")) {
686  normType = convertStringToNormType(params->get<std::string> ("Residual Norm"));
687  }
688  }
689  convTest_->defineResForm(StatusTestResNorm_t::Implicit, normType);
690  convTest_->defineScaleForm (resScaleType, Belos::TwoNorm);
691  }
692  catch (std::exception& e) {
693  // Make sure the convergence test gets constructed again.
694  newResTest = true;
695  }
696  }
697  }
698  }
699 
700  // Create status tests if we need to.
701 
702  // Basic test checks maximum iterations and native residual.
703  if (maxIterTest_ == Teuchos::null)
704  maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
705 
706  // Implicit residual test, using the native residual to determine if convergence was achieved.
707  if (convTest_.is_null () || newResTest) {
708 
709  NormType normType = Belos::TwoNorm;
710  if (params->isParameter("Residual Norm")) {
711  if (params->isType<std::string> ("Residual Norm")) {
712  normType = convertStringToNormType(params->get<std::string> ("Residual Norm"));
713  }
714  }
715 
716  convTest_ = rcp (new StatusTestResNorm_t (convtol_, 1, showMaxResNormOnly_));
717  convTest_->defineResForm(StatusTestResNorm_t::Implicit, normType);
718  convTest_->defineScaleForm (convertStringToScaleType (resScale_), Belos::TwoNorm);
719  }
720 
721  if (sTest_.is_null () || newResTest)
722  sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
723 
724  if (outputTest_.is_null () || newResTest) {
725 
726  // Create the status test output class.
727  // This class manages and formats the output from the status test.
728  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
729  outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
730 
731  // Set the solver string for the output test
732  std::string solverDesc = " Block CG ";
733  outputTest_->setSolverDesc( solverDesc );
734 
735  }
736 
737  // BelosCgIter accepts a parameter specifying whether to assert for the positivity of p^H*A*p in the CG iteration
738  if (params->isParameter("Assert Positive Definiteness")) {
739  assertPositiveDefiniteness_ = Teuchos::getParameter<bool>(*params,"Assert Positive Definiteness");
740  params_->set("Assert Positive Definiteness", assertPositiveDefiniteness_);
741  }
742 
743  // Create the timer if we need to.
744  if (timerSolve_ == Teuchos::null) {
745  std::string solveLabel = label_ + ": BlockCGSolMgr total solve time";
746 #ifdef BELOS_TEUCHOS_TIME_MONITOR
747  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
748 #endif
749  }
750 
751  // Inform the solver manager that the current parameters were set.
752  isSet_ = true;
753 }
754 
755 
756 template<class ScalarType, class MV, class OP>
759 {
761 
762  // Set all the valid parameters and their default values.
763  if(is_null(validPL)) {
764  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
765  pl->set("Convergence Tolerance", static_cast<MagnitudeType>(DefaultSolverParameters::convTol),
766  "The relative residual tolerance that needs to be achieved by the\n"
767  "iterative solver in order for the linear system to be declared converged.");
768  pl->set("Maximum Iterations", static_cast<int>(maxIters_default_),
769  "The maximum number of block iterations allowed for each\n"
770  "set of RHS solved.");
771  pl->set("Block Size", static_cast<int>(blockSize_default_),
772  "The number of vectors in each block.");
773  pl->set("Adaptive Block Size", static_cast<bool>(adaptiveBlockSize_default_),
774  "Whether the solver manager should adapt to the block size\n"
775  "based on the number of RHS to solve.");
776  pl->set("Verbosity", static_cast<int>(verbosity_default_),
777  "What type(s) of solver information should be outputted\n"
778  "to the output stream.");
779  pl->set("Output Style", static_cast<int>(outputStyle_default_),
780  "What style is used for the solver information outputted\n"
781  "to the output stream.");
782  pl->set("Output Frequency", static_cast<int>(outputFreq_default_),
783  "How often convergence information should be outputted\n"
784  "to the output stream.");
785  pl->set("Output Stream", Teuchos::rcpFromRef(std::cout),
786  "A reference-counted pointer to the output stream where all\n"
787  "solver output is sent.");
788  pl->set("Show Maximum Residual Norm Only", static_cast<bool>(showMaxResNormOnly_default_),
789  "When convergence information is printed, only show the maximum\n"
790  "relative residual norm when the block size is greater than one.");
791  pl->set("Use Single Reduction", static_cast<bool>(useSingleReduction_default_),
792  "Use single reduction iteration when the block size is one.");
793  pl->set("Implicit Residual Scaling", resScale_default_,
794  "The type of scaling used in the residual convergence test.");
795  pl->set("Timer Label", static_cast<const char *>(label_default_),
796  "The string to use as a prefix for the timer labels.");
797  pl->set("Orthogonalization", static_cast<const char *>(orthoType_default_),
798  "The type of orthogonalization to use: DGKS, ICGS, or IMGS.");
799  pl->set("Assert Positive Definiteness",static_cast<bool>(assertPositiveDefiniteness_default_),
800  "Assert for positivity of p^H*A*p in CG iteration.");
801  pl->set("Orthogonalization Constant",static_cast<MagnitudeType>(DefaultSolverParameters::orthoKappa),
802  "The constant used by DGKS orthogonalization to determine\n"
803  "whether another step of classical Gram-Schmidt is necessary.");
804  pl->set("Residual Norm",static_cast<const char *>(resNorm_default_),
805  "Norm used for the convergence check on the residual.");
806  pl->set("Fold Convergence Detection Into Allreduce",static_cast<bool>(foldConvergenceDetectionIntoAllreduce_default_),
807  "Merge the allreduce for convergence detection with the one for CG.\n"
808  "This saves one all-reduce, but incurs more computation.");
809  validPL = pl;
810  }
811  return validPL;
812 }
813 
814 
815 // solve()
816 template<class ScalarType, class MV, class OP>
818  using Teuchos::RCP;
819  using Teuchos::rcp;
820  using Teuchos::rcp_const_cast;
821  using Teuchos::rcp_dynamic_cast;
822 
823  // Set the current parameters if they were not set before. NOTE:
824  // This may occur if the user generated the solver manager with the
825  // default constructor and then didn't set any parameters using
826  // setParameters().
827  if (!isSet_) {
828  setParameters(Teuchos::parameterList(*getValidParameters()));
829  }
830 
832 
833  TEUCHOS_TEST_FOR_EXCEPTION( !problem_->isProblemSet(),
835  "Belos::BlockCGSolMgr::solve(): Linear problem is not ready, setProblem() "
836  "has not been called.");
837 
838  // Create indices for the linear systems to be solved.
839  int startPtr = 0;
840  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
841  int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
842 
843  std::vector<int> currIdx, currIdx2;
844  // If an adaptive block size is allowed then only the linear
845  // systems that need to be solved are solved. Otherwise, the index
846  // set is generated that informs the linear problem that some
847  // linear systems are augmented.
848  if ( adaptiveBlockSize_ ) {
849  blockSize_ = numCurrRHS;
850  currIdx.resize( numCurrRHS );
851  currIdx2.resize( numCurrRHS );
852  for (int i=0; i<numCurrRHS; ++i)
853  { currIdx[i] = startPtr+i; currIdx2[i]=i; }
854 
855  }
856  else {
857  currIdx.resize( blockSize_ );
858  currIdx2.resize( blockSize_ );
859  for (int i=0; i<numCurrRHS; ++i)
860  { currIdx[i] = startPtr+i; currIdx2[i]=i; }
861  for (int i=numCurrRHS; i<blockSize_; ++i)
862  { currIdx[i] = -1; currIdx2[i] = i; }
863  }
864 
865  // Inform the linear problem of the current linear system to solve.
866  problem_->setLSIndex( currIdx );
867 
869  // Set up the parameter list for the Iteration subclass.
871  plist.set("Block Size",blockSize_);
872 
873  // Reset the output status test (controls all the other status tests).
874  outputTest_->reset();
875 
876  // Assume convergence is achieved, then let any failed convergence
877  // set this to false. "Innocent until proven guilty."
878  bool isConverged = true;
879 
881  // Set up the BlockCG Iteration subclass.
882 
883  plist.set("Assert Positive Definiteness", assertPositiveDefiniteness_);
884 
885  RCP<CGIteration<ScalarType,MV,OP> > block_cg_iter;
886  if (blockSize_ == 1) {
887  // Standard (nonblock) CG is faster for the special case of a
888  // block size of 1. A single reduction iteration can also be used
889  // if collectives are more expensive than vector updates.
890  plist.set("Fold Convergence Detection Into Allreduce",
891  foldConvergenceDetectionIntoAllreduce_);
892  if (useSingleReduction_) {
893  block_cg_iter =
894  rcp (new CGSingleRedIter<ScalarType,MV,OP> (problem_, printer_,
895  outputTest_, convTest_, plist));
896  }
897  else {
898  block_cg_iter =
899  rcp (new CGIter<ScalarType,MV,OP> (problem_, printer_,
900  outputTest_, convTest_, plist));
901  }
902  } else {
903  block_cg_iter =
904  rcp (new BlockCGIter<ScalarType,MV,OP> (problem_, printer_, outputTest_,
905  ortho_, plist));
906  }
907 
908 
909  // Enter solve() iterations
910  {
911 #ifdef BELOS_TEUCHOS_TIME_MONITOR
912  Teuchos::TimeMonitor slvtimer(*timerSolve_);
913 #endif
914 
915  while ( numRHS2Solve > 0 ) {
916  //
917  // Reset the active / converged vectors from this block
918  std::vector<int> convRHSIdx;
919  std::vector<int> currRHSIdx( currIdx );
920  currRHSIdx.resize(numCurrRHS);
921 
922  // Reset the number of iterations.
923  block_cg_iter->resetNumIters();
924 
925  // Reset the number of calls that the status test output knows about.
926  outputTest_->resetNumCalls();
927 
928  // Get the current residual for this block of linear systems.
929  RCP<MV> R_0 = MVT::CloneViewNonConst( *(rcp_const_cast<MV>(problem_->getInitResVec())), currIdx );
930 
931  // Set the new state and initialize the solver.
933  newstate.R = R_0;
934  block_cg_iter->initializeCG(newstate);
935 
936  while(1) {
937 
938  // tell block_cg_iter to iterate
939  try {
940  block_cg_iter->iterate();
941  //
942  // Check whether any of the linear systems converged.
943  //
944  if (convTest_->getStatus() == Passed) {
945  // At least one of the linear system(s) converged.
946  //
947  // Get the column indices of the linear systems that converged.
948  typedef StatusTestGenResNorm<ScalarType,MV,OP> conv_test_type;
949  std::vector<int> convIdx =
950  rcp_dynamic_cast<conv_test_type>(convTest_)->convIndices();
951 
952  // If the number of converged linear systems equals the
953  // number of linear systems currently being solved, then
954  // we are done with this block.
955  if (convIdx.size() == currRHSIdx.size())
956  break; // break from while(1){block_cg_iter->iterate()}
957 
958  // Inform the linear problem that we are finished with
959  // this current linear system.
960  problem_->setCurrLS();
961 
962  // Reset currRHSIdx to contain the right-hand sides that
963  // are left to converge for this block.
964  int have = 0;
965  std::vector<int> unconvIdx(currRHSIdx.size());
966  for (unsigned int i=0; i<currRHSIdx.size(); ++i) {
967  bool found = false;
968  for (unsigned int j=0; j<convIdx.size(); ++j) {
969  if (currRHSIdx[i] == convIdx[j]) {
970  found = true;
971  break;
972  }
973  }
974  if (!found) {
975  currIdx2[have] = currIdx2[i];
976  currRHSIdx[have++] = currRHSIdx[i];
977  }
978  else {
979  }
980  }
981  currRHSIdx.resize(have);
982  currIdx2.resize(have);
983 
984  // Set the remaining indices after deflation.
985  problem_->setLSIndex( currRHSIdx );
986 
987  // Get the current residual vector.
988  std::vector<MagnitudeType> norms;
989  R_0 = MVT::CloneCopy( *(block_cg_iter->getNativeResiduals(&norms)),currIdx2 );
990  for (int i=0; i<have; ++i) { currIdx2[i] = i; }
991 
992  // Set the new blocksize for the solver.
993  block_cg_iter->setBlockSize( have );
994 
995  // Set the new state and initialize the solver.
997  defstate.R = R_0;
998  block_cg_iter->initializeCG(defstate);
999  }
1000  //
1001  // None of the linear systems converged. Check whether the
1002  // maximum iteration count was reached.
1003  //
1004  else if (maxIterTest_->getStatus() == Passed) {
1005  isConverged = false; // None of the linear systems converged.
1006  break; // break from while(1){block_cg_iter->iterate()}
1007  }
1008  //
1009  // iterate() returned, but none of our status tests Passed.
1010  // This indicates a bug.
1011  //
1012  else {
1013  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
1014  "Belos::BlockCGSolMgr::solve(): Neither the convergence test nor "
1015  "the maximum iteration count test passed. Please report this bug "
1016  "to the Belos developers.");
1017  }
1018  }
1019  catch (const StatusTestNaNError& e) {
1020  // A NaN was detected in the solver. Set the solution to zero and return unconverged.
1021  achievedTol_ = MT::one();
1022  Teuchos::RCP<MV> X = problem_->getLHS();
1023  MVT::MvInit( *X, SCT::zero() );
1024  printer_->stream(Warnings) << "Belos::BlockCGSolMgr::solve(): Warning! NaN has been detected!"
1025  << std::endl;
1026  return Unconverged;
1027  }
1028  catch (const std::exception &e) {
1029  std::ostream& err = printer_->stream (Errors);
1030  err << "Error! Caught std::exception in CGIteration::iterate() at "
1031  << "iteration " << block_cg_iter->getNumIters() << std::endl
1032  << e.what() << std::endl;
1033  throw;
1034  }
1035  }
1036 
1037  // Inform the linear problem that we are finished with this
1038  // block linear system.
1039  problem_->setCurrLS();
1040 
1041  // Update indices for the linear systems to be solved.
1042  startPtr += numCurrRHS;
1043  numRHS2Solve -= numCurrRHS;
1044  if ( numRHS2Solve > 0 ) {
1045  numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1046 
1047  if ( adaptiveBlockSize_ ) {
1048  blockSize_ = numCurrRHS;
1049  currIdx.resize( numCurrRHS );
1050  currIdx2.resize( numCurrRHS );
1051  for (int i=0; i<numCurrRHS; ++i)
1052  { currIdx[i] = startPtr+i; currIdx2[i] = i; }
1053  }
1054  else {
1055  currIdx.resize( blockSize_ );
1056  currIdx2.resize( blockSize_ );
1057  for (int i=0; i<numCurrRHS; ++i)
1058  { currIdx[i] = startPtr+i; currIdx2[i] = i; }
1059  for (int i=numCurrRHS; i<blockSize_; ++i)
1060  { currIdx[i] = -1; currIdx2[i] = i; }
1061  }
1062  // Set the next indices.
1063  problem_->setLSIndex( currIdx );
1064 
1065  // Set the new blocksize for the solver.
1066  block_cg_iter->setBlockSize( blockSize_ );
1067  }
1068  else {
1069  currIdx.resize( numRHS2Solve );
1070  }
1071 
1072  }// while ( numRHS2Solve > 0 )
1073 
1074  }
1075 
1076  // print final summary
1077  sTest_->print( printer_->stream(FinalSummary) );
1078 
1079  // print timing information
1080 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1081  // Calling summarize() requires communication in general, so don't
1082  // call it unless the user wants to print out timing details.
1083  // summarize() will do all the work even if it's passed a "black
1084  // hole" output stream.
1085  if (verbosity_ & TimingDetails) {
1086  Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
1087  }
1088 #endif
1089 
1090  // Save the iteration count for this solve.
1091  numIters_ = maxIterTest_->getNumIters();
1092 
1093  // Save the convergence test value ("achieved tolerance") for this solve.
1094  {
1095  typedef StatusTestGenResNorm<ScalarType,MV,OP> conv_test_type;
1096  // testValues is nonnull and not persistent.
1097  const std::vector<MagnitudeType>* pTestValues =
1098  rcp_dynamic_cast<conv_test_type>(convTest_)->getTestValue();
1099 
1100  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
1101  "Belos::BlockCGSolMgr::solve(): The convergence test's getTestValue() "
1102  "method returned NULL. Please report this bug to the Belos developers.");
1103 
1104  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
1105  "Belos::BlockCGSolMgr::solve(): The convergence test's getTestValue() "
1106  "method returned a vector of length zero. Please report this bug to the "
1107  "Belos developers.");
1108 
1109  // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
1110  // achieved tolerances for all vectors in the current solve(), or
1111  // just for the vectors from the last deflation?
1112  achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1113  }
1114 
1115  if (!isConverged) {
1116  return Unconverged; // return from BlockCGSolMgr::solve()
1117  }
1118  return Converged; // return from BlockCGSolMgr::solve()
1119 }
1120 
1121 // This method requires the solver manager to return a std::string that describes itself.
1122 template<class ScalarType, class MV, class OP>
1124 {
1125  std::ostringstream oss;
1126  oss << "Belos::BlockCGSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
1127  oss << "{";
1128  oss << "Ortho Type='"<<orthoType_<<"\', Block Size=" << blockSize_;
1129  oss << "}";
1130  return oss.str();
1131 }
1132 
1133 } // end Belos namespace
1134 
1135 #endif /* BELOS_BLOCK_CG_SOLMGR_HPP */
static const double orthoKappa
DGKS orthogonalization constant.
Definition: BelosTypes.hpp:299
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Definition: BelosTypes.cpp:106
Teuchos::RCP< const MV > R
The current residual.
bool isLOADetected() const override
Return whether a loss of accuracy was detected by this solver during the most current solve...
Collection of types and exceptions used within the Belos solvers.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
int maxIters_
Maximum iteration count (read from parameter list).
Class which manages the output and verbosity of the Belos solvers.
bool is_null(const boost::shared_ptr< T > &p)
bool isSet_
Whether or not the parameters have been set (via setParameters()).
Teuchos::RCP< Teuchos::ParameterList > params_
Current parameter list.
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > outputTest_
Output &quot;status test&quot; that controls all the other status tests.
ScaleType
The type of scaling to use on the residual norm value.
Definition: BelosTypes.hpp:120
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
This class implements the preconditioned Conjugate Gradient (CG) iteration.
Definition: BelosCGIter.hpp:79
T & get(ParameterList &l, const std::string &name)
Structure to contain pointers to CGIteration state variables.
Belos concrete class for performing the conjugate-gradient (CG) iteration.
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)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
A factory class for generating StatusTestOutput objects.
Belos concrete class for performing the block conjugate-gradient (CG) iteration.
An implementation of StatusTestResNorm using a family of residual norms.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
std::string label_
Prefix label for all the timers.
The Belos::BlockCGSolMgr provides a powerful and fully-featured solver manager over the CG and BlockC...
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()
BlockCGSolMgrLinearProblemFailure(const std::string &what_arg)
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
BlockCGSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i.e.
A factory class for generating StatusTestOutput objects.
Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > ortho_
Orthogonalization manager.
BlockCGSolMgr(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl)
Teuchos::RCP< StatusTestMaxIters< ScalarType, MV, OP > > maxIterTest_
Maximum iteration count stopping criterion.
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > sTest_
Aggregate stopping criterion.
Details::SolverManagerRequiresLapack< ScalarType, MV, OP, requiresLapack > base_type
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
void reset(const ResetType type) override
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
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. ...
Belos concrete class for performing a single-reduction conjugate-gradient (CG) iteration.
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)
bool is_null(const RCP< T > &p)
NormType convertStringToNormType(const std::string &normType)
Convert the given string to its NormType enum value.
Definition: BelosTypes.cpp:88
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
MagnitudeType achievedTol_
Tolerance achieved by the last solve() invocation.
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
Type traits class that says whether Teuchos::LAPACK has a valid implementation for the given ScalarTy...
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:155
Teuchos::RCP< StatusTestGenResNorm< ScalarType, MV, OP > > convTest_
Convergence stopping criterion.
int getNumIters() const override
Get the iteration count for the most recent call to solve().
Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > problem_
The linear problem to solve.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
NormType
The type of vector norm to compute.
Definition: BelosTypes.hpp:97
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.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
MagnitudeType convtol_
Convergence tolerance (read from parameter list).
bool isType(const std::string &name) const
Teuchos::RCP< Teuchos::Time > timerSolve_
Solve timer.
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
A class for extending the status testing capabilities of Belos via logical combinations.
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Return a reference to the linear problem being solved by this solver manager.
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
This class implements the preconditioned single-reduction Conjugate Gradient (CG) iteration...
Teuchos::RCP< OutputManager< ScalarType > > printer_
Output manager, that handles printing of different kinds of messages.
Base class for Belos::SolverManager subclasses which normally can only compile with ScalarType types ...
Belos header file which uses auto-configuration information to include necessary C++ headers...
int numIters_
Number of iterations taken by the last solve() invocation.
Teuchos::ScalarTraits< MagnitudeType > MT
static const bool requiresLapack
Teuchos::RCP< std::ostream > outputStream_
Output stream to which the output manager prints.
MagnitudeType orthoKappa_
Orthogonalization parameter (read from parameter list).
Stub implementation of BlockCGIter, for ScalarType types for which Teuchos::LAPACK does NOT have a va...
Teuchos::RCP< Belos::MatOrthoManager< Scalar, MV, OP > > makeMatOrthoManager(const std::string &ortho, const Teuchos::RCP< const OP > &M, const Teuchos::RCP< OutputManager< Scalar > > &, const std::string &label, const Teuchos::RCP< Teuchos::ParameterList > &params)
Return an instance of the specified MatOrthoManager subclass.