45#ifndef THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP
46#define THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP
49#include "Thyra_BelosLinearOpWithSolveFactory_decl.hpp"
50#include "Thyra_BelosLinearOpWithSolve.hpp"
51#include "Thyra_ScaledAdjointLinearOpBase.hpp"
66#include "Thyra_BelosTpetrasSolverAdapter.hpp"
68#include "Teuchos_VerboseObjectParameterListHelpers.hpp"
69#include "Teuchos_StandardParameterEntryValidators.hpp"
70#include "Teuchos_ParameterList.hpp"
71#include "Teuchos_dyn_cast.hpp"
72#include "Teuchos_ValidatorXMLConverterDB.hpp"
73#include "Teuchos_StandardValidatorXMLConverters.hpp"
101template<
class Scalar>
103template<
class Scalar>
105template<
class Scalar>
107template<
class Scalar>
109template<
class Scalar>
111template<
class Scalar>
113template<
class Scalar>
115template<
class Scalar>
118template<
class Scalar>
122const std::string LeftPreconditionerIfUnspecified_name =
"Left Preconditioner If Unspecified";
128template<
class Scalar>
130 :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES),
131 convergenceTestFrequency_(1)
133 updateThisValidParamList();
137template<
class Scalar>
139 const RCP<PreconditionerFactoryBase<Scalar> > &precFactory
141 :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES)
150template<
class Scalar>
157template<
class Scalar>
159 const RCP<PreconditionerFactoryBase<Scalar> > &precFactory,
160 const std::string &precFactoryName
163 TEUCHOS_TEST_FOR_EXCEPT(!precFactory.get());
164 RCP<const Teuchos::ParameterList>
165 precFactoryValidPL = precFactory->getValidParameters();
166 const std::string _precFactoryName =
167 ( precFactoryName !=
""
169 : ( precFactoryValidPL.get() ? precFactoryValidPL->name() :
"GENERIC PRECONDITIONER FACTORY" )
171 precFactory_ = precFactory;
172 precFactoryName_ = _precFactoryName;
173 updateThisValidParamList();
177template<
class Scalar>
178RCP<PreconditionerFactoryBase<Scalar> >
185template<
class Scalar>
187 RCP<PreconditionerFactoryBase<Scalar> > *precFactory,
188 std::string *precFactoryName
191 if(precFactory) *precFactory = precFactory_;
192 if(precFactoryName) *precFactoryName = precFactoryName_;
193 precFactory_ = Teuchos::null;
194 precFactoryName_ =
"";
195 updateThisValidParamList();
199template<
class Scalar>
201 const LinearOpSourceBase<Scalar> &fwdOpSrc
204 if(precFactory_.get())
205 return precFactory_->isCompatible(fwdOpSrc);
210template<
class Scalar>
211RCP<LinearOpWithSolveBase<Scalar> >
218template<
class Scalar>
220 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
221 LinearOpWithSolveBase<Scalar> *Op,
222 const ESupportSolveUse supportSolveUse
226 initializeOpImpl(fwdOpSrc,null,null,
false,Op,supportSolveUse);
230template<
class Scalar>
232 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
233 LinearOpWithSolveBase<Scalar> *Op
237 initializeOpImpl(fwdOpSrc,null,null,
true,Op,SUPPORT_SOLVE_UNSPECIFIED);
241template<
class Scalar>
243 const EPreconditionerInputType precOpType
246 if(precFactory_.get())
248 return (precOpType==PRECONDITIONER_INPUT_TYPE_AS_OPERATOR);
252template<
class Scalar>
254 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
255 const RCP<
const PreconditionerBase<Scalar> > &prec,
256 LinearOpWithSolveBase<Scalar> *Op,
257 const ESupportSolveUse supportSolveUse
261 initializeOpImpl(fwdOpSrc,null,prec,
false,Op,supportSolveUse);
265template<
class Scalar>
267 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
268 const RCP<
const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
269 LinearOpWithSolveBase<Scalar> *Op,
270 const ESupportSolveUse supportSolveUse
274 initializeOpImpl(fwdOpSrc,approxFwdOpSrc,null,
false,Op,supportSolveUse);
278template<
class Scalar>
280 LinearOpWithSolveBase<Scalar> *Op,
281 RCP<
const LinearOpSourceBase<Scalar> > *fwdOpSrc,
282 RCP<
const PreconditionerBase<Scalar> > *prec,
283 RCP<
const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
284 ESupportSolveUse *supportSolveUse
288 TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
291 &belosOp = Teuchos::dyn_cast<BelosLinearOpWithSolve<Scalar> >(*Op);
292 RCP<const LinearOpSourceBase<Scalar> >
294 RCP<const PreconditionerBase<Scalar> >
299 RCP<const LinearOpSourceBase<Scalar> >
303 if(fwdOpSrc) *fwdOpSrc = _fwdOpSrc;
304 if(prec) *prec = _prec;
305 if(approxFwdOpSrc) *approxFwdOpSrc = _approxFwdOpSrc;
306 if(supportSolveUse) *supportSolveUse = _supportSolveUse;
313template<
class Scalar>
315 RCP<Teuchos::ParameterList>
const& paramList
318 TEUCHOS_TEST_FOR_EXCEPT(paramList.get()==NULL);
320 paramList_ = paramList;
322 Teuchos::getIntegralValue<EBelosSolverType>(*paramList_,
SolverType_name);
323 convergenceTestFrequency_ =
325 Teuchos::readVerboseObjectSublist(&*paramList_,
this);
329template<
class Scalar>
330RCP<Teuchos::ParameterList>
337template<
class Scalar>
338RCP<Teuchos::ParameterList>
341 RCP<Teuchos::ParameterList> _paramList = paramList_;
342 paramList_ = Teuchos::null;
347template<
class Scalar>
348RCP<const Teuchos::ParameterList>
355template<
class Scalar>
356RCP<const Teuchos::ParameterList>
359 return thisValidParamList_;
366template<
class Scalar>
369 std::ostringstream oss;
370 oss << Teuchos::Describable::description();
372 if (nonnull(precFactory_)) {
373 oss << precFactory_->description();
383template<
class Scalar>
384RCP<const Teuchos::ParameterList>
385BelosLinearOpWithSolveFactory<Scalar>::generateAndGetValidParameters()
388 using Teuchos::tuple;
389 using Teuchos::setStringToIntegralParameter;
390Teuchos::ValidatorXMLConverterDB::addConverter(
391 Teuchos::DummyObjectGetter<
392 Teuchos::StringToIntegralParameterEntryValidator<EBelosSolverType>
394 Teuchos::DummyObjectGetter<Teuchos::StringToIntegralValidatorXMLConverter<
395 EBelosSolverType> >::getDummyObject());
397 typedef MultiVectorBase<Scalar> MV_t;
398 typedef LinearOpBase<Scalar> LO_t;
399 static RCP<Teuchos::ParameterList> validParamList;
400 if(validParamList.get()==NULL) {
401 validParamList = Teuchos::rcp(
new Teuchos::ParameterList(
"BelosLinearOpWithSolveFactory"));
402 setStringToIntegralParameter<EBelosSolverType>(
403 SolverType_name, SolverType_default,
404 "Type of linear solver algorithm to use.",
407 "Pseudo Block GMRES",
410 "Pseudo Block Stochastic CG",
418 "TPETRA GMRES PIPELINE",
419 "TPETRA GMRES SINGLE REDUCE",
420 "TPETRA GMRES S-STEP"
423 "Block GMRES solver for nonsymmetric linear systems. It can also solve "
424 "single right-hand side systems, and can also perform Flexible GMRES "
425 "(where the preconditioner may change at every iteration, for example "
426 "for inner-outer iterations), by setting options in the \"Block GMRES\" "
429 "GMRES solver for nonsymmetric linear systems, that performs single "
430 "right-hand side solves on multiple right-hand sides at once. It "
431 "exploits operator multivector multiplication in order to amortize "
432 "global communication costs. Individual linear systems are deflated "
433 "out as they are solved.",
435 "Block CG solver for symmetric (Hermitian in complex arithmetic) "
436 "positive definite linear systems. It can also solve single "
437 "right-hand-side systems.",
439 "CG solver that performs single right-hand side CG on multiple right-hand "
440 "sides at once. It exploits operator multivector multiplication in order "
441 "to amortize global communication costs. Individual linear systems are "
442 "deflated out as they are solved.",
444 "stochastic CG solver that performs single right-hand side CG on multiple right-hand "
445 "sides at once. It exploits operator multivector multiplication in order "
446 "to amortize global communication costs. Individual linear systems are "
447 "deflated out as they are solved. [EXPERIMENTAL]",
449 "Variant of GMRES that performs subspace recycling to accelerate "
450 "convergence for sequences of solves with related linear systems. "
451 "Individual linear systems are deflated out as they are solved. "
452 "The current implementation only supports real-valued Scalar types.",
454 "CG solver for symmetric (Hermitian in complex arithmetic) positive "
455 "definite linear systems, that performs subspace recycling to "
456 "accelerate convergence for sequences of related linear systems.",
458 "MINRES solver for symmetric indefinite linear systems. It performs "
459 "single-right-hand-side solves on multiple right-hand sides sequentially.",
461 "TFQMR (Transpose-Free QMR) solver for nonsymmetric linear systems.",
463 "BiCGStab solver for nonsymmetric linear systems.",
465 "Fixed point iteration",
467 "Native Tpetra implementation of GMRES",
469 "Native Tpetra implementation of pipeline GMRES",
471 "Native Tpetra implementation of single-reduce GMRES",
473 "Native Tpetra implementation of s-step GMRES"
475 tuple<EBelosSolverType>(
476 SOLVER_TYPE_BLOCK_GMRES,
477 SOLVER_TYPE_PSEUDO_BLOCK_GMRES,
478 SOLVER_TYPE_BLOCK_CG,
479 SOLVER_TYPE_PSEUDO_BLOCK_CG,
480 SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG,
485 SOLVER_TYPE_BICGSTAB,
486 SOLVER_TYPE_FIXEDPOINT,
487 SOLVER_TYPE_TPETRA_GMRES,
488 SOLVER_TYPE_TPETRA_GMRES_PIPELINE,
489 SOLVER_TYPE_TPETRA_GMRES_SINGLE_REDUCE,
490 SOLVER_TYPE_TPETRA_GMRES_SSTEP
494 validParamList->set(ConvergenceTestFrequency_name, as<int>(1),
495 "Number of linear solver iterations to skip between applying"
496 " user-defined convergence test.");
498 LeftPreconditionerIfUnspecified_name,
false,
499 "If the preconditioner does not specify if it is left or right, and this\n"
500 "option is set to true, put the preconditioner on the left side.\n"
501 "Historically, preconditioning is on the right. Some solvers may not\n"
502 "support left preconditioning.");
503 Teuchos::ParameterList
504 &solverTypesSL = validParamList->sublist(SolverTypes_name);
507 const bool scalarIsReal = !Teuchos::ScalarTraits<Scalar>::isComplex;
521 if (lapackSupportsScalar) {
522 Belos::BlockCGSolMgr<Scalar,MV_t,LO_t> mgr;
524 *mgr.getValidParameters()
527 if (lapackSupportsScalar) {
528 Belos::PseudoBlockCGSolMgr<Scalar,MV_t,LO_t> mgr;
529 solverTypesSL.sublist(PseudoBlockCG_name).setParameters(
530 *mgr.getValidParameters()
534 Belos::PseudoBlockStochasticCGSolMgr<Scalar,MV_t,LO_t> mgr;
535 solverTypesSL.sublist(PseudoBlockStochasticCG_name).setParameters(
539 if (lapackSupportsScalar) {
540 Belos::GCRODRSolMgr<Scalar,MV_t,LO_t> mgr;
542 *mgr.getValidParameters()
545 if (lapackSupportsScalar && scalarIsReal) {
546 Belos::RCGSolMgr<Scalar,MV_t,LO_t> mgr;
547 solverTypesSL.sublist(RCG_name).setParameters(
548 *mgr.getValidParameters()
552 Belos::MinresSolMgr<Scalar,MV_t,LO_t> mgr;
553 solverTypesSL.sublist(MINRES_name).setParameters(
558 Belos::TFQMRSolMgr<Scalar,MV_t,LO_t> mgr;
564 Belos::BiCGStabSolMgr<Scalar,MV_t,LO_t> mgr;
570 Belos::FixedPointSolMgr<Scalar,MV_t,LO_t> mgr;
576 Thyra::BelosTpetraGmres<Scalar,MV_t,LO_t> mgr;
578 *mgr.getValidParameters()
582 Thyra::BelosTpetraGmresPipeline<Scalar,MV_t,LO_t> mgr;
583 solverTypesSL.sublist(TpetraGmresPipeline_name).setParameters(
584 *mgr.getValidParameters()
588 Thyra::BelosTpetraGmresSingleReduce<Scalar,MV_t,LO_t> mgr;
589 solverTypesSL.sublist(TpetraGmresSingleReduce_name).setParameters(
590 *mgr.getValidParameters()
594 Thyra::BelosTpetraGmresSstep<Scalar,MV_t,LO_t> mgr;
595 solverTypesSL.sublist(TpetraGmresSstep_name).setParameters(
596 *mgr.getValidParameters()
600 return validParamList;
604template<
class Scalar>
605void BelosLinearOpWithSolveFactory<Scalar>::updateThisValidParamList()
607 thisValidParamList_ = Teuchos::rcp(
608 new Teuchos::ParameterList(*generateAndGetValidParameters())
610 Teuchos::setupVerboseObjectSublist(&*thisValidParamList_);
614template<
class Scalar>
615void BelosLinearOpWithSolveFactory<Scalar>::initializeOpImpl(
616 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
617 const RCP<
const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
618 const RCP<
const PreconditionerBase<Scalar> > &prec_in,
619 const bool reusePrec,
620 LinearOpWithSolveBase<Scalar> *Op,
621 const ESupportSolveUse supportSolveUse
626 using Teuchos::set_extra_data;
627 typedef Teuchos::ScalarTraits<Scalar> ST;
628 typedef MultiVectorBase<Scalar> MV_t;
629 typedef LinearOpBase<Scalar> LO_t;
631 const RCP<Teuchos::FancyOStream> out = this->getOStream();
632 const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
633 Teuchos::OSTab tab(out);
634 if(out.get() &&
static_cast<int>(verbLevel) >
static_cast<int>(Teuchos::VERB_LOW))
635 *out <<
"\nEntering Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<
">::initializeOpImpl(...) ...\n";
642 TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
643 TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc.get()==NULL);
644 TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc->getOp().get()==NULL);
645 RCP<const LinearOpBase<Scalar> >
646 fwdOp = fwdOpSrc->getOp(),
647 approxFwdOp = ( approxFwdOpSrc.get() ? approxFwdOpSrc->getOp() : Teuchos::null );
654 *belosOp = &Teuchos::dyn_cast<BelosLinearOpWithSolve<Scalar> >(*Op);
660 RCP<PreconditionerBase<Scalar> > myPrec = Teuchos::null;
661 RCP<const PreconditionerBase<Scalar> > prec = Teuchos::null;
668 if(precFactory_.get()) {
670 ( !belosOp->isExternalPrec()
671 ? Teuchos::rcp_const_cast<PreconditionerBase<Scalar> >(belosOp->extract_prec())
674 bool hasExistingPrec =
false;
676 hasExistingPrec =
true;
681 hasExistingPrec =
false;
682 myPrec = precFactory_->createPrec();
684 if( hasExistingPrec && reusePrec ) {
689 if(approxFwdOp.get())
690 precFactory_->initializePrec(approxFwdOpSrc,&*myPrec);
692 precFactory_->initializePrec(fwdOpSrc,&*myPrec);
702 bool oldIsExternalPrec =
false;
703 RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > oldLP = Teuchos::null;
704 RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > oldIterSolver = Teuchos::null;
705 RCP<const LinearOpSourceBase<Scalar> > oldFwdOpSrc = Teuchos::null;
706 RCP<const LinearOpSourceBase<Scalar> > oldApproxFwdOpSrc = Teuchos::null;
707 ESupportSolveUse oldSupportSolveUse = SUPPORT_SOLVE_UNSPECIFIED;
709 belosOp->
uninitialize( &oldLP, NULL, &oldIterSolver, &oldFwdOpSrc,
710 NULL, &oldIsExternalPrec, &oldApproxFwdOpSrc, &oldSupportSolveUse );
717 typedef Belos::LinearProblem<Scalar,MV_t,LO_t> LP_t;
719 if (oldLP != Teuchos::null) {
723 lp = rcp(
new LP_t());
730 lp->setOperator(fwdOp);
737 RCP<const LinearOpBase<Scalar> > unspecified = prec->getUnspecifiedPrecOp();
738 RCP<const LinearOpBase<Scalar> > left = prec->getLeftPrecOp();
739 RCP<const LinearOpBase<Scalar> > right = prec->getRightPrecOp();
740 TEUCHOS_TEST_FOR_EXCEPTION(
741 !( left.get() || right.get() || unspecified.get() ), std::logic_error
742 ,
"Error, at least one preconditoner linear operator objects must be set!"
744 if (nonnull(unspecified)) {
745 if (paramList_->get<
bool>(LeftPreconditionerIfUnspecified_name,
false))
746 lp->setLeftPrec(unspecified);
748 lp->setRightPrec(unspecified);
750 else if (nonnull(left)) {
751 lp->setLeftPrec(left);
753 else if (nonnull(right)) {
754 lp->setRightPrec(right);
758 TEUCHOS_TEST_FOR_EXCEPTION(
759 nonnull(left) && nonnull(right),std::logic_error
760 ,
"Error, we can not currently handle split preconditioners!"
765 set_extra_data<RCP<PreconditionerBase<Scalar> > >(myPrec,
"Belos::InternalPrec",
766 Teuchos::inOutArg(lp), Teuchos::POST_DESTROY,
false);
768 else if(prec.get()) {
769 set_extra_data<RCP<const PreconditionerBase<Scalar> > >(prec,
"Belos::ExternalPrec",
770 Teuchos::inOutArg(lp), Teuchos::POST_DESTROY,
false);
777 typedef Belos::SolverManager<Scalar,MV_t,LO_t> IterativeSolver_t;
778 RCP<IterativeSolver_t> iterativeSolver = Teuchos::null;
779 RCP<Teuchos::ParameterList> solverPL = Teuchos::rcp(
new Teuchos::ParameterList() );
781 switch(solverType_) {
782 case SOLVER_TYPE_BLOCK_GMRES:
785 if(paramList_.get()) {
786 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
787 Teuchos::ParameterList &gmresPL = solverTypesPL.sublist(BlockGMRES_name);
788 solverPL = Teuchos::rcp( &gmresPL,
false );
791 if (oldIterSolver != Teuchos::null) {
792 iterativeSolver = oldIterSolver;
793 iterativeSolver->setProblem( lp );
794 iterativeSolver->setParameters( solverPL );
797 iterativeSolver = rcp(
new Belos::BlockGmresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
801 case SOLVER_TYPE_PSEUDO_BLOCK_GMRES:
804 if(paramList_.get()) {
805 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
806 Teuchos::ParameterList &pbgmresPL = solverTypesPL.sublist(PseudoBlockGMRES_name);
807 solverPL = Teuchos::rcp( &pbgmresPL,
false );
812 if (oldIterSolver != Teuchos::null) {
813 iterativeSolver = oldIterSolver;
814 iterativeSolver->setProblem( lp );
815 iterativeSolver->setParameters( solverPL );
818 iterativeSolver = rcp(
new Belos::PseudoBlockGmresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
822 case SOLVER_TYPE_BLOCK_CG:
825 if(paramList_.get()) {
826 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
827 Teuchos::ParameterList &cgPL = solverTypesPL.sublist(BlockCG_name);
828 solverPL = Teuchos::rcp( &cgPL,
false );
831 if (oldIterSolver != Teuchos::null) {
832 iterativeSolver = oldIterSolver;
833 iterativeSolver->setProblem( lp );
834 iterativeSolver->setParameters( solverPL );
837 iterativeSolver = rcp(
new Belos::BlockCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
841 case SOLVER_TYPE_PSEUDO_BLOCK_CG:
844 if(paramList_.get()) {
845 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
846 Teuchos::ParameterList &pbcgPL = solverTypesPL.sublist(PseudoBlockCG_name);
847 solverPL = Teuchos::rcp( &pbcgPL,
false );
852 if (oldIterSolver != Teuchos::null) {
853 iterativeSolver = oldIterSolver;
854 iterativeSolver->setProblem( lp );
855 iterativeSolver->setParameters( solverPL );
858 iterativeSolver = rcp(
new Belos::PseudoBlockCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
862 case SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG:
865 if(paramList_.get()) {
866 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
867 Teuchos::ParameterList &pbcgPL = solverTypesPL.sublist(PseudoBlockStochasticCG_name);
868 solverPL = Teuchos::rcp( &pbcgPL,
false );
873 if (oldIterSolver != Teuchos::null) {
874 iterativeSolver = oldIterSolver;
875 iterativeSolver->setProblem( lp );
876 iterativeSolver->setParameters( solverPL );
879 iterativeSolver = rcp(
new Belos::PseudoBlockStochasticCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
883 case SOLVER_TYPE_GCRODR:
886 if(paramList_.get()) {
887 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
888 Teuchos::ParameterList &gcrodrPL = solverTypesPL.sublist(GCRODR_name);
889 solverPL = Teuchos::rcp( &gcrodrPL,
false );
892 if (oldIterSolver != Teuchos::null) {
893 iterativeSolver = oldIterSolver;
894 iterativeSolver->setProblem( lp );
895 iterativeSolver->setParameters( solverPL );
898 iterativeSolver = rcp(
new Belos::GCRODRSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
902 case SOLVER_TYPE_RCG:
905 if(paramList_.get()) {
906 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
907 Teuchos::ParameterList &rcgPL = solverTypesPL.sublist(RCG_name);
908 solverPL = Teuchos::rcp( &rcgPL,
false );
911 if (oldIterSolver != Teuchos::null) {
912 iterativeSolver = oldIterSolver;
913 iterativeSolver->setProblem( lp );
914 iterativeSolver->setParameters( solverPL );
917 iterativeSolver = rcp(
new Belos::RCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
921 case SOLVER_TYPE_MINRES:
924 if(paramList_.get()) {
925 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
926 Teuchos::ParameterList &minresPL = solverTypesPL.sublist(MINRES_name);
927 solverPL = Teuchos::rcp( &minresPL,
false );
930 if (oldIterSolver != Teuchos::null) {
931 iterativeSolver = oldIterSolver;
932 iterativeSolver->setProblem( lp );
933 iterativeSolver->setParameters( solverPL );
936 iterativeSolver = rcp(
new Belos::MinresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
940 case SOLVER_TYPE_TFQMR:
943 if(paramList_.get()) {
944 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
945 Teuchos::ParameterList &minresPL = solverTypesPL.sublist(TFQMR_name);
946 solverPL = Teuchos::rcp( &minresPL,
false );
949 if (oldIterSolver != Teuchos::null) {
950 iterativeSolver = oldIterSolver;
951 iterativeSolver->setProblem( lp );
952 iterativeSolver->setParameters( solverPL );
955 iterativeSolver = rcp(
new Belos::TFQMRSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
959 case SOLVER_TYPE_BICGSTAB:
962 if(paramList_.get()) {
963 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
964 Teuchos::ParameterList &bicgstabPL = solverTypesPL.sublist(BiCGStab_name);
965 solverPL = Teuchos::rcp( &bicgstabPL,
false );
968 if (oldIterSolver != Teuchos::null) {
969 iterativeSolver = oldIterSolver;
970 iterativeSolver->setProblem( lp );
971 iterativeSolver->setParameters( solverPL );
974 iterativeSolver = rcp(
new Belos::BiCGStabSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
978 case SOLVER_TYPE_FIXEDPOINT:
981 if(paramList_.get()) {
982 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
983 Teuchos::ParameterList &fixedPointPL = solverTypesPL.sublist(FixedPoint_name);
984 solverPL = Teuchos::rcp( &fixedPointPL,
false );
987 if (oldIterSolver != Teuchos::null) {
988 iterativeSolver = oldIterSolver;
989 iterativeSolver->setProblem( lp );
990 iterativeSolver->setParameters( solverPL );
993 iterativeSolver = rcp(
new Belos::FixedPointSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
997 case SOLVER_TYPE_TPETRA_GMRES:
1000 if(paramList_.get()) {
1001 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
1002 Teuchos::ParameterList &tpetraGmresPL = solverTypesPL.sublist(TpetraGmres_name);
1003 solverPL = Teuchos::rcp( &tpetraGmresPL,
false );
1006 iterativeSolver = rcp(
new Thyra::BelosTpetraGmres<Scalar,MV_t,LO_t>());
1008 iterativeSolver->setProblem( lp );
1009 iterativeSolver->setParameters( solverPL );
1012 case SOLVER_TYPE_TPETRA_GMRES_PIPELINE:
1015 if(paramList_.get()) {
1016 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
1017 Teuchos::ParameterList &tpetraGmresPipelinePL = solverTypesPL.sublist(TpetraGmresPipeline_name);
1018 solverPL = Teuchos::rcp( &tpetraGmresPipelinePL,
false );
1021 iterativeSolver = rcp(
new Thyra::BelosTpetraGmresPipeline<Scalar,MV_t,LO_t>());
1023 iterativeSolver->setProblem( lp );
1024 iterativeSolver->setParameters( solverPL );
1027 case SOLVER_TYPE_TPETRA_GMRES_SINGLE_REDUCE:
1030 if(paramList_.get()) {
1031 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
1032 Teuchos::ParameterList &tpetraGmresSingleReducePL = solverTypesPL.sublist(TpetraGmresSingleReduce_name);
1033 solverPL = Teuchos::rcp( &tpetraGmresSingleReducePL,
false );
1036 iterativeSolver = rcp(
new Thyra::BelosTpetraGmresSingleReduce<Scalar,MV_t,LO_t>());
1038 iterativeSolver->setProblem( lp );
1039 iterativeSolver->setParameters( solverPL );
1042 case SOLVER_TYPE_TPETRA_GMRES_SSTEP:
1045 if(paramList_.get()) {
1046 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
1047 Teuchos::ParameterList &tpetraGmresSstepPL = solverTypesPL.sublist(TpetraGmresSstep_name);
1048 solverPL = Teuchos::rcp( &tpetraGmresSstepPL,
false );
1051 iterativeSolver = rcp(
new Thyra::BelosTpetraGmresSstep<Scalar,MV_t,LO_t>());
1053 iterativeSolver->setProblem( lp );
1054 iterativeSolver->setParameters( solverPL );
1060 TEUCHOS_TEST_FOR_EXCEPT(
true);
1068 belosOp->initialize(
1069 lp, solverPL, iterativeSolver,
1070 fwdOpSrc, prec, myPrec.get()==NULL, approxFwdOpSrc,
1071 supportSolveUse, convergenceTestFrequency_
1073 belosOp->setOStream(out);
1074 belosOp->setVerbLevel(verbLevel);
1076 if(paramList_.get()) {
1078 paramList_->validateParameters(*this->getValidParameters(),1);
1081 if(out.get() &&
static_cast<int>(verbLevel) >
static_cast<int>(Teuchos::VERB_LOW))
1082 *out <<
"\nLeaving Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<
">::initializeOpImpl(...) ...\n";
Thyra specializations of MultiVecTraits and OperatorTraits.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms) override
static const std::string TpetraGmresSstep_name
void unsetPreconditionerFactory(Teuchos::RCP< PreconditionerFactoryBase< Scalar > > *precFactory, std::string *precFactoryName)
bool supportsPreconditionerInputType(const EPreconditionerInputType precOpType) const
bool isCompatible(const LinearOpSourceBase< Scalar > &fwdOpSrc) const
void setParameterList(Teuchos::RCP< Teuchos::ParameterList > const ¶mList)
static const std::string TFQMR_name
static const std::string ConvergenceTestFrequency_name
static const std::string TpetraGmresPipeline_name
static const std::string TpetraGmresSingleReduce_name
static const std::string BlockGMRES_name
static const std::string SolverType_name
void initializePreconditionedOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const Teuchos::RCP< const PreconditionerBase< Scalar > > &prec, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
std::string description() const
static const std::string SolverTypes_name
static const std::string MINRES_name
static const std::string SolverType_default
static const std::string BiCGStab_name
void setPreconditionerFactory(const Teuchos::RCP< PreconditionerFactoryBase< Scalar > > &precFactory, const std::string &precFactoryName)
static const std::string GCRODR_name
static const std::string TpetraGmres_name
void initializeAndReuseOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, LinearOpWithSolveBase< Scalar > *Op) const
static const std::string PseudoBlockCG_name
Teuchos::RCP< const Teuchos::ParameterList > getParameterList() const
static const std::string RCG_name
static const std::string BlockCG_name
void initializeOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
bool acceptsPreconditionerFactory() const
Returns true .
static const std::string FixedPoint_name
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
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
BelosLinearOpWithSolveFactory()
Construct without preconditioner factory.
void initializeApproxPreconditionedOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &approxFwdOpSrc, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
Teuchos::RCP< LinearOpWithSolveBase< Scalar > > createOp() const
static const std::string PseudoBlockGMRES_name
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
static const std::string PseudoBlockStochasticCG_name
Teuchos::RCP< PreconditionerFactoryBase< Scalar > > getPreconditionerFactory() const
Concrete LinearOpWithSolveBase subclass in terms of Belos.
void uninitialize(RCP< Belos::LinearProblem< Scalar, MV_t, LO_t > > *lp=NULL, RCP< Teuchos::ParameterList > *solverPL=NULL, RCP< Belos::SolverManager< Scalar, MV_t, LO_t > > *iterativeSolver=NULL, RCP< const LinearOpSourceBase< Scalar > > *fwdOpSrc=NULL, RCP< const PreconditionerBase< Scalar > > *prec=NULL, bool *isExternalPrec=NULL, RCP< const LinearOpSourceBase< Scalar > > *approxFwdOpSrc=NULL, ESupportSolveUse *supportSolveUse=NULL)
Uninitializes and returns stored quantities.
ESupportSolveUse supportSolveUse() const
RCP< const LinearOpSourceBase< Scalar > > extract_approxFwdOpSrc()
RCP< const PreconditionerBase< Scalar > > extract_prec()
bool isExternalPrec() const
RCP< const LinearOpSourceBase< Scalar > > extract_fwdOpSrc()