45 #ifndef THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP
46 #define THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP
49 #include "Thyra_GeneralSolveCriteriaBelosStatusTest.hpp"
50 #include "Thyra_LinearOpWithSolveHelpers.hpp"
103 const std::string& residualScalingType)
132 if (solverValidParams->
isParameter (
"Implicit Residual Scaling")) {
133 solverParams->
set (
"Implicit Residual Scaling", residualScalingType);
135 if (solverValidParams->
isParameter (
"Explicit Residual Scaling")) {
136 solverParams->
set (
"Explicit Residual Scaling", residualScalingType);
149 template<
class Scalar>
151 :convergenceTestFrequency_(-1),
152 isExternalPrec_(false),
153 supportSolveUse_(SUPPORT_SOLVE_UNSPECIFIED),
159 template<
class Scalar>
161 const RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > &lp,
163 const RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > &iterativeSolver,
164 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
165 const RCP<
const PreconditionerBase<Scalar> > &prec,
166 const bool isExternalPrec_in,
167 const RCP<
const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
168 const ESupportSolveUse &supportSolveUse_in,
169 const int convergenceTestFrequency
177 this->setLinePrefix(
"BELOS/T");
180 solverPL_ = solverPL;
181 iterativeSolver_ = iterativeSolver;
182 fwdOpSrc_ = fwdOpSrc;
184 isExternalPrec_ = isExternalPrec_in;
185 approxFwdOpSrc_ = approxFwdOpSrc;
186 supportSolveUse_ = supportSolveUse_in;
187 convergenceTestFrequency_ = convergenceTestFrequency;
191 if (solverPL_->isParameter(
"Convergence Tolerance")) {
197 if (solverPL_->isType<
double> (
"Convergence Tolerance")) {
199 as<magnitude_type> (solverPL_->get<
double> (
"Convergence Tolerance"));
205 true, std::invalid_argument,
"BelosLinearOpWithSolve::initialize: "
206 "The \"Convergence Tolerance\" parameter, which you provided, must "
207 "have type double (the type of the magnitude of Scalar = double).");
209 else if (solverPL_->isType<magnitude_type> (
"Convergence Tolerance")) {
210 defaultTol_ = solverPL_->get<magnitude_type> (
"Convergence Tolerance");
219 true, InvalidParameterType,
"BelosLinearOpWithSolve::initialize: "
220 "The \"Convergence Tolerance\" parameter, which you provided, must "
221 "have type double (preferred) or the type of the magnitude of Scalar "
222 "= " << TypeNameTraits<Scalar>::name () <<
", which is " <<
223 TypeNameTraits<magnitude_type>::name () <<
" in this case. You can "
224 "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
228 if (solverPL_->isParameter(
"Timer Label") && solverPL_->isType<std::string>(
"Timer Label")) {
229 label_ = solverPL_->get<std::string>(
"Timer Label");
230 lp_->setLabel(label_);
235 iterativeSolver->getValidParameters();
241 if (defaultPL->
isType<
double> (
"Convergence Tolerance")) {
243 as<magnitude_type> (defaultPL->
get<
double> (
"Convergence Tolerance"));
249 true, std::invalid_argument,
"BelosLinearOpWithSolve::initialize: "
250 "The \"Convergence Tolerance\" parameter, which you provided, must "
251 "have type double (the type of the magnitude of Scalar = double).");
253 else if (defaultPL->
isType<magnitude_type> (
"Convergence Tolerance")) {
254 defaultTol_ = defaultPL->
get<magnitude_type> (
"Convergence Tolerance");
263 true, InvalidParameterType,
"BelosLinearOpWithSolve::initialize: "
264 "The \"Convergence Tolerance\" parameter, which you provided, must "
265 "have type double (preferred) or the type of the magnitude of Scalar "
266 "= " << TypeNameTraits<Scalar>::name () <<
", which is " <<
267 TypeNameTraits<magnitude_type>::name () <<
" in this case. You can "
268 "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
274 template<
class Scalar>
279 _fwdOpSrc = fwdOpSrc_;
285 template<
class Scalar>
296 template<
class Scalar>
299 return isExternalPrec_;
303 template<
class Scalar>
308 _approxFwdOpSrc = approxFwdOpSrc_;
310 return _approxFwdOpSrc;
314 template<
class Scalar>
317 return supportSolveUse_;
321 template<
class Scalar>
323 RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > *lp,
325 RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > *iterativeSolver,
326 RCP<
const LinearOpSourceBase<Scalar> > *fwdOpSrc,
327 RCP<
const PreconditionerBase<Scalar> > *prec,
328 bool *isExternalPrec_in,
329 RCP<
const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
330 ESupportSolveUse *supportSolveUse_in
334 if (solverPL) *solverPL = solverPL_;
335 if (iterativeSolver) *iterativeSolver = iterativeSolver_;
336 if (fwdOpSrc) *fwdOpSrc = fwdOpSrc_;
337 if (prec) *prec = prec_;
338 if (isExternalPrec_in) *isExternalPrec_in = isExternalPrec_;
339 if (approxFwdOpSrc) *approxFwdOpSrc = approxFwdOpSrc_;
340 if (supportSolveUse_in) *supportSolveUse_in = supportSolveUse_;
347 isExternalPrec_ =
false;
349 supportSolveUse_ = SUPPORT_SOLVE_UNSPECIFIED;
356 template<
class Scalar>
361 return lp_->getOperator()->range();
366 template<
class Scalar>
371 return lp_->getOperator()->domain();
376 template<
class Scalar>
387 template<
class Scalar>
390 std::ostringstream oss;
394 oss <<
"iterativeSolver=\'"<<iterativeSolver_->description()<<
"\'";
395 oss <<
",fwdOp=\'"<<lp_->getOperator()->description()<<
"\'";
396 if (lp_->getLeftPrec().get())
397 oss <<
",leftPrecOp=\'"<<lp_->getLeftPrec()->description()<<
"\'";
398 if (lp_->getRightPrec().get())
399 oss <<
",rightPrecOp=\'"<<lp_->getRightPrec()->description()<<
"\'";
408 template<
class Scalar>
416 using Teuchos::describe;
424 *out << this->description() << std::endl;
431 <<
"rangeDim=" << this->range()->dim()
432 <<
",domainDim=" << this->domain()->dim() <<
"}\n";
433 if (lp_->getOperator().get()) {
436 <<
"iterativeSolver = "<<describe(*iterativeSolver_,verbLevel)
437 <<
"fwdOp = " << describe(*lp_->getOperator(),verbLevel);
438 if (lp_->getLeftPrec().get())
439 *out <<
"leftPrecOp = "<<describe(*lp_->getLeftPrec(),verbLevel);
440 if (lp_->getRightPrec().get())
441 *out <<
"rightPrecOp = "<<describe(*lp_->getRightPrec(),verbLevel);
457 template<
class Scalar>
460 return ::Thyra::opSupported(*lp_->getOperator(),M_trans);
464 template<
class Scalar>
466 const EOpTransp M_trans,
467 const MultiVectorBase<Scalar> &X,
468 const Ptr<MultiVectorBase<Scalar> > &Y,
473 ::Thyra::apply<Scalar>(*lp_->getOperator(), M_trans, X, Y, alpha, beta);
480 template<
class Scalar>
488 template<
class Scalar>
491 const Ptr<
const SolveCriteria<Scalar> > )
const
494 if (real_trans(transp)==
NOTRANS)
return true;
512 template<
class Scalar>
515 EOpTransp M_trans,
const SolveMeasureType& solveMeasureType)
const
517 SolveCriteria<Scalar> solveCriteria(solveMeasureType, SolveCriteria<Scalar>::unspecifiedTolerance());
518 return solveSupportsNewImpl(M_trans, Teuchos::constOptInArg(solveCriteria));
522 template<
class Scalar>
525 const EOpTransp M_trans,
526 const MultiVectorBase<Scalar> &
B,
527 const Ptr<MultiVectorBase<Scalar> > &X,
528 const Ptr<
const SolveCriteria<Scalar> > solveCriteria
532 THYRA_FUNC_TIME_MONITOR(
"Stratimikos: BelosLOWS");
535 using Teuchos::rcpFromRef;
536 using Teuchos::rcpFromPtr;
540 using Teuchos::parameterList;
541 using Teuchos::describe;
543 typedef typename ST::magnitudeType ScalarMag;
545 totalTimer.start(
true);
547 assertSolveSupports(*
this, M_trans, solveCriteria);
553 OSTab tab = this->getOSTab();
555 *out <<
"\nStarting iterations with Belos:\n";
557 *out <<
"Using forward operator = " << describe(*fwdOpSrc_->getOp(),verbLevel);
558 *out <<
"Using iterative solver = " << describe(*iterativeSolver_,verbLevel);
559 *out <<
"With #Eqns="<<B.range()->dim()<<
", #RHSs="<<B.domain()->dim()<<
" ...\n";
566 bool ret = lp_->setProblem( rcpFromPtr(X), rcpFromRef(B) );
568 ret ==
false, CatastrophicSolveFailure
569 ,
"Error, the Belos::LinearProblem could not be set for the current solve!"
582 SolveMeasureType solveMeasureType;
585 solveMeasureType = solveCriteria->solveMeasureType;
586 const ScalarMag requestedTol = solveCriteria->requestedTol;
587 if (solveMeasureType.useDefault()) {
588 tmpPL->
set(
"Convergence Tolerance", defaultTol_);
590 else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_RHS)) {
591 if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
592 tmpPL->
set(
"Convergence Tolerance", requestedTol);
595 tmpPL->
set(
"Convergence Tolerance", defaultTol_);
597 setResidualScalingType (tmpPL, validPL,
"Norm of RHS");
599 else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_INIT_RESIDUAL)) {
600 if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
601 tmpPL->
set(
"Convergence Tolerance", requestedTol);
604 tmpPL->
set(
"Convergence Tolerance", defaultTol_);
606 setResidualScalingType (tmpPL, validPL,
"Norm of Initial Residual");
610 generalSolveCriteriaBelosStatusTest = createGeneralSolveCriteriaBelosStatusTest(
611 *solveCriteria, convergenceTestFrequency_);
613 generalSolveCriteriaBelosStatusTest->setOStream(out);
614 generalSolveCriteriaBelosStatusTest->setVerbLevel(
incrVerbLevel(verbLevel, -1));
617 tmpPL->
set(
"Convergence Tolerance", 1.0);
620 if (
nonnull(solveCriteria->extraParameters)) {
621 if (Teuchos::isParameterType<int>(*solveCriteria->extraParameters,
"Maximum Iterations")) {
622 tmpPL->
set(
"Maximum Iterations", Teuchos::get<int>(*solveCriteria->extraParameters,
"Maximum Iterations"));
628 validPL->
isParameter (
"Implicit Residual Scaling"))
629 tmpPL->
set(
"Implicit Residual Scaling",
630 "Norm of Preconditioned Initial Residual");
634 tmpPL->
set(
"Convergence Tolerance", defaultTol_);
641 Belos::ReturnType belosSolveStatus;
651 tmpPL->
set(
"Output Stream", outUsed);
653 if (
nonnull(generalSolveCriteriaBelosStatusTest)) {
654 iterativeSolver_->setUserConvStatusTest(generalSolveCriteriaBelosStatusTest);
657 belosSolveStatus = iterativeSolver_->solve();
659 catch (Belos::BelosError& ex) {
661 CatastrophicSolveFailure,
672 SolveStatus<Scalar> solveStatus;
674 switch (belosSolveStatus) {
675 case Belos::Unconverged: {
676 solveStatus.solveStatus = SOLVE_STATUS_UNCONVERGED;
687 solveStatus.achievedTol = iterativeSolver_->achievedTol();
688 }
catch (std::runtime_error&) {
693 case Belos::Converged: {
694 solveStatus.solveStatus = SOLVE_STATUS_CONVERGED;
695 if (
nonnull(generalSolveCriteriaBelosStatusTest)) {
700 const ArrayView<const ScalarMag> achievedTol =
701 generalSolveCriteriaBelosStatusTest->achievedTol();
702 solveStatus.achievedTol = ST::zero();
703 for (
Ordinal i = 0; i < achievedTol.size(); ++i) {
704 solveStatus.achievedTol = std::max(solveStatus.achievedTol, achievedTol[i]);
711 solveStatus.achievedTol = iterativeSolver_->achievedTol();
712 }
catch (std::runtime_error&) {
715 solveStatus.achievedTol = tmpPL->
get(
"Convergence Tolerance", defaultTol_);
723 std::ostringstream ossmessage;
725 <<
"The Belos solver " << (label_ !=
"" ? (
"\"" + label_ +
"\" ") :
"")
726 <<
"of type \""<<iterativeSolver_->description()
727 <<
"\" returned a solve status of \""<<
toString(solveStatus.solveStatus) <<
"\""
728 <<
" in " << iterativeSolver_->getNumIters() <<
" iterations"
729 <<
" with total CPU time of " << totalTimer.totalElapsedTime() <<
" sec" ;
731 *out <<
"\n" << ossmessage.str() <<
"\n";
733 solveStatus.message = ossmessage.str();
738 if (solveStatus.extraParameters.is_null()) {
739 solveStatus.extraParameters = parameterList ();
741 solveStatus.extraParameters->set (
"Belos/Iteration Count",
742 iterativeSolver_->getNumIters());\
744 solveStatus.extraParameters->set (
"Iteration Count",
745 iterativeSolver_->getNumIters());\
752 solveStatus.extraParameters->set (
"Belos/Achieved Tolerance",
753 solveStatus.achievedTol);
768 #endif // THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP
virtual bool solveSupportsImpl(EOpTransp M_trans) const
RCP< const VectorSpaceBase< Scalar > > domain() const
const std::string toString(const EAdjointEpetraOp adjointEpetraOp)
bool is_null(const boost::shared_ptr< T > &p)
virtual SolveStatus< Scalar > solveImpl(const EOpTransp transp, const MultiVectorBase< Scalar > &B, const Ptr< MultiVectorBase< Scalar > > &X, const Ptr< const SolveCriteria< Scalar > > solveCriteria) const
RCP< const LinearOpSourceBase< Scalar > > extract_approxFwdOpSrc()
T & get(ParameterList &l, const std::string &name)
void initialize(const RCP< Belos::LinearProblem< Scalar, MV_t, LO_t > > &lp, const RCP< Teuchos::ParameterList > &solverPL, const RCP< Belos::SolverManager< Scalar, MV_t, LO_t > > &iterativeSolver, const RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const RCP< const PreconditionerBase< Scalar > > &prec, const bool isExternalPrec, const RCP< const LinearOpSourceBase< Scalar > > &approxFwdOpSrc, const ESupportSolveUse &supportSolveUse, const int convergenceTestFrequency)
Initializes given precreated solver objects.
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
bool nonnull(const std::shared_ptr< T > &p)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
RCP< const PreconditionerBase< Scalar > > extract_prec()
virtual bool solveSupportsNewImpl(EOpTransp transp, const Ptr< const SolveCriteria< Scalar > > solveCriteria) const
RCP< const VectorSpaceBase< Scalar > > range() const
virtual void applyImpl(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
bool isExternalPrec() const
basic_OSTab< char > OSTab
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
bool isParameter(const std::string &name) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
#define TEUCHOS_SWITCH_DEFAULT_DEBUG_ASSERT()
virtual std::string description() const
std::string description() const
TypeTo as(const TypeFrom &t)
ESupportSolveUse supportSolveUse() const
ParameterList & setParameters(const ParameterList &source)
RCP< const LinearOpSourceBase< Scalar > > extract_fwdOpSrc()
TEUCHOSCORE_LIB_DLL_EXPORT EVerbosityLevel incrVerbLevel(const EVerbosityLevel inputVerbLevel, const int numLevels)
bool nonnull(const boost::shared_ptr< T > &p)
basic_FancyOStream< char > FancyOStream
bool isType(const std::string &name) const
BelosLinearOpWithSolve()
Construct to unintialize.
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.
RCP< const LinearOpBase< Scalar > > clone() const
virtual bool opSupportedImpl(EOpTransp M_trans) const
virtual bool solveSupportsSolveMeasureTypeImpl(EOpTransp M_trans, const SolveMeasureType &solveMeasureType) const
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)