Stratimikos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Thyra_BelosLinearOpWithSolveFactory_def.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Stratimikos: Thyra-based strategies for linear solvers
4 //
5 // Copyright 2006 NTESS and the Stratimikos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP
11 #define THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP
12 
13 
14 #include "Thyra_BelosLinearOpWithSolveFactory_decl.hpp"
15 #include "Thyra_BelosLinearOpWithSolve.hpp"
16 #include "Thyra_ScaledAdjointLinearOpBase.hpp"
17 
20 #include "BelosBlockCGSolMgr.hpp"
23 #include "BelosGCRODRSolMgr.hpp"
24 #include "BelosRCGSolMgr.hpp"
25 #include "BelosMinresSolMgr.hpp"
26 #include "BelosTFQMRSolMgr.hpp"
27 #include "BelosBiCGStabSolMgr.hpp"
29 #include "BelosThyraAdapter.hpp"
30 
31 #if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
32 #include "Thyra_BelosTpetrasSolverAdapter.hpp"
33 #endif
34 
35 #include "Teuchos_VerboseObjectParameterListHelpers.hpp"
36 #include "Teuchos_StandardParameterEntryValidators.hpp"
37 #include "Teuchos_ParameterList.hpp"
38 #include "Teuchos_dyn_cast.hpp"
39 #include "Teuchos_ValidatorXMLConverterDB.hpp"
40 #include "Teuchos_StandardValidatorXMLConverters.hpp"
41 
42 
43 namespace Thyra {
44 
45 
46 // Parameter names for Parameter List
47 
48 template<class Scalar>
49 const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverType_name = "Solver Type";
50 template<class Scalar>
51 const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverType_default = "Pseudo Block GMRES";
52 template<class Scalar>
53 const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverTypes_name = "Solver Types";
54 template<class Scalar>
55 const std::string BelosLinearOpWithSolveFactory<Scalar>::BlockGMRES_name = "Block GMRES";
56 template<class Scalar>
57 const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockGMRES_name = "Pseudo Block GMRES";
58 template<class Scalar>
59 const std::string BelosLinearOpWithSolveFactory<Scalar>::BlockCG_name = "Block CG";
60 template<class Scalar>
61 const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockCG_name = "Pseudo Block CG";
62 template<class Scalar>
63 const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockStochasticCG_name = "Pseudo Block Stochastic CG";
64 template<class Scalar>
65 const std::string BelosLinearOpWithSolveFactory<Scalar>::GCRODR_name = "GCRODR";
66 template<class Scalar>
67 const std::string BelosLinearOpWithSolveFactory<Scalar>::RCG_name = "RCG";
68 template<class Scalar>
69 const std::string BelosLinearOpWithSolveFactory<Scalar>::MINRES_name = "MINRES";
70 template<class Scalar>
71 const std::string BelosLinearOpWithSolveFactory<Scalar>::TFQMR_name = "TFQMR";
72 template<class Scalar>
73 const std::string BelosLinearOpWithSolveFactory<Scalar>::BiCGStab_name = "BiCGStab";
74 template<class Scalar>
75 const std::string BelosLinearOpWithSolveFactory<Scalar>::FixedPoint_name = "Fixed Point";
76 
77 #if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
78 template<class Scalar>
79 const std::string BelosLinearOpWithSolveFactory<Scalar>::TpetraGmres_name = "TPETRA GMRES";
80 template<class Scalar>
81 const std::string BelosLinearOpWithSolveFactory<Scalar>::TpetraGmresPipeline_name = "TPETRA GMRES PIPELINE";
82 template<class Scalar>
83 const std::string BelosLinearOpWithSolveFactory<Scalar>::TpetraGmresSingleReduce_name = "TPETRA GMRES SINGLE REDUCE";
84 template<class Scalar>
85 const std::string BelosLinearOpWithSolveFactory<Scalar>::TpetraGmresSstep_name = "TPETRA GMRES S-STEP";
86 #endif
87 
88 template<class Scalar>
89 const std::string BelosLinearOpWithSolveFactory<Scalar>::ConvergenceTestFrequency_name = "Convergence Test Frequency";
90 
91 namespace {
92 const std::string LeftPreconditionerIfUnspecified_name = "Left Preconditioner If Unspecified";
93 }
94 
95 // Constructors/initializers/accessors
96 
97 
98 template<class Scalar>
100  :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES),
101  convergenceTestFrequency_(1)
102 {
103  updateThisValidParamList();
104 }
105 
106 
107 template<class Scalar>
109  const RCP<PreconditionerFactoryBase<Scalar> > &precFactory
110  )
111  :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES)
112 {
113  this->setPreconditionerFactory(precFactory, "");
114 }
115 
116 
117 // Overridden from LinearOpWithSolveFactoryBase
118 
119 
120 template<class Scalar>
122 {
123  return true;
124 }
125 
126 
127 template<class Scalar>
129  const RCP<PreconditionerFactoryBase<Scalar> > &precFactory,
130  const std::string &precFactoryName
131  )
132 {
133  TEUCHOS_TEST_FOR_EXCEPT(!precFactory.get());
135  precFactoryValidPL = precFactory->getValidParameters();
136  const std::string _precFactoryName =
137  ( precFactoryName != ""
138  ? precFactoryName
139  : ( precFactoryValidPL.get() ? precFactoryValidPL->name() : "GENERIC PRECONDITIONER FACTORY" )
140  );
141  precFactory_ = precFactory;
142  precFactoryName_ = _precFactoryName;
143  updateThisValidParamList();
144 }
145 
146 
147 template<class Scalar>
150 {
151  return precFactory_;
152 }
153 
154 
155 template<class Scalar>
157  RCP<PreconditionerFactoryBase<Scalar> > *precFactory,
158  std::string *precFactoryName
159  )
160 {
161  if(precFactory) *precFactory = precFactory_;
162  if(precFactoryName) *precFactoryName = precFactoryName_;
163  precFactory_ = Teuchos::null;
164  precFactoryName_ = "";
165  updateThisValidParamList();
166 }
167 
168 
169 template<class Scalar>
171  const LinearOpSourceBase<Scalar> &fwdOpSrc
172  ) const
173 {
174  if(precFactory_.get())
175  return precFactory_->isCompatible(fwdOpSrc);
176  return true; // Without a preconditioner, we are compatible with all linear operators!
177 }
178 
179 
180 template<class Scalar>
183 {
185 }
186 
187 
188 template<class Scalar>
190  const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
192  const ESupportSolveUse supportSolveUse
193  ) const
194 {
195  using Teuchos::null;
196  initializeOpImpl(fwdOpSrc,null,null,false,Op,supportSolveUse);
197 }
198 
199 
200 template<class Scalar>
202  const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
204  ) const
205 {
206  using Teuchos::null;
207  initializeOpImpl(fwdOpSrc,null,null,true,Op,SUPPORT_SOLVE_UNSPECIFIED);
208 }
209 
210 
211 template<class Scalar>
213  const EPreconditionerInputType precOpType
214  ) const
215 {
216  if(precFactory_.get())
217  return true;
218  return (precOpType==PRECONDITIONER_INPUT_TYPE_AS_OPERATOR);
219 }
220 
221 
222 template<class Scalar>
224  const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
225  const RCP<const PreconditionerBase<Scalar> > &prec,
227  const ESupportSolveUse supportSolveUse
228  ) const
229 {
230  using Teuchos::null;
231  initializeOpImpl(fwdOpSrc,null,prec,false,Op,supportSolveUse);
232 }
233 
234 
235 template<class Scalar>
237  const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
238  const RCP<const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
240  const ESupportSolveUse supportSolveUse
241  ) const
242 {
243  using Teuchos::null;
244  initializeOpImpl(fwdOpSrc,approxFwdOpSrc,null,false,Op,supportSolveUse);
245 }
246 
247 
248 template<class Scalar>
251  RCP<const LinearOpSourceBase<Scalar> > *fwdOpSrc,
252  RCP<const PreconditionerBase<Scalar> > *prec,
253  RCP<const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
254  ESupportSolveUse *supportSolveUse
255  ) const
256 {
257 #ifdef TEUCHOS_DEBUG
258  TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
259 #endif
263  _fwdOpSrc = belosOp.extract_fwdOpSrc();
265  _prec = ( belosOp.isExternalPrec() ? belosOp.extract_prec() : Teuchos::null );
266  // Note: above we only extract the preconditioner if it was passed in
267  // externally. Otherwise, we need to hold on to it so that we can reuse it
268  // in the next initialization.
270  _approxFwdOpSrc = belosOp.extract_approxFwdOpSrc();
272  _supportSolveUse = belosOp.supportSolveUse();
273  if(fwdOpSrc) *fwdOpSrc = _fwdOpSrc;
274  if(prec) *prec = _prec;
275  if(approxFwdOpSrc) *approxFwdOpSrc = _approxFwdOpSrc;
276  if(supportSolveUse) *supportSolveUse = _supportSolveUse;
277 }
278 
279 
280 // Overridden from ParameterListAcceptor
281 
282 
283 template<class Scalar>
285  RCP<Teuchos::ParameterList> const& paramList
286  )
287 {
288  TEUCHOS_TEST_FOR_EXCEPT(paramList.get()==NULL);
289  paramList->validateParametersAndSetDefaults(*this->getValidParameters(), 1);
290  paramList_ = paramList;
291  solverType_ =
292  Teuchos::getIntegralValue<EBelosSolverType>(*paramList_, SolverType_name);
293  convergenceTestFrequency_ =
294  Teuchos::getParameter<int>(*paramList_, ConvergenceTestFrequency_name);
295  Teuchos::readVerboseObjectSublist(&*paramList_,this);
296 }
297 
298 
299 template<class Scalar>
302 {
303  return paramList_;
304 }
305 
306 
307 template<class Scalar>
310 {
311  RCP<Teuchos::ParameterList> _paramList = paramList_;
312  paramList_ = Teuchos::null;
313  return _paramList;
314 }
315 
316 
317 template<class Scalar>
320 {
321  return paramList_;
322 }
323 
324 
325 template<class Scalar>
328 {
329  return thisValidParamList_;
330 }
331 
332 
333 // Public functions overridden from Teuchos::Describable
334 
335 
336 template<class Scalar>
338 {
339  std::ostringstream oss;
341  oss << "{";
342  if (nonnull(precFactory_)) {
343  oss << precFactory_->description();
344  }
345  oss << "}";
346  return oss.str();
347 }
348 
349 
350 // private
351 
352 
353 template<class Scalar>
356 {
357  using Teuchos::as;
358  using Teuchos::tuple;
359  using Teuchos::setStringToIntegralParameter;
363  >::getDummyObject(),
365  EBelosSolverType> >::getDummyObject());
366 
367  typedef MultiVectorBase<Scalar> MV_t;
368  typedef LinearOpBase<Scalar> LO_t;
369  static RCP<Teuchos::ParameterList> validParamList;
370  if(validParamList.get()==NULL) {
371  validParamList = Teuchos::rcp(new Teuchos::ParameterList("BelosLinearOpWithSolveFactory"));
372  setStringToIntegralParameter<EBelosSolverType>(
373  SolverType_name, SolverType_default,
374  "Type of linear solver algorithm to use.",
375  tuple<std::string>(
376  "Block GMRES",
377  "Pseudo Block GMRES",
378  "Block CG",
379  "Pseudo Block CG",
380  "Pseudo Block Stochastic CG",
381  "GCRODR",
382  "RCG",
383  "MINRES",
384  "TFQMR",
385  "BiCGStab",
386  "Fixed Point"
387 #if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
388  ,"TPETRA GMRES",
389  "TPETRA GMRES PIPELINE",
390  "TPETRA GMRES SINGLE REDUCE",
391  "TPETRA GMRES S-STEP"
392 #endif
393  ),
394  tuple<std::string>(
395  "Block GMRES solver for nonsymmetric linear systems. It can also solve "
396  "single right-hand side systems, and can also perform Flexible GMRES "
397  "(where the preconditioner may change at every iteration, for example "
398  "for inner-outer iterations), by setting options in the \"Block GMRES\" "
399  "sublist.",
400 
401  "GMRES solver for nonsymmetric linear systems, that performs single "
402  "right-hand side solves on multiple right-hand sides at once. It "
403  "exploits operator multivector multiplication in order to amortize "
404  "global communication costs. Individual linear systems are deflated "
405  "out as they are solved.",
406 
407  "Block CG solver for symmetric (Hermitian in complex arithmetic) "
408  "positive definite linear systems. It can also solve single "
409  "right-hand-side systems.",
410 
411  "CG solver that performs single right-hand side CG on multiple right-hand "
412  "sides at once. It exploits operator multivector multiplication in order "
413  "to amortize global communication costs. Individual linear systems are "
414  "deflated out as they are solved.",
415 
416  "stochastic CG solver that performs single right-hand side CG on multiple right-hand "
417  "sides at once. It exploits operator multivector multiplication in order "
418  "to amortize global communication costs. Individual linear systems are "
419  "deflated out as they are solved. [EXPERIMENTAL]",
420 
421  "Variant of GMRES that performs subspace recycling to accelerate "
422  "convergence for sequences of solves with related linear systems. "
423  "Individual linear systems are deflated out as they are solved. "
424  "The current implementation only supports real-valued Scalar types.",
425 
426  "CG solver for symmetric (Hermitian in complex arithmetic) positive "
427  "definite linear systems, that performs subspace recycling to "
428  "accelerate convergence for sequences of related linear systems.",
429 
430  "MINRES solver for symmetric indefinite linear systems. It performs "
431  "single-right-hand-side solves on multiple right-hand sides sequentially.",
432 
433  "TFQMR (Transpose-Free QMR) solver for nonsymmetric linear systems.",
434 
435  "BiCGStab solver for nonsymmetric linear systems.",
436 
437  "Fixed point iteration"
438 
439 #if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
440  ,"Native Tpetra implementation of GMRES",
441 
442  "Native Tpetra implementation of pipeline GMRES",
443 
444  "Native Tpetra implementation of single-reduce GMRES",
445 
446  "Native Tpetra implementation of s-step GMRES"
447 #endif
448  ),
449  tuple<EBelosSolverType>(
450  SOLVER_TYPE_BLOCK_GMRES,
451  SOLVER_TYPE_PSEUDO_BLOCK_GMRES,
452  SOLVER_TYPE_BLOCK_CG,
453  SOLVER_TYPE_PSEUDO_BLOCK_CG,
454  SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG,
455  SOLVER_TYPE_GCRODR,
456  SOLVER_TYPE_RCG,
457  SOLVER_TYPE_MINRES,
458  SOLVER_TYPE_TFQMR,
459  SOLVER_TYPE_BICGSTAB,
460  SOLVER_TYPE_FIXEDPOINT
461 #if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
462  ,SOLVER_TYPE_TPETRA_GMRES,
463  SOLVER_TYPE_TPETRA_GMRES_PIPELINE,
464  SOLVER_TYPE_TPETRA_GMRES_SINGLE_REDUCE,
465  SOLVER_TYPE_TPETRA_GMRES_SSTEP
466 #endif
467  ),
468  &*validParamList
469  );
470  validParamList->set(ConvergenceTestFrequency_name, as<int>(1),
471  "Number of linear solver iterations to skip between applying"
472  " user-defined convergence test.");
473  validParamList->set(
474  LeftPreconditionerIfUnspecified_name, false,
475  "If the preconditioner does not specify if it is left or right, and this\n"
476  "option is set to true, put the preconditioner on the left side.\n"
477  "Historically, preconditioning is on the right. Some solvers may not\n"
478  "support left preconditioning.");
480  &solverTypesSL = validParamList->sublist(SolverTypes_name);
481 
482  const bool lapackSupportsScalar = Belos::Details::LapackSupportsScalar<Scalar>::value;
483  const bool scalarIsReal = !Teuchos::ScalarTraits<Scalar>::isComplex;
484 
485  {
487  solverTypesSL.sublist(BlockGMRES_name).setParameters(
488  *mgr.getValidParameters()
489  );
490  }
491  {
493  solverTypesSL.sublist(PseudoBlockGMRES_name).setParameters(
494  *mgr.getValidParameters()
495  );
496  }
497  if (lapackSupportsScalar) {
499  solverTypesSL.sublist(BlockCG_name).setParameters(
500  *mgr.getValidParameters()
501  );
502  }
503  if (lapackSupportsScalar) {
505  solverTypesSL.sublist(PseudoBlockCG_name).setParameters(
506  *mgr.getValidParameters()
507  );
508  }
509  {
511  solverTypesSL.sublist(PseudoBlockStochasticCG_name).setParameters(
512  *mgr.getValidParameters()
513  );
514  }
515  if (lapackSupportsScalar) {
517  solverTypesSL.sublist(GCRODR_name).setParameters(
518  *mgr.getValidParameters()
519  );
520  }
521  if (lapackSupportsScalar && scalarIsReal) {
523  solverTypesSL.sublist(RCG_name).setParameters(
524  *mgr.getValidParameters()
525  );
526  }
527  {
529  solverTypesSL.sublist(MINRES_name).setParameters(
530  *mgr.getValidParameters()
531  );
532  }
533  {
535  solverTypesSL.sublist(TFQMR_name).setParameters(
536  *mgr.getValidParameters()
537  );
538  }
539  {
541  solverTypesSL.sublist(BiCGStab_name).setParameters(
542  *mgr.getValidParameters()
543  );
544  }
545  {
547  solverTypesSL.sublist(FixedPoint_name).setParameters(
548  *mgr.getValidParameters()
549  );
550  }
551 #if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
552  {
553  Thyra::BelosTpetraGmres<Scalar,MV_t,LO_t> mgr;
554  solverTypesSL.sublist(TpetraGmres_name).setParameters(
555  *mgr.getValidParameters()
556  );
557  }
558  {
559  Thyra::BelosTpetraGmresPipeline<Scalar,MV_t,LO_t> mgr;
560  solverTypesSL.sublist(TpetraGmresPipeline_name).setParameters(
561  *mgr.getValidParameters()
562  );
563  }
564  {
565  Thyra::BelosTpetraGmresSingleReduce<Scalar,MV_t,LO_t> mgr;
566  solverTypesSL.sublist(TpetraGmresSingleReduce_name).setParameters(
567  *mgr.getValidParameters()
568  );
569  }
570  {
571  Thyra::BelosTpetraGmresSstep<Scalar,MV_t,LO_t> mgr;
572  solverTypesSL.sublist(TpetraGmresSstep_name).setParameters(
573  *mgr.getValidParameters()
574  );
575  }
576 #endif
577  }
578  return validParamList;
579 }
580 
581 
582 template<class Scalar>
583 void BelosLinearOpWithSolveFactory<Scalar>::updateThisValidParamList()
584 {
585  thisValidParamList_ = Teuchos::rcp(
586  new Teuchos::ParameterList(*generateAndGetValidParameters())
587  );
588  Teuchos::setupVerboseObjectSublist(&*thisValidParamList_);
589 }
590 
591 
592 template<class Scalar>
593 void BelosLinearOpWithSolveFactory<Scalar>::initializeOpImpl(
594  const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
595  const RCP<const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
596  const RCP<const PreconditionerBase<Scalar> > &prec_in,
597  const bool reusePrec,
598  LinearOpWithSolveBase<Scalar> *Op,
599  const ESupportSolveUse supportSolveUse
600  ) const
601 {
602 
603  using Teuchos::rcp;
604  using Teuchos::set_extra_data;
606  typedef MultiVectorBase<Scalar> MV_t;
607  typedef LinearOpBase<Scalar> LO_t;
608 
609  const RCP<Teuchos::FancyOStream> out = this->getOStream();
610  const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
611  Teuchos::OSTab tab(out);
612  if(out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
613  *out << "\nEntering Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<">::initializeOpImpl(...) ...\n";
614 
615  // These lines are changing the verbosity of the preconditioner, which has its own verbose object list,
616  // so I am commenting these out, as it is not the job of the linear solver to dictate preconditioner verbosity.
617  //typedef Teuchos::VerboseObjectTempState<PreconditionerFactoryBase<Scalar> > VOTSPF;
618  //VOTSPF precFactoryOutputTempState(precFactory_,out,verbLevel);
619 
620  TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
621  TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc.get()==NULL);
622  TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc->getOp().get()==NULL);
623  RCP<const LinearOpBase<Scalar> >
624  fwdOp = fwdOpSrc->getOp(),
625  approxFwdOp = ( approxFwdOpSrc.get() ? approxFwdOpSrc->getOp() : Teuchos::null );
626 
627  //
628  // Get the BelosLinearOpWithSolve interface
629  //
630 
631  BelosLinearOpWithSolve<Scalar>
632  *belosOp = &Teuchos::dyn_cast<BelosLinearOpWithSolve<Scalar> >(*Op);
633 
634  //
635  // Get/Create the preconditioner
636  //
637 
638  RCP<PreconditionerBase<Scalar> > myPrec = Teuchos::null;
639  RCP<const PreconditionerBase<Scalar> > prec = Teuchos::null;
640  if(prec_in.get()) {
641  // Use an externally defined preconditioner
642  prec = prec_in;
643  }
644  else {
645  // Try and generate a preconditioner on our own
646  if(precFactory_.get()) {
647  myPrec =
648  ( !belosOp->isExternalPrec()
649  ? Teuchos::rcp_const_cast<PreconditionerBase<Scalar> >(belosOp->extract_prec())
650  : Teuchos::null
651  );
652  bool hasExistingPrec = false;
653  if(myPrec.get()) {
654  hasExistingPrec = true;
655  // ToDo: Get the forward operator and validate that it is the same
656  // operator that is used here!
657  }
658  else {
659  hasExistingPrec = false;
660  myPrec = precFactory_->createPrec();
661  }
662  if( hasExistingPrec && reusePrec ) {
663  // Just reuse the existing preconditioner again!
664  }
665  else {
666  // Update the preconditioner
667  if(approxFwdOp.get())
668  precFactory_->initializePrec(approxFwdOpSrc,&*myPrec);
669  else
670  precFactory_->initializePrec(fwdOpSrc,&*myPrec);
671  }
672  prec = myPrec;
673  }
674  }
675 
676  //
677  // Uninitialize the current solver object
678  //
679 
680  bool oldIsExternalPrec = false;
681  RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > oldLP = Teuchos::null;
682  RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > oldIterSolver = Teuchos::null;
683  RCP<const LinearOpSourceBase<Scalar> > oldFwdOpSrc = Teuchos::null;
684  RCP<const LinearOpSourceBase<Scalar> > oldApproxFwdOpSrc = Teuchos::null;
685  ESupportSolveUse oldSupportSolveUse = SUPPORT_SOLVE_UNSPECIFIED;
686 
687  belosOp->uninitialize( &oldLP, NULL, &oldIterSolver, &oldFwdOpSrc,
688  NULL, &oldIsExternalPrec, &oldApproxFwdOpSrc, &oldSupportSolveUse );
689 
690  //
691  // Create the Belos linear problem
692  // NOTE: If one exists already, reuse it.
693  //
694 
696  RCP<LP_t> lp;
697  if (oldLP != Teuchos::null) {
698  lp = oldLP;
699  }
700  else {
701  lp = rcp(new LP_t());
702  }
703 
704  //
705  // Set the operator
706  //
707 
708  lp->setOperator(fwdOp);
709 
710  //
711  // Set the preconditioner
712  //
713 
714  if(prec.get()) {
715  RCP<const LinearOpBase<Scalar> > unspecified = prec->getUnspecifiedPrecOp();
716  RCP<const LinearOpBase<Scalar> > left = prec->getLeftPrecOp();
717  RCP<const LinearOpBase<Scalar> > right = prec->getRightPrecOp();
719  !( left.get() || right.get() || unspecified.get() ), std::logic_error
720  ,"Error, at least one preconditoner linear operator objects must be set!"
721  );
722  if (nonnull(unspecified)) {
723  if (paramList_->get<bool>(LeftPreconditionerIfUnspecified_name, false))
724  lp->setLeftPrec(unspecified);
725  else
726  lp->setRightPrec(unspecified);
727  }
728  else if (nonnull(left)) {
729  lp->setLeftPrec(left);
730  }
731  else if (nonnull(right)) {
732  lp->setRightPrec(right);
733  }
734  else {
735  // Set a left, right or split preconditioner
737  nonnull(left) && nonnull(right),std::logic_error
738  ,"Error, we can not currently handle split preconditioners!"
739  );
740  }
741  }
742  if(myPrec.get()) {
743  set_extra_data<RCP<PreconditionerBase<Scalar> > >(myPrec,"Belos::InternalPrec",
744  Teuchos::inOutArg(lp), Teuchos::POST_DESTROY, false);
745  }
746  else if(prec.get()) {
747  set_extra_data<RCP<const PreconditionerBase<Scalar> > >(prec,"Belos::ExternalPrec",
748  Teuchos::inOutArg(lp), Teuchos::POST_DESTROY, false);
749  }
750 
751  //
752  // Generate the parameter list.
753  //
754 
755  typedef Belos::SolverManager<Scalar,MV_t,LO_t> IterativeSolver_t;
756  RCP<IterativeSolver_t> iterativeSolver = Teuchos::null;
757  RCP<Teuchos::ParameterList> solverPL = Teuchos::rcp( new Teuchos::ParameterList() );
758 
759  switch(solverType_) {
760  case SOLVER_TYPE_BLOCK_GMRES:
761  {
762  // Set the PL
763  if(paramList_.get()) {
764  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
765  Teuchos::ParameterList &gmresPL = solverTypesPL.sublist(BlockGMRES_name);
766  solverPL = Teuchos::rcp( &gmresPL, false );
767  }
768  // Create the solver
769  if (oldIterSolver != Teuchos::null) {
770  iterativeSolver = oldIterSolver;
771  iterativeSolver->setProblem( lp );
772  iterativeSolver->setParameters( solverPL );
773  }
774  else {
775  iterativeSolver = rcp(new Belos::BlockGmresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
776  }
777  break;
778  }
779  case SOLVER_TYPE_PSEUDO_BLOCK_GMRES:
780  {
781  // Set the PL
782  if(paramList_.get()) {
783  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
784  Teuchos::ParameterList &pbgmresPL = solverTypesPL.sublist(PseudoBlockGMRES_name);
785  solverPL = Teuchos::rcp( &pbgmresPL, false );
786  }
787  //
788  // Create the solver
789  //
790  if (oldIterSolver != Teuchos::null) {
791  iterativeSolver = oldIterSolver;
792  iterativeSolver->setProblem( lp );
793  iterativeSolver->setParameters( solverPL );
794  }
795  else {
796  iterativeSolver = rcp(new Belos::PseudoBlockGmresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
797  }
798  break;
799  }
800  case SOLVER_TYPE_BLOCK_CG:
801  {
802  // Set the PL
803  if(paramList_.get()) {
804  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
805  Teuchos::ParameterList &cgPL = solverTypesPL.sublist(BlockCG_name);
806  solverPL = Teuchos::rcp( &cgPL, false );
807  }
808  // Create the solver
809  if (oldIterSolver != Teuchos::null) {
810  iterativeSolver = oldIterSolver;
811  iterativeSolver->setProblem( lp );
812  iterativeSolver->setParameters( solverPL );
813  }
814  else {
815  iterativeSolver = rcp(new Belos::BlockCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
816  }
817  break;
818  }
819  case SOLVER_TYPE_PSEUDO_BLOCK_CG:
820  {
821  // Set the PL
822  if(paramList_.get()) {
823  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
824  Teuchos::ParameterList &pbcgPL = solverTypesPL.sublist(PseudoBlockCG_name);
825  solverPL = Teuchos::rcp( &pbcgPL, false );
826  }
827  //
828  // Create the solver
829  //
830  if (oldIterSolver != Teuchos::null) {
831  iterativeSolver = oldIterSolver;
832  iterativeSolver->setProblem( lp );
833  iterativeSolver->setParameters( solverPL );
834  }
835  else {
836  iterativeSolver = rcp(new Belos::PseudoBlockCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
837  }
838  break;
839  }
840  case SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG:
841  {
842  // Set the PL
843  if(paramList_.get()) {
844  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
845  Teuchos::ParameterList &pbcgPL = solverTypesPL.sublist(PseudoBlockStochasticCG_name);
846  solverPL = Teuchos::rcp( &pbcgPL, false );
847  }
848  //
849  // Create the solver
850  //
851  if (oldIterSolver != Teuchos::null) {
852  iterativeSolver = oldIterSolver;
853  iterativeSolver->setProblem( lp );
854  iterativeSolver->setParameters( solverPL );
855  }
856  else {
857  iterativeSolver = rcp(new Belos::PseudoBlockStochasticCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
858  }
859  break;
860  }
861  case SOLVER_TYPE_GCRODR:
862  {
863  // Set the PL
864  if(paramList_.get()) {
865  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
866  Teuchos::ParameterList &gcrodrPL = solverTypesPL.sublist(GCRODR_name);
867  solverPL = Teuchos::rcp( &gcrodrPL, false );
868  }
869  // Create the solver
870  if (oldIterSolver != Teuchos::null) {
871  iterativeSolver = oldIterSolver;
872  iterativeSolver->setProblem( lp );
873  iterativeSolver->setParameters( solverPL );
874  }
875  else {
876  iterativeSolver = rcp(new Belos::GCRODRSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
877  }
878  break;
879  }
880  case SOLVER_TYPE_RCG:
881  {
882  // Set the PL
883  if(paramList_.get()) {
884  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
885  Teuchos::ParameterList &rcgPL = solverTypesPL.sublist(RCG_name);
886  solverPL = Teuchos::rcp( &rcgPL, false );
887  }
888  // Create the solver
889  if (oldIterSolver != Teuchos::null) {
890  iterativeSolver = oldIterSolver;
891  iterativeSolver->setProblem( lp );
892  iterativeSolver->setParameters( solverPL );
893  }
894  else {
895  iterativeSolver = rcp(new Belos::RCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
896  }
897  break;
898  }
899  case SOLVER_TYPE_MINRES:
900  {
901  // Set the PL
902  if(paramList_.get()) {
903  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
904  Teuchos::ParameterList &minresPL = solverTypesPL.sublist(MINRES_name);
905  solverPL = Teuchos::rcp( &minresPL, false );
906  }
907  // Create the solver
908  if (oldIterSolver != Teuchos::null) {
909  iterativeSolver = oldIterSolver;
910  iterativeSolver->setProblem( lp );
911  iterativeSolver->setParameters( solverPL );
912  }
913  else {
914  iterativeSolver = rcp(new Belos::MinresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
915  }
916  break;
917  }
918  case SOLVER_TYPE_TFQMR:
919  {
920  // Set the PL
921  if(paramList_.get()) {
922  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
923  Teuchos::ParameterList &minresPL = solverTypesPL.sublist(TFQMR_name);
924  solverPL = Teuchos::rcp( &minresPL, false );
925  }
926  // Create the solver
927  if (oldIterSolver != Teuchos::null) {
928  iterativeSolver = oldIterSolver;
929  iterativeSolver->setProblem( lp );
930  iterativeSolver->setParameters( solverPL );
931  }
932  else {
933  iterativeSolver = rcp(new Belos::TFQMRSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
934  }
935  break;
936  }
937  case SOLVER_TYPE_BICGSTAB:
938  {
939  // Set the PL
940  if(paramList_.get()) {
941  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
942  Teuchos::ParameterList &bicgstabPL = solverTypesPL.sublist(BiCGStab_name);
943  solverPL = Teuchos::rcp( &bicgstabPL, false );
944  }
945  // Create the solver
946  if (oldIterSolver != Teuchos::null) {
947  iterativeSolver = oldIterSolver;
948  iterativeSolver->setProblem( lp );
949  iterativeSolver->setParameters( solverPL );
950  }
951  else {
952  iterativeSolver = rcp(new Belos::BiCGStabSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
953  }
954  break;
955  }
956  case SOLVER_TYPE_FIXEDPOINT:
957  {
958  // Set the PL
959  if(paramList_.get()) {
960  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
961  Teuchos::ParameterList &fixedPointPL = solverTypesPL.sublist(FixedPoint_name);
962  solverPL = Teuchos::rcp( &fixedPointPL, false );
963  }
964  // Create the solver
965  if (oldIterSolver != Teuchos::null) {
966  iterativeSolver = oldIterSolver;
967  iterativeSolver->setProblem( lp );
968  iterativeSolver->setParameters( solverPL );
969  }
970  else {
971  iterativeSolver = rcp(new Belos::FixedPointSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
972  }
973  break;
974  }
975 #if defined(HAVE_BELOS_TPETRA) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS)
976  case SOLVER_TYPE_TPETRA_GMRES:
977  {
978  // Get the PL
979  if(paramList_.get()) {
980  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
981  Teuchos::ParameterList &tpetraGmresPL = solverTypesPL.sublist(TpetraGmres_name);
982  solverPL = Teuchos::rcp( &tpetraGmresPL, false );
983  }
984  // Create the solver, always
985  iterativeSolver = rcp(new Thyra::BelosTpetraGmres<Scalar,MV_t,LO_t>());
986  // Set the Problem & PL
987  iterativeSolver->setProblem( lp );
988  iterativeSolver->setParameters( solverPL );
989  break;
990  }
991  case SOLVER_TYPE_TPETRA_GMRES_PIPELINE:
992  {
993  // Get the PL
994  if(paramList_.get()) {
995  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
996  Teuchos::ParameterList &tpetraGmresPipelinePL = solverTypesPL.sublist(TpetraGmresPipeline_name);
997  solverPL = Teuchos::rcp( &tpetraGmresPipelinePL, false );
998  }
999  // Create the solver, always
1000  iterativeSolver = rcp(new Thyra::BelosTpetraGmresPipeline<Scalar,MV_t,LO_t>());
1001  // Set the Problem & PL
1002  iterativeSolver->setProblem( lp );
1003  iterativeSolver->setParameters( solverPL );
1004  break;
1005  }
1006  case SOLVER_TYPE_TPETRA_GMRES_SINGLE_REDUCE:
1007  {
1008  // Get the PL
1009  if(paramList_.get()) {
1010  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
1011  Teuchos::ParameterList &tpetraGmresSingleReducePL = solverTypesPL.sublist(TpetraGmresSingleReduce_name);
1012  solverPL = Teuchos::rcp( &tpetraGmresSingleReducePL, false );
1013  }
1014  // Create the solver, always
1015  iterativeSolver = rcp(new Thyra::BelosTpetraGmresSingleReduce<Scalar,MV_t,LO_t>());
1016  // Set the Problem & PL
1017  iterativeSolver->setProblem( lp );
1018  iterativeSolver->setParameters( solverPL );
1019  break;
1020  }
1021  case SOLVER_TYPE_TPETRA_GMRES_SSTEP:
1022  {
1023  // Get the PL
1024  if(paramList_.get()) {
1025  Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
1026  Teuchos::ParameterList &tpetraGmresSstepPL = solverTypesPL.sublist(TpetraGmresSstep_name);
1027  solverPL = Teuchos::rcp( &tpetraGmresSstepPL, false );
1028  }
1029  // Create the solver, always
1030  iterativeSolver = rcp(new Thyra::BelosTpetraGmresSstep<Scalar,MV_t,LO_t>());
1031  // Set the Problem & PL
1032  iterativeSolver->setProblem( lp );
1033  iterativeSolver->setParameters( solverPL );
1034  break;
1035  }
1036 #endif
1037  default:
1038  {
1040  }
1041  }
1042 
1043  //
1044  // Initialize the LOWS object
1045  //
1046 
1047  belosOp->initialize(
1048  lp, solverPL, iterativeSolver,
1049  fwdOpSrc, prec, myPrec.get()==NULL, approxFwdOpSrc,
1050  supportSolveUse, convergenceTestFrequency_
1051  );
1052  belosOp->setOStream(out);
1053  belosOp->setVerbLevel(verbLevel);
1054 #ifdef TEUCHOS_DEBUG
1055  if(paramList_.get()) {
1056  // Make sure we read the list correctly
1057  paramList_->validateParameters(*this->getValidParameters(),1); // Validate 0th and 1st level deep
1058  }
1059 #endif
1060  if(out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
1061  *out << "\nLeaving Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<">::initializeOpImpl(...) ...\n";
1062 
1063 }
1064 
1065 
1066 } // namespace Thyra
1067 
1068 
1069 #endif // THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
void uninitializeOp(LinearOpWithSolveBase< Scalar > *Op, Teuchos::RCP< const LinearOpSourceBase< Scalar > > *fwdOpSrc, Teuchos::RCP< const PreconditionerBase< Scalar > > *prec, Teuchos::RCP< const LinearOpSourceBase< Scalar > > *approxFwdOpSrc, ESupportSolveUse *supportSolveUse) const
RCP< const LinearOpSourceBase< Scalar > > extract_approxFwdOpSrc()
static void addConverter(RCP< const ParameterEntryValidator > validator, RCP< ValidatorXMLConverter > converterToAdd)
Concrete LinearOpWithSolveBase subclass in terms of Belos.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
RCP< const PreconditionerBase< Scalar > > extract_prec()
T_To & dyn_cast(T_From &from)
void initializeAndReuseOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, LinearOpWithSolveBase< Scalar > *Op) const
T * get() const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
void unsetPreconditionerFactory(Teuchos::RCP< PreconditionerFactoryBase< Scalar > > *precFactory, std::string *precFactoryName)
Thyra specializations of MultiVecTraits and OperatorTraits.
void setParameterList(Teuchos::RCP< Teuchos::ParameterList > const &paramList)
bool isCompatible(const LinearOpSourceBase< Scalar > &fwdOpSrc) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
virtual std::string description() const
void initializeOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
void validateParametersAndSetDefaults(ParameterList const &validParamList, int const depth=1000)
void setPreconditionerFactory(const Teuchos::RCP< PreconditionerFactoryBase< Scalar > > &precFactory, const std::string &precFactoryName)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
LinearOpWithSolveFactoryBase subclass implemented in terms of Belos.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
ParameterList & setParameters(const ParameterList &source)
ESupportSolveUse
bool supportsPreconditionerInputType(const EPreconditionerInputType precOpType) const
Teuchos::RCP< const Teuchos::ParameterList > getParameterList() const
RCP< const LinearOpSourceBase< Scalar > > extract_fwdOpSrc()
bool nonnull(const boost::shared_ptr< T > &p)
Teuchos::RCP< PreconditionerFactoryBase< Scalar > > getPreconditionerFactory() const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
void initializeApproxPreconditionedOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &approxFwdOpSrc, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
TypeTo as(const TypeFrom &t)
void initializePreconditionedOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const Teuchos::RCP< const PreconditionerBase< Scalar > > &prec, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
EPreconditionerInputType
Teuchos::RCP< LinearOpWithSolveBase< Scalar > > createOp() const
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
BelosLinearOpWithSolveFactory()
Construct without preconditioner factory.
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)

Generated on Thu Nov 21 2024 09:22:16 for Stratimikos by doxygen 1.8.5