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

Generated on Mon Feb 5 2018 15:01:51 for Belos by doxygen 1.8.13