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"
64 #include "Teuchos_VerboseObjectParameterListHelpers.hpp"
65 #include "Teuchos_StandardParameterEntryValidators.hpp"
66 #include "Teuchos_ParameterList.hpp"
67 #include "Teuchos_dyn_cast.hpp"
68 #include "Teuchos_ValidatorXMLConverterDB.hpp"
69 #include "Teuchos_StandardValidatorXMLConverters.hpp"
77 template<
class Scalar>
79 template<
class Scalar>
81 template<
class Scalar>
83 template<
class Scalar>
85 template<
class Scalar>
87 template<
class Scalar>
89 template<
class Scalar>
91 template<
class Scalar>
93 template<
class Scalar>
95 template<
class Scalar>
97 template<
class Scalar>
99 template<
class Scalar>
101 template<
class Scalar>
105 const std::string LeftPreconditionerIfUnspecified_name =
"Left Preconditioner If Unspecified";
111 template<
class Scalar>
113 :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES),
114 convergenceTestFrequency_(1)
116 updateThisValidParamList();
120 template<
class Scalar>
122 const RCP<PreconditionerFactoryBase<Scalar> > &precFactory
124 :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES)
133 template<
class Scalar>
140 template<
class Scalar>
142 const RCP<PreconditionerFactoryBase<Scalar> > &precFactory,
143 const std::string &precFactoryName
146 TEUCHOS_TEST_FOR_EXCEPT(!precFactory.get());
147 RCP<const Teuchos::ParameterList>
148 precFactoryValidPL = precFactory->getValidParameters();
149 const std::string _precFactoryName =
150 ( precFactoryName !=
""
152 : ( precFactoryValidPL.get() ? precFactoryValidPL->name() :
"GENERIC PRECONDITIONER FACTORY" )
154 precFactory_ = precFactory;
155 precFactoryName_ = _precFactoryName;
156 updateThisValidParamList();
160 template<
class Scalar>
161 RCP<PreconditionerFactoryBase<Scalar> >
168 template<
class Scalar>
170 RCP<PreconditionerFactoryBase<Scalar> > *precFactory,
171 std::string *precFactoryName
174 if(precFactory) *precFactory = precFactory_;
175 if(precFactoryName) *precFactoryName = precFactoryName_;
176 precFactory_ = Teuchos::null;
177 precFactoryName_ =
"";
178 updateThisValidParamList();
182 template<
class Scalar>
184 const LinearOpSourceBase<Scalar> &fwdOpSrc
187 if(precFactory_.get())
188 return precFactory_->isCompatible(fwdOpSrc);
193 template<
class Scalar>
194 RCP<LinearOpWithSolveBase<Scalar> >
201 template<
class Scalar>
203 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
204 LinearOpWithSolveBase<Scalar> *Op,
205 const ESupportSolveUse supportSolveUse
209 initializeOpImpl(fwdOpSrc,
null,
null,
false,Op,supportSolveUse);
213 template<
class Scalar>
215 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
216 LinearOpWithSolveBase<Scalar> *Op
220 initializeOpImpl(fwdOpSrc,
null,
null,
true,Op,SUPPORT_SOLVE_UNSPECIFIED);
224 template<
class Scalar>
226 const EPreconditionerInputType precOpType
229 if(precFactory_.get())
231 return (precOpType==PRECONDITIONER_INPUT_TYPE_AS_OPERATOR);
235 template<
class Scalar>
237 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
238 const RCP<
const PreconditionerBase<Scalar> > &prec,
239 LinearOpWithSolveBase<Scalar> *Op,
240 const ESupportSolveUse supportSolveUse
244 initializeOpImpl(fwdOpSrc,
null,prec,
false,Op,supportSolveUse);
248 template<
class Scalar>
250 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
251 const RCP<
const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
252 LinearOpWithSolveBase<Scalar> *Op,
253 const ESupportSolveUse supportSolveUse
257 initializeOpImpl(fwdOpSrc,approxFwdOpSrc,
null,
false,Op,supportSolveUse);
261 template<
class Scalar>
263 LinearOpWithSolveBase<Scalar> *Op,
264 RCP<
const LinearOpSourceBase<Scalar> > *fwdOpSrc,
265 RCP<
const PreconditionerBase<Scalar> > *prec,
266 RCP<
const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
267 ESupportSolveUse *supportSolveUse
271 TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
275 RCP<const LinearOpSourceBase<Scalar> >
277 RCP<const PreconditionerBase<Scalar> >
282 RCP<const LinearOpSourceBase<Scalar> >
286 if(fwdOpSrc) *fwdOpSrc = _fwdOpSrc;
287 if(prec) *prec = _prec;
288 if(approxFwdOpSrc) *approxFwdOpSrc = _approxFwdOpSrc;
289 if(supportSolveUse) *supportSolveUse = _supportSolveUse;
296 template<
class Scalar>
298 RCP<Teuchos::ParameterList>
const& paramList
301 TEUCHOS_TEST_FOR_EXCEPT(paramList.get()==NULL);
302 paramList->validateParametersAndSetDefaults(*this->getValidParameters(), 1);
303 paramList_ = paramList;
305 Teuchos::getIntegralValue<EBelosSolverType>(*paramList_, SolverType_name);
306 convergenceTestFrequency_ =
307 Teuchos::getParameter<int>(*paramList_, ConvergenceTestFrequency_name);
308 Teuchos::readVerboseObjectSublist(&*paramList_,
this);
312 template<
class Scalar>
313 RCP<Teuchos::ParameterList>
320 template<
class Scalar>
321 RCP<Teuchos::ParameterList>
324 RCP<Teuchos::ParameterList> _paramList = paramList_;
325 paramList_ = Teuchos::null;
330 template<
class Scalar>
331 RCP<const Teuchos::ParameterList>
338 template<
class Scalar>
339 RCP<const Teuchos::ParameterList>
342 return thisValidParamList_;
349 template<
class Scalar>
352 std::ostringstream oss;
353 oss <<
"Thyra::BelosLinearOpWithSolveFactory";
364 template<
class Scalar>
365 RCP<const Teuchos::ParameterList>
369 using Teuchos::tuple;
370 using Teuchos::setStringToIntegralParameter;
371 Teuchos::ValidatorXMLConverterDB::addConverter(
372 Teuchos::DummyObjectGetter<
373 Teuchos::StringToIntegralParameterEntryValidator<EBelosSolverType>
375 Teuchos::DummyObjectGetter<Teuchos::StringToIntegralValidatorXMLConverter<
376 EBelosSolverType> >::getDummyObject());
378 typedef MultiVectorBase<Scalar> MV_t;
379 typedef LinearOpBase<Scalar> LO_t;
380 static RCP<Teuchos::ParameterList> validParamList;
381 if(validParamList.get()==NULL) {
382 validParamList = Teuchos::rcp(
new Teuchos::ParameterList(
"BelosLinearOpWithSolveFactory"));
383 setStringToIntegralParameter<EBelosSolverType>(
384 SolverType_name, SolverType_default,
385 "Type of linear solver algorithm to use.",
388 "Pseudo Block GMRES",
391 "Pseudo Block Stochastic CG",
398 "Block GMRES solver for nonsymmetric linear systems. It can also solve "
399 "single right-hand side systems, and can also perform Flexible GMRES "
400 "(where the preconditioner may change at every iteration, for example "
401 "for inner-outer iterations), by setting options in the \"Block GMRES\" "
404 "GMRES solver for nonsymmetric linear systems, that performs single "
405 "right-hand side solves on multiple right-hand sides at once. It "
406 "exploits operator multivector multiplication in order to amortize "
407 "global communication costs. Individual linear systems are deflated "
408 "out as they are solved.",
410 "Block CG solver for symmetric (Hermitian in complex arithmetic) "
411 "positive definite linear systems. It can also solve single "
412 "right-hand-side systems.",
414 "CG solver that performs single right-hand side CG on multiple right-hand "
415 "sides at once. It exploits operator multivector multiplication in order "
416 "to amortize global communication costs. Individual linear systems are "
417 "deflated out as they are solved.",
419 "stochastic CG solver that performs single right-hand side CG on multiple right-hand "
420 "sides at once. It exploits operator multivector multiplication in order "
421 "to amortize global communication costs. Individual linear systems are "
422 "deflated out as they are solved. [EXPERIMENTAL]",
424 "Variant of GMRES that performs subspace recycling to accelerate "
425 "convergence for sequences of solves with related linear systems. "
426 "Individual linear systems are deflated out as they are solved. "
427 "The current implementation only supports real-valued Scalar types.",
429 "CG solver for symmetric (Hermitian in complex arithmetic) positive "
430 "definite linear systems, that performs subspace recycling to "
431 "accelerate convergence for sequences of related linear systems.",
433 "MINRES solver for symmetric indefinite linear systems. It performs "
434 "single-right-hand-side solves on multiple right-hand sides sequentially.",
436 "TFQMR (Transpose-Free QMR) solver for nonsymmetric linear systems."
438 tuple<EBelosSolverType>(
439 SOLVER_TYPE_BLOCK_GMRES,
440 SOLVER_TYPE_PSEUDO_BLOCK_GMRES,
441 SOLVER_TYPE_BLOCK_CG,
442 SOLVER_TYPE_PSEUDO_BLOCK_CG,
443 SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG,
451 validParamList->set(ConvergenceTestFrequency_name, as<int>(1),
452 "Number of linear solver iterations to skip between applying"
453 " user-defined convergence test.");
455 LeftPreconditionerIfUnspecified_name,
false,
456 "If the preconditioner does not specify if it is left or right, and this\n"
457 "option is set to true, put the preconditioner on the left side.\n"
458 "Historically, preconditioning is on the right. Some solvers may not\n"
459 "support left preconditioning.");
460 Teuchos::ParameterList
461 &solverTypesSL = validParamList->sublist(SolverTypes_name);
477 *mgr.getValidParameters()
482 solverTypesSL.sublist(PseudoBlockCG_name).setParameters(
483 *mgr.getValidParameters()
488 solverTypesSL.sublist(PseudoBlockStochasticCG_name).setParameters(
495 *mgr.getValidParameters()
500 solverTypesSL.sublist(RCG_name).setParameters(
501 *mgr.getValidParameters()
506 solverTypesSL.sublist(MINRES_name).setParameters(
517 return validParamList;
521 template<
class Scalar>
522 void BelosLinearOpWithSolveFactory<Scalar>::updateThisValidParamList()
524 thisValidParamList_ = Teuchos::rcp(
525 new Teuchos::ParameterList(*generateAndGetValidParameters())
527 Teuchos::setupVerboseObjectSublist(&*thisValidParamList_);
531 template<
class Scalar>
532 void BelosLinearOpWithSolveFactory<Scalar>::initializeOpImpl(
533 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
534 const RCP<
const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
535 const RCP<
const PreconditionerBase<Scalar> > &prec_in,
536 const bool reusePrec,
537 LinearOpWithSolveBase<Scalar> *Op,
538 const ESupportSolveUse supportSolveUse
543 using Teuchos::set_extra_data;
544 typedef Teuchos::ScalarTraits<Scalar> ST;
545 typedef MultiVectorBase<Scalar> MV_t;
546 typedef LinearOpBase<Scalar> LO_t;
548 const RCP<Teuchos::FancyOStream> out = this->getOStream();
549 const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
550 Teuchos::OSTab tab(out);
551 if(out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
552 *out <<
"\nEntering Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<
">::initializeOpImpl(...) ...\n";
559 TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
560 TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc.get()==NULL);
561 TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc->getOp().get()==NULL);
562 RCP<const LinearOpBase<Scalar> >
563 fwdOp = fwdOpSrc->getOp(),
564 approxFwdOp = ( approxFwdOpSrc.get() ? approxFwdOpSrc->getOp() : Teuchos::null );
570 BelosLinearOpWithSolve<Scalar>
571 *belosOp = &Teuchos::dyn_cast<BelosLinearOpWithSolve<Scalar> >(*Op);
577 RCP<PreconditionerBase<Scalar> > myPrec = Teuchos::null;
578 RCP<const PreconditionerBase<Scalar> > prec = Teuchos::null;
585 if(precFactory_.get()) {
587 ( !belosOp->isExternalPrec()
588 ? Teuchos::rcp_const_cast<PreconditionerBase<Scalar> >(belosOp->extract_prec())
591 bool hasExistingPrec =
false;
593 hasExistingPrec =
true;
598 hasExistingPrec =
false;
599 myPrec = precFactory_->createPrec();
601 if( hasExistingPrec && reusePrec ) {
606 if(approxFwdOp.get())
607 precFactory_->initializePrec(approxFwdOpSrc,&*myPrec);
609 precFactory_->initializePrec(fwdOpSrc,&*myPrec);
619 bool oldIsExternalPrec =
false;
620 RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > oldLP = Teuchos::null;
621 RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > oldIterSolver = Teuchos::null;
622 RCP<const LinearOpSourceBase<Scalar> > oldFwdOpSrc = Teuchos::null;
623 RCP<const LinearOpSourceBase<Scalar> > oldApproxFwdOpSrc = Teuchos::null;
624 ESupportSolveUse oldSupportSolveUse = SUPPORT_SOLVE_UNSPECIFIED;
626 belosOp->uninitialize( &oldLP, NULL, &oldIterSolver, &oldFwdOpSrc,
627 NULL, &oldIsExternalPrec, &oldApproxFwdOpSrc, &oldSupportSolveUse );
636 if (oldLP != Teuchos::null) {
640 lp = rcp(
new LP_t());
647 lp->setOperator(fwdOp);
654 RCP<const LinearOpBase<Scalar> > unspecified = prec->getUnspecifiedPrecOp();
655 RCP<const LinearOpBase<Scalar> > left = prec->getLeftPrecOp();
656 RCP<const LinearOpBase<Scalar> > right = prec->getRightPrecOp();
657 TEUCHOS_TEST_FOR_EXCEPTION(
658 !( left.get() || right.get() || unspecified.get() ), std::logic_error
659 ,
"Error, at least one preconditoner linear operator objects must be set!"
661 if (nonnull(unspecified)) {
662 if (paramList_->get<
bool>(LeftPreconditionerIfUnspecified_name,
false))
663 lp->setLeftPrec(unspecified);
665 lp->setRightPrec(unspecified);
667 else if (nonnull(left)) {
668 lp->setLeftPrec(left);
670 else if (nonnull(right)) {
671 lp->setRightPrec(right);
675 TEUCHOS_TEST_FOR_EXCEPTION(
676 nonnull(left) && nonnull(right),std::logic_error
677 ,
"Error, we can not currently handle split preconditioners!"
682 set_extra_data<RCP<PreconditionerBase<Scalar> > >(myPrec,
"Belos::InternalPrec",
683 Teuchos::inOutArg(lp), Teuchos::POST_DESTROY,
false);
685 else if(prec.get()) {
686 set_extra_data<RCP<const PreconditionerBase<Scalar> > >(prec,
"Belos::ExternalPrec",
687 Teuchos::inOutArg(lp), Teuchos::POST_DESTROY,
false);
695 RCP<IterativeSolver_t> iterativeSolver = Teuchos::null;
696 RCP<Teuchos::ParameterList> solverPL = Teuchos::rcp(
new Teuchos::ParameterList() );
698 switch(solverType_) {
699 case SOLVER_TYPE_BLOCK_GMRES:
702 if(paramList_.get()) {
703 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
704 Teuchos::ParameterList &gmresPL = solverTypesPL.sublist(BlockGMRES_name);
705 solverPL = Teuchos::rcp( &gmresPL,
false );
708 if (oldIterSolver != Teuchos::null) {
709 iterativeSolver = oldIterSolver;
710 iterativeSolver->setProblem( lp );
711 iterativeSolver->setParameters( solverPL );
718 case SOLVER_TYPE_PSEUDO_BLOCK_GMRES:
721 if(paramList_.get()) {
722 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
723 Teuchos::ParameterList &pbgmresPL = solverTypesPL.sublist(PseudoBlockGMRES_name);
724 solverPL = Teuchos::rcp( &pbgmresPL,
false );
729 if (oldIterSolver != Teuchos::null) {
730 iterativeSolver = oldIterSolver;
731 iterativeSolver->setProblem( lp );
732 iterativeSolver->setParameters( solverPL );
739 case SOLVER_TYPE_BLOCK_CG:
742 if(paramList_.get()) {
743 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
744 Teuchos::ParameterList &cgPL = solverTypesPL.sublist(BlockCG_name);
745 solverPL = Teuchos::rcp( &cgPL,
false );
748 if (oldIterSolver != Teuchos::null) {
749 iterativeSolver = oldIterSolver;
750 iterativeSolver->setProblem( lp );
751 iterativeSolver->setParameters( solverPL );
758 case SOLVER_TYPE_PSEUDO_BLOCK_CG:
761 if(paramList_.get()) {
762 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
763 Teuchos::ParameterList &pbcgPL = solverTypesPL.sublist(PseudoBlockCG_name);
764 solverPL = Teuchos::rcp( &pbcgPL,
false );
769 if (oldIterSolver != Teuchos::null) {
770 iterativeSolver = oldIterSolver;
771 iterativeSolver->setProblem( lp );
772 iterativeSolver->setParameters( solverPL );
779 case SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG:
782 if(paramList_.get()) {
783 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
784 Teuchos::ParameterList &pbcgPL = solverTypesPL.sublist(PseudoBlockStochasticCG_name);
785 solverPL = Teuchos::rcp( &pbcgPL,
false );
790 if (oldIterSolver != Teuchos::null) {
791 iterativeSolver = oldIterSolver;
792 iterativeSolver->setProblem( lp );
793 iterativeSolver->setParameters( solverPL );
800 case SOLVER_TYPE_GCRODR:
803 if(paramList_.get()) {
804 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
805 Teuchos::ParameterList &gcrodrPL = solverTypesPL.sublist(GCRODR_name);
806 solverPL = Teuchos::rcp( &gcrodrPL,
false );
809 if (oldIterSolver != Teuchos::null) {
810 iterativeSolver = oldIterSolver;
811 iterativeSolver->setProblem( lp );
812 iterativeSolver->setParameters( solverPL );
819 case SOLVER_TYPE_RCG:
822 if(paramList_.get()) {
823 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
824 Teuchos::ParameterList &rcgPL = solverTypesPL.sublist(RCG_name);
825 solverPL = Teuchos::rcp( &rcgPL,
false );
828 if (oldIterSolver != Teuchos::null) {
829 iterativeSolver = oldIterSolver;
830 iterativeSolver->setProblem( lp );
831 iterativeSolver->setParameters( solverPL );
838 case SOLVER_TYPE_MINRES:
841 if(paramList_.get()) {
842 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
843 Teuchos::ParameterList &minresPL = solverTypesPL.sublist(MINRES_name);
844 solverPL = Teuchos::rcp( &minresPL,
false );
847 if (oldIterSolver != Teuchos::null) {
848 iterativeSolver = oldIterSolver;
849 iterativeSolver->setProblem( lp );
850 iterativeSolver->setParameters( solverPL );
857 case SOLVER_TYPE_TFQMR:
860 if(paramList_.get()) {
861 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
862 Teuchos::ParameterList &minresPL = solverTypesPL.sublist(TFQMR_name);
863 solverPL = Teuchos::rcp( &minresPL,
false );
866 if (oldIterSolver != Teuchos::null) {
867 iterativeSolver = oldIterSolver;
868 iterativeSolver->setProblem( lp );
869 iterativeSolver->setParameters( solverPL );
879 TEUCHOS_TEST_FOR_EXCEPT(
true);
888 lp, solverPL, iterativeSolver,
889 fwdOpSrc, prec, myPrec.get()==NULL, approxFwdOpSrc,
890 supportSolveUse, convergenceTestFrequency_
892 belosOp->setOStream(out);
893 belosOp->setVerbLevel(verbLevel);
895 if(paramList_.get()) {
897 paramList_->validateParameters(*this->getValidParameters(),1);
900 if(out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
901 *out <<
"\nLeaving Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<
">::initializeOpImpl(...) ...\n";
909 #endif // THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP