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>
78 const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverType_name =
"Solver Type";
79 template<
class Scalar>
80 const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverType_default =
"Pseudo Block GMRES";
81 template<
class Scalar>
82 const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverTypes_name =
"Solver Types";
83 template<
class Scalar>
84 const std::string BelosLinearOpWithSolveFactory<Scalar>::BlockGMRES_name =
"Block GMRES";
85 template<
class Scalar>
86 const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockGMRES_name =
"Pseudo Block GMRES";
87 template<
class Scalar>
88 const std::string BelosLinearOpWithSolveFactory<Scalar>::BlockCG_name =
"Block CG";
89 template<
class Scalar>
90 const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockCG_name =
"Pseudo Block CG";
91 template<
class Scalar>
92 const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockStochasticCG_name =
"Pseudo Block Stochastic CG";
93 template<
class Scalar>
94 const std::string BelosLinearOpWithSolveFactory<Scalar>::GCRODR_name =
"GCRODR";
95 template<
class Scalar>
96 const std::string BelosLinearOpWithSolveFactory<Scalar>::RCG_name =
"RCG";
97 template<
class Scalar>
98 const std::string BelosLinearOpWithSolveFactory<Scalar>::MINRES_name =
"MINRES";
99 template<
class Scalar>
100 const std::string BelosLinearOpWithSolveFactory<Scalar>::TFQMR_name =
"TFQMR";
101 template<
class Scalar>
102 const std::string BelosLinearOpWithSolveFactory<Scalar>::ConvergenceTestFrequency_name =
"Convergence Test Frequency";
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>
124 :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES)
133 template<
class Scalar>
140 template<
class Scalar>
143 const std::string &precFactoryName
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>
168 template<
class Scalar>
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>
187 if(precFactory_.get())
188 return precFactory_->isCompatible(fwdOpSrc);
193 template<
class Scalar>
201 template<
class Scalar>
209 initializeOpImpl(fwdOpSrc,null,null,
false,Op,supportSolveUse);
213 template<
class Scalar>
224 template<
class Scalar>
229 if(precFactory_.get())
235 template<
class Scalar>
244 initializeOpImpl(fwdOpSrc,null,prec,
false,Op,supportSolveUse);
248 template<
class Scalar>
257 initializeOpImpl(fwdOpSrc,approxFwdOpSrc,null,
false,Op,supportSolveUse);
261 template<
class 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>
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>
320 template<
class Scalar>
325 paramList_ = Teuchos::null;
330 template<
class Scalar>
338 template<
class Scalar>
342 return thisValidParamList_;
349 template<
class Scalar>
352 std::ostringstream oss;
353 oss <<
"Thyra::BelosLinearOpWithSolveFactory";
364 template<
class Scalar>
369 using Teuchos::tuple;
370 using Teuchos::setStringToIntegralParameter;
376 EBelosSolverType> >::getDummyObject());
381 if(validParamList.
get()==NULL) {
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.");
461 &solverTypesSL = validParamList->
sublist(SolverTypes_name);
477 *mgr.getValidParameters()
483 *mgr.getValidParameters()
495 *mgr.getValidParameters()
501 *mgr.getValidParameters()
517 return validParamList;
521 template<
class Scalar>
522 void BelosLinearOpWithSolveFactory<Scalar>::updateThisValidParamList()
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,
543 using Teuchos::set_extra_data;
545 typedef MultiVectorBase<Scalar> MV_t;
546 typedef LinearOpBase<Scalar> LO_t;
548 const RCP<Teuchos::FancyOStream> out = this->getOStream();
551 if(out.get() &&
static_cast<int>(verbLevel) > static_cast<int>(
Teuchos::VERB_LOW))
552 *out <<
"\nEntering Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<
">::initializeOpImpl(...) ...\n";
562 RCP<const LinearOpBase<Scalar> >
563 fwdOp = fwdOpSrc->getOp(),
564 approxFwdOp = ( approxFwdOpSrc.get() ? approxFwdOpSrc->getOp() : Teuchos::null );
570 BelosLinearOpWithSolve<Scalar>
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;
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();
658 !( left.get() || right.get() || unspecified.get() ), std::logic_error
659 ,
"Error, at least one preconditoner linear operator objects must be set!"
662 if (paramList_->get<
bool>(LeftPreconditionerIfUnspecified_name,
false))
663 lp->setLeftPrec(unspecified);
665 lp->setRightPrec(unspecified);
668 lp->setLeftPrec(left);
671 lp->setRightPrec(right);
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;
698 switch(solverType_) {
699 case SOLVER_TYPE_BLOCK_GMRES:
702 if(paramList_.get()) {
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()) {
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()) {
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()) {
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()) {
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()) {
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()) {
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()) {
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()) {
866 if (oldIterSolver != Teuchos::null) {
867 iterativeSolver = oldIterSolver;
868 iterativeSolver->setProblem( lp );
869 iterativeSolver->setParameters( solverPL );
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
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.
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#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
bool isExternalPrec() const
void unsetPreconditionerFactory(Teuchos::RCP< PreconditionerFactoryBase< Scalar > > *precFactory, std::string *precFactoryName)
Thyra specializations of MultiVecTraits and OperatorTraits.
void setParameterList(Teuchos::RCP< Teuchos::ParameterList > const ¶mList)
bool isCompatible(const LinearOpSourceBase< Scalar > &fwdOpSrc) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
bool acceptsPreconditionerFactory() const
Returns true .
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
void initializeOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
std::string description() 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.
ESupportSolveUse supportSolveUse() const
ParameterList & setParameters(const ParameterList &source)
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
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)