136 const RCP<Teuchos::ParameterList> &solverPL,
138 const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
139 const RCP<
const PreconditionerBase<Scalar> > &prec,
140 const bool isExternalPrec_in,
141 const RCP<
const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
142 const ESupportSolveUse &supportSolveUse_in,
143 const int convergenceTestFrequency
147 using Teuchos::TypeNameTraits;
148 using Teuchos::Exceptions::InvalidParameterType;
149 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitude_type;
151 this->setLinePrefix(
"BELOS/T");
154 solverPL_ = solverPL;
155 iterativeSolver_ = iterativeSolver;
156 fwdOpSrc_ = fwdOpSrc;
158 isExternalPrec_ = isExternalPrec_in;
159 approxFwdOpSrc_ = approxFwdOpSrc;
160 supportSolveUse_ = supportSolveUse_in;
161 convergenceTestFrequency_ = convergenceTestFrequency;
165 if ( !is_null(solverPL_) ) {
166 if (solverPL_->isParameter(
"Convergence Tolerance")) {
172 if (solverPL_->isType<
double> (
"Convergence Tolerance")) {
174 as<magnitude_type> (solverPL_->get<
double> (
"Convergence Tolerance"));
176 else if (std::is_same_v<double, magnitude_type>) {
179 TEUCHOS_TEST_FOR_EXCEPTION(
180 true, std::invalid_argument,
"BelosLinearOpWithSolve::initialize: "
181 "The \"Convergence Tolerance\" parameter, which you provided, must "
182 "have type double (the type of the magnitude of Scalar = double).");
184 else if (solverPL_->isType<magnitude_type> (
"Convergence Tolerance")) {
185 defaultTol_ = solverPL_->get<magnitude_type> (
"Convergence Tolerance");
193 TEUCHOS_TEST_FOR_EXCEPTION(
194 true, InvalidParameterType,
"BelosLinearOpWithSolve::initialize: "
195 "The \"Convergence Tolerance\" parameter, which you provided, must "
196 "have type double (preferred) or the type of the magnitude of Scalar "
197 "= " << TypeNameTraits<Scalar>::name () <<
", which is " <<
198 TypeNameTraits<magnitude_type>::name () <<
" in this case. You can "
199 "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
203 if (solverPL_->isParameter(
"Timer Label") && solverPL_->isType<std::string>(
"Timer Label")) {
204 label_ = solverPL_->get<std::string>(
"Timer Label");
205 lp_->setLabel(label_);
207 if (solverPL_->isParameter(
"Filename LHS") && solverPL_->isType<std::string>(
"Filename LHS")) {
208 filenameLHS_ = solverPL_->get<std::string>(
"Filename LHS");
210 if (solverPL_->isParameter(
"Filename RHS") && solverPL_->isType<std::string>(
"Filename RHS")) {
211 filenameRHS_ = solverPL_->get<std::string>(
"Filename RHS");
215 RCP<const Teuchos::ParameterList> defaultPL =
216 iterativeSolver->getValidParameters();
222 if (defaultPL->isType<
double> (
"Convergence Tolerance")) {
224 as<magnitude_type> (defaultPL->get<
double> (
"Convergence Tolerance"));
226 else if (std::is_same_v<double, magnitude_type>) {
229 TEUCHOS_TEST_FOR_EXCEPTION(
230 true, std::invalid_argument,
"BelosLinearOpWithSolve::initialize: "
231 "The \"Convergence Tolerance\" parameter, which you provided, must "
232 "have type double (the type of the magnitude of Scalar = double).");
234 else if (defaultPL->isType<magnitude_type> (
"Convergence Tolerance")) {
235 defaultTol_ = defaultPL->get<magnitude_type> (
"Convergence Tolerance");
243 TEUCHOS_TEST_FOR_EXCEPTION(
244 true, InvalidParameterType,
"BelosLinearOpWithSolve::initialize: "
245 "The \"Convergence Tolerance\" parameter, which you provided, must "
246 "have type double (preferred) or the type of the magnitude of Scalar "
247 "= " << TypeNameTraits<Scalar>::name () <<
", which is " <<
248 TypeNameTraits<magnitude_type>::name () <<
" in this case. You can "
249 "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
305 RCP<Teuchos::ParameterList> *solverPL,
307 RCP<
const LinearOpSourceBase<Scalar> > *fwdOpSrc,
308 RCP<
const PreconditionerBase<Scalar> > *prec,
309 bool *isExternalPrec_in,
310 RCP<
const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
311 ESupportSolveUse *supportSolveUse_in
315 if (solverPL) *solverPL = solverPL_;
316 if (iterativeSolver) *iterativeSolver = iterativeSolver_;
317 if (fwdOpSrc) *fwdOpSrc = fwdOpSrc_;
318 if (prec) *prec = prec_;
319 if (isExternalPrec_in) *isExternalPrec_in = isExternalPrec_;
320 if (approxFwdOpSrc) *approxFwdOpSrc = approxFwdOpSrc_;
321 if (supportSolveUse_in) *supportSolveUse_in = supportSolveUse_;
324 solverPL_ = Teuchos::null;
325 iterativeSolver_ = Teuchos::null;
326 fwdOpSrc_ = Teuchos::null;
327 prec_ = Teuchos::null;
328 isExternalPrec_ =
false;
329 approxFwdOpSrc_ = Teuchos::null;
330 supportSolveUse_ = SUPPORT_SOLVE_UNSPECIFIED;
506 const EOpTransp M_trans,
507 const MultiVectorBase<Scalar> &B,
508 const Ptr<MultiVectorBase<Scalar> > &X,
509 const Ptr<
const SolveCriteria<Scalar> > solveCriteria
513 THYRA_FUNC_TIME_MONITOR(
"Stratimikos: BelosLOWS");
516 using Teuchos::rcpFromRef;
517 using Teuchos::rcpFromPtr;
518 using Teuchos::FancyOStream;
519 using Teuchos::OSTab;
520 using Teuchos::ParameterList;
521 using Teuchos::parameterList;
522 using Teuchos::describe;
523 typedef Teuchos::ScalarTraits<Scalar> ST;
524 typedef typename ST::magnitudeType ScalarMag;
525 Teuchos::Time totalTimer(
"Stratimikos::BelosLinearOpWithSolve::totalTime");
526 totalTimer.start(
true);
528 assertSolveSupports(*
this, M_trans, solveCriteria);
532 const RCP<FancyOStream> out = this->getOStream();
533 const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
534 OSTab tab = this->getOSTab();
535 if (out.get() &&
static_cast<int>(verbLevel) >
static_cast<int>(Teuchos::VERB_LOW)) {
536 *out <<
"\nStarting iterations with Belos:\n";
538 *out <<
"Using forward operator = " << describe(*fwdOpSrc_->getOp(),verbLevel);
539 *out <<
"Using iterative solver = " << describe(*iterativeSolver_,verbLevel);
540 *out <<
"With #Eqns="<<B.range()->dim()<<
", #RHSs="<<B.domain()->dim()<<
" ...\n";
543#ifdef HAVE_STRATIMIKOS_THYRATPETRAADAPTERS
547 if (filenameLHS_ !=
"") {
549 auto tmv = Thyra::TpetraOperatorVectorExtraction<Scalar>::getTpetraMultiVector(Teuchos::rcpFromPtr(X));
550 Tpetra::MatrixMarket::Writer<::Tpetra::CrsMatrix<Scalar> >::writeDenseFile(filenameLHS_+
"." + label_ +
"." + std::to_string(counter_), *tmv,
"",
"");
551 }
catch (
const std::logic_error&) {
552 *out <<
"Warning: Cannot write LHS multivector to file.\n";
555 if (filenameRHS_ !=
"") {
557 auto tmv = Thyra::TpetraOperatorVectorExtraction<Scalar>::getConstTpetraMultiVector(Teuchos::rcpFromRef(B));
558 Tpetra::MatrixMarket::Writer<::Tpetra::CrsMatrix<Scalar> >::writeDenseFile(filenameRHS_+
"." + label_ +
"." + std::to_string(counter_), *tmv,
"",
"");
559 }
catch (
const std::logic_error&) {
560 *out <<
"Warning: Cannot write RHS multivector to file.\n";
570 bool ret = lp_->setProblem( rcpFromPtr(X), rcpFromRef(B) );
571 TEUCHOS_TEST_FOR_EXCEPTION(
572 ret ==
false, CatastrophicSolveFailure
573 ,
"Error, the Belos::LinearProblem could not be set for the current solve!"
581 const RCP<ParameterList> tmpPL = Teuchos::parameterList();
584 RCP<const ParameterList> validPL = iterativeSolver_->getValidParameters();
586 SolveMeasureType solveMeasureType;
587 RCP<GeneralSolveCriteriaBelosStatusTest<Scalar> > generalSolveCriteriaBelosStatusTest;
588 if (nonnull(solveCriteria)) {
589 solveMeasureType = solveCriteria->solveMeasureType;
590 const ScalarMag requestedTol = solveCriteria->requestedTol;
591 if (solveMeasureType.useDefault()) {
592 tmpPL->set(
"Convergence Tolerance", defaultTol_);
594 else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_RHS)) {
595 if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
596 tmpPL->set(
"Convergence Tolerance", requestedTol);
599 tmpPL->set(
"Convergence Tolerance", defaultTol_);
601 setResidualScalingType (tmpPL, validPL,
"Norm of RHS");
603 else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_INIT_RESIDUAL)) {
604 if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
605 tmpPL->set(
"Convergence Tolerance", requestedTol);
608 tmpPL->set(
"Convergence Tolerance", defaultTol_);
610 setResidualScalingType (tmpPL, validPL,
"Norm of Initial Residual");
614 generalSolveCriteriaBelosStatusTest = createGeneralSolveCriteriaBelosStatusTest(
615 *solveCriteria, convergenceTestFrequency_);
617 generalSolveCriteriaBelosStatusTest->setOStream(out);
618 generalSolveCriteriaBelosStatusTest->setVerbLevel(incrVerbLevel(verbLevel, -1));
621 tmpPL->set(
"Convergence Tolerance", 1.0);
624 if (nonnull(solveCriteria->extraParameters)) {
625 if (Teuchos::isParameterType<int>(*solveCriteria->extraParameters,
"Maximum Iterations")) {
626 tmpPL->set(
"Maximum Iterations", Teuchos::get<int>(*solveCriteria->extraParameters,
"Maximum Iterations"));
631 if (Teuchos::nonnull(lp_->getLeftPrec()) &&
632 validPL->isParameter (
"Implicit Residual Scaling"))
633 tmpPL->set(
"Implicit Residual Scaling",
634 "Norm of Preconditioned Initial Residual");
638 tmpPL->set(
"Convergence Tolerance", defaultTol_);
650 (
static_cast<int>(verbLevel) >=
static_cast<int>(Teuchos::VERB_LOW)
652 : rcp(
new FancyOStream(rcp(
new Teuchos::oblackholestream())))
654 Teuchos::OSTab tab1(outUsed,1,
"BELOS");
655 tmpPL->set(
"Output Stream", outUsed);
657 iterativeSolver_->setParameters(tmpPL);
658 init_ = !nonnull(solveCriteria);
660 if (nonnull(generalSolveCriteriaBelosStatusTest)) {
661 iterativeSolver_->setUserConvStatusTest(generalSolveCriteriaBelosStatusTest);
664 belosSolveStatus = iterativeSolver_->solve();
667 TEUCHOS_TEST_FOR_EXCEPTION(
true,
668 CatastrophicSolveFailure,
679 SolveStatus<Scalar> solveStatus;
681 switch (belosSolveStatus) {
683 solveStatus.solveStatus = SOLVE_STATUS_UNCONVERGED;
694 solveStatus.achievedTol = iterativeSolver_->achievedTol();
695 }
catch (std::runtime_error&) {
701 solveStatus.solveStatus = SOLVE_STATUS_CONVERGED;
702 if (nonnull(generalSolveCriteriaBelosStatusTest)) {
707 const ArrayView<const ScalarMag> achievedTol =
708 generalSolveCriteriaBelosStatusTest->achievedTol();
709 solveStatus.achievedTol = Teuchos::ScalarTraits<ScalarMag>::zero();
710 for (Ordinal i = 0; i < achievedTol.size(); ++i) {
711 solveStatus.achievedTol = std::max(solveStatus.achievedTol, achievedTol[i]);
718 solveStatus.achievedTol = iterativeSolver_->achievedTol();
719 }
catch (std::runtime_error&) {
722 solveStatus.achievedTol = tmpPL->get(
"Convergence Tolerance", defaultTol_);
727 TEUCHOS_SWITCH_DEFAULT_DEBUG_ASSERT();
730 std::ostringstream ossmessage;
732 <<
"The Belos solver " << (label_ !=
"" ? (
"\"" + label_ +
"\" ") :
"")
733 <<
"of type \""<<iterativeSolver_->description()
734 <<
"\" returned a solve status of \""<< toString(solveStatus.solveStatus) <<
"\""
735 <<
" in " << iterativeSolver_->getNumIters() <<
" iterations"
736 <<
" with total CPU time of " << totalTimer.totalElapsedTime() <<
" sec" ;
737 if (out.get() &&
static_cast<int>(verbLevel) >
static_cast<int>(Teuchos::VERB_LOW))
738 *out <<
"\n" << ossmessage.str() <<
"\n";
740 solveStatus.message = ossmessage.str();
745 if (solveStatus.extraParameters.is_null()) {
746 solveStatus.extraParameters = parameterList ();
748 solveStatus.extraParameters->set (
"Belos/Iteration Count",
749 iterativeSolver_->getNumIters());\
751 solveStatus.extraParameters->set (
"Iteration Count",
752 iterativeSolver_->getNumIters());\
759 solveStatus.extraParameters->set (
"Belos/Achieved Tolerance",
760 solveStatus.achievedTol);