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),
158 template<
class Scalar>
160 const RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > &lp,
162 const RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > &iterativeSolver,
163 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
164 const RCP<
const PreconditionerBase<Scalar> > &prec,
165 const bool isExternalPrec_in,
166 const RCP<
const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
167 const ESupportSolveUse &supportSolveUse_in,
168 const int convergenceTestFrequency
176 this->setLinePrefix(
"BELOS/T");
179 solverPL_ = solverPL;
180 iterativeSolver_ = iterativeSolver;
181 fwdOpSrc_ = fwdOpSrc;
183 isExternalPrec_ = isExternalPrec_in;
184 approxFwdOpSrc_ = approxFwdOpSrc;
185 supportSolveUse_ = supportSolveUse_in;
186 convergenceTestFrequency_ = convergenceTestFrequency;
190 if (solverPL_->isParameter(
"Convergence Tolerance")) {
196 if (solverPL_->isType<
double> (
"Convergence Tolerance")) {
198 as<magnitude_type> (solverPL_->get<
double> (
"Convergence Tolerance"));
204 true, std::invalid_argument,
"BelosLinearOpWithSolve::initialize: "
205 "The \"Convergence Tolerance\" parameter, which you provided, must "
206 "have type double (the type of the magnitude of Scalar = double).");
208 else if (solverPL_->isType<magnitude_type> (
"Convergence Tolerance")) {
209 defaultTol_ = solverPL_->get<magnitude_type> (
"Convergence Tolerance");
218 true, InvalidParameterType,
"BelosLinearOpWithSolve::initialize: "
219 "The \"Convergence Tolerance\" parameter, which you provided, must "
220 "have type double (preferred) or the type of the magnitude of Scalar "
221 "= " << TypeNameTraits<Scalar>::name () <<
", which is " <<
222 TypeNameTraits<magnitude_type>::name () <<
" in this case. You can "
223 "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
227 if (solverPL_->isParameter(
"Timer Label") && solverPL_->isType<std::string>(
"Timer Label"))
228 lp_->setLabel(solverPL_->get<std::string>(
"Timer Label"));
232 iterativeSolver->getValidParameters();
238 if (defaultPL->
isType<
double> (
"Convergence Tolerance")) {
240 as<magnitude_type> (defaultPL->
get<
double> (
"Convergence Tolerance"));
246 true, std::invalid_argument,
"BelosLinearOpWithSolve::initialize: "
247 "The \"Convergence Tolerance\" parameter, which you provided, must "
248 "have type double (the type of the magnitude of Scalar = double).");
250 else if (defaultPL->
isType<magnitude_type> (
"Convergence Tolerance")) {
251 defaultTol_ = defaultPL->
get<magnitude_type> (
"Convergence Tolerance");
260 true, InvalidParameterType,
"BelosLinearOpWithSolve::initialize: "
261 "The \"Convergence Tolerance\" parameter, which you provided, must "
262 "have type double (preferred) or the type of the magnitude of Scalar "
263 "= " << TypeNameTraits<Scalar>::name () <<
", which is " <<
264 TypeNameTraits<magnitude_type>::name () <<
" in this case. You can "
265 "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
271 template<
class Scalar>
276 _fwdOpSrc = fwdOpSrc_;
277 fwdOpSrc_ = Teuchos::null;
282 template<
class Scalar>
288 prec_ = Teuchos::null;
293 template<
class Scalar>
296 return isExternalPrec_;
300 template<
class Scalar>
305 _approxFwdOpSrc = approxFwdOpSrc_;
306 approxFwdOpSrc_ = Teuchos::null;
307 return _approxFwdOpSrc;
311 template<
class Scalar>
314 return supportSolveUse_;
318 template<
class Scalar>
320 RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > *lp,
322 RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > *iterativeSolver,
323 RCP<
const LinearOpSourceBase<Scalar> > *fwdOpSrc,
324 RCP<
const PreconditionerBase<Scalar> > *prec,
325 bool *isExternalPrec_in,
326 RCP<
const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
327 ESupportSolveUse *supportSolveUse_in
331 if (solverPL) *solverPL = solverPL_;
332 if (iterativeSolver) *iterativeSolver = iterativeSolver_;
333 if (fwdOpSrc) *fwdOpSrc = fwdOpSrc_;
334 if (prec) *prec = prec_;
335 if (isExternalPrec_in) *isExternalPrec_in = isExternalPrec_;
336 if (approxFwdOpSrc) *approxFwdOpSrc = approxFwdOpSrc_;
337 if (supportSolveUse_in) *supportSolveUse_in = supportSolveUse_;
340 solverPL_ = Teuchos::null;
341 iterativeSolver_ = Teuchos::null;
342 fwdOpSrc_ = Teuchos::null;
343 prec_ = Teuchos::null;
344 isExternalPrec_ =
false;
345 approxFwdOpSrc_ = Teuchos::null;
346 supportSolveUse_ = SUPPORT_SOLVE_UNSPECIFIED;
353 template<
class Scalar>
358 return lp_->getOperator()->range();
359 return Teuchos::null;
363 template<
class Scalar>
368 return lp_->getOperator()->domain();
369 return Teuchos::null;
373 template<
class Scalar>
377 return Teuchos::null;
384 template<
class Scalar>
387 std::ostringstream oss;
391 oss <<
"iterativeSolver=\'"<<iterativeSolver_->description()<<
"\'";
392 oss <<
",fwdOp=\'"<<lp_->getOperator()->description()<<
"\'";
393 if (lp_->getLeftPrec().get())
394 oss <<
",leftPrecOp=\'"<<lp_->getLeftPrec()->description()<<
"\'";
395 if (lp_->getRightPrec().get())
396 oss <<
",rightPrecOp=\'"<<lp_->getRightPrec()->description()<<
"\'";
405 template<
class Scalar>
413 using Teuchos::describe;
421 *out << this->description() << std::endl;
428 <<
"rangeDim=" << this->range()->dim()
429 <<
",domainDim=" << this->domain()->dim() <<
"}\n";
430 if (lp_->getOperator().get()) {
433 <<
"iterativeSolver = "<<describe(*iterativeSolver_,verbLevel)
434 <<
"fwdOp = " << describe(*lp_->getOperator(),verbLevel);
435 if (lp_->getLeftPrec().get())
436 *out <<
"leftPrecOp = "<<describe(*lp_->getLeftPrec(),verbLevel);
437 if (lp_->getRightPrec().get())
438 *out <<
"rightPrecOp = "<<describe(*lp_->getRightPrec(),verbLevel);
454 template<
class Scalar>
457 return ::Thyra::opSupported(*lp_->getOperator(),M_trans);
461 template<
class Scalar>
463 const EOpTransp M_trans,
464 const MultiVectorBase<Scalar> &X,
465 const Ptr<MultiVectorBase<Scalar> > &Y,
470 ::Thyra::apply<Scalar>(*lp_->getOperator(), M_trans, X, Y, alpha, beta);
477 template<
class Scalar>
481 return solveSupportsNewImpl(M_trans, Teuchos::null);
485 template<
class Scalar>
488 const Ptr<
const SolveCriteria<Scalar> > solveCriteria)
const
491 if (real_trans(transp)==
NOTRANS)
return true;
509 template<
class Scalar>
512 EOpTransp M_trans,
const SolveMeasureType& solveMeasureType)
const
514 SolveCriteria<Scalar> solveCriteria(solveMeasureType, SolveCriteria<Scalar>::unspecifiedTolerance());
515 return solveSupportsNewImpl(M_trans, Teuchos::constOptInArg(solveCriteria));
519 template<
class Scalar>
522 const EOpTransp M_trans,
523 const MultiVectorBase<Scalar> &
B,
524 const Ptr<MultiVectorBase<Scalar> > &X,
525 const Ptr<
const SolveCriteria<Scalar> > solveCriteria
529 THYRA_FUNC_TIME_MONITOR(
"Stratimikos: BelosLOWS");
532 using Teuchos::rcpFromRef;
533 using Teuchos::rcpFromPtr;
537 using Teuchos::parameterList;
538 using Teuchos::describe;
540 typedef typename ST::magnitudeType ScalarMag;
542 totalTimer.start(
true);
544 assertSolveSupports(*
this, M_trans, solveCriteria);
550 OSTab tab = this->getOSTab();
552 *out <<
"\nStarting iterations with Belos:\n";
554 *out <<
"Using forward operator = " << describe(*fwdOpSrc_->getOp(),verbLevel);
555 *out <<
"Using iterative solver = " << describe(*iterativeSolver_,verbLevel);
556 *out <<
"With #Eqns="<<
B.range()->dim()<<
", #RHSs="<<
B.domain()->dim()<<
" ...\n";
563 bool ret = lp_->setProblem( rcpFromPtr(X), rcpFromRef(
B) );
565 ret ==
false, CatastrophicSolveFailure
566 ,
"Error, the Belos::LinearProblem could not be set for the current solve!"
579 SolveMeasureType solveMeasureType;
582 solveMeasureType = solveCriteria->solveMeasureType;
583 const ScalarMag requestedTol = solveCriteria->requestedTol;
584 if (solveMeasureType.useDefault()) {
585 tmpPL->
set(
"Convergence Tolerance", defaultTol_);
587 else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_RHS)) {
588 if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
589 tmpPL->
set(
"Convergence Tolerance", requestedTol);
592 tmpPL->
set(
"Convergence Tolerance", defaultTol_);
594 setResidualScalingType (tmpPL, validPL,
"Norm of RHS");
596 else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_INIT_RESIDUAL)) {
597 if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
598 tmpPL->
set(
"Convergence Tolerance", requestedTol);
601 tmpPL->
set(
"Convergence Tolerance", defaultTol_);
603 setResidualScalingType (tmpPL, validPL,
"Norm of Initial Residual");
607 generalSolveCriteriaBelosStatusTest = createGeneralSolveCriteriaBelosStatusTest(
608 *solveCriteria, convergenceTestFrequency_);
610 generalSolveCriteriaBelosStatusTest->setOStream(out);
611 generalSolveCriteriaBelosStatusTest->setVerbLevel(
incrVerbLevel(verbLevel, -1));
614 tmpPL->
set(
"Convergence Tolerance", 1.0);
617 if (
nonnull(solveCriteria->extraParameters)) {
618 if (Teuchos::isParameterType<int>(*solveCriteria->extraParameters,
"Maximum Iterations")) {
619 tmpPL->
set(
"Maximum Iterations", Teuchos::get<int>(*solveCriteria->extraParameters,
"Maximum Iterations"));
625 validPL->
isParameter (
"Implicit Residual Scaling"))
626 tmpPL->
set(
"Implicit Residual Scaling",
627 "Norm of Preconditioned Initial Residual");
631 tmpPL->
set(
"Convergence Tolerance", defaultTol_);
638 Belos::ReturnType belosSolveStatus;
648 tmpPL->
set(
"Output Stream", outUsed);
650 if (
nonnull(generalSolveCriteriaBelosStatusTest)) {
651 iterativeSolver_->setUserConvStatusTest(generalSolveCriteriaBelosStatusTest);
654 belosSolveStatus = iterativeSolver_->solve();
656 catch (Belos::BelosError&) {
657 belosSolveStatus = Belos::Unconverged;
667 SolveStatus<Scalar> solveStatus;
669 switch (belosSolveStatus) {
670 case Belos::Unconverged: {
671 solveStatus.solveStatus = SOLVE_STATUS_UNCONVERGED;
682 solveStatus.achievedTol = iterativeSolver_->achievedTol();
683 }
catch (std::runtime_error&) {
688 case Belos::Converged: {
689 solveStatus.solveStatus = SOLVE_STATUS_CONVERGED;
690 if (
nonnull(generalSolveCriteriaBelosStatusTest)) {
695 const ArrayView<const ScalarMag> achievedTol =
696 generalSolveCriteriaBelosStatusTest->achievedTol();
697 solveStatus.achievedTol = ST::zero();
698 for (
Ordinal i = 0; i < achievedTol.size(); ++i) {
699 solveStatus.achievedTol = std::max(solveStatus.achievedTol, achievedTol[i]);
706 solveStatus.achievedTol = iterativeSolver_->achievedTol();
707 }
catch (std::runtime_error&) {
710 solveStatus.achievedTol = tmpPL->
get(
"Convergence Tolerance", defaultTol_);
718 std::ostringstream ossmessage;
720 <<
"The Belos solver of type \""<<iterativeSolver_->description()
721 <<
"\" returned a solve status of \""<<
toString(solveStatus.solveStatus) <<
"\""
722 <<
" in " << iterativeSolver_->getNumIters() <<
" iterations"
723 <<
" with total CPU time of " << totalTimer.totalElapsedTime() <<
" sec" ;
725 *out <<
"\n" << ossmessage.str() <<
"\n";
727 solveStatus.message = ossmessage.str();
732 if (solveStatus.extraParameters.is_null()) {
733 solveStatus.extraParameters = parameterList ();
735 solveStatus.extraParameters->set (
"Belos/Iteration Count",
736 iterativeSolver_->getNumIters());\
738 solveStatus.extraParameters->set (
"Iteration Count",
739 iterativeSolver_->getNumIters());\
746 solveStatus.extraParameters->set (
"Belos/Achieved Tolerance",
747 solveStatus.achievedTol);
762 #endif // THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP