10 #ifndef THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP 
   11 #define THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP 
   14 #include "Thyra_GeneralSolveCriteriaBelosStatusTest.hpp" 
   15 #include "Thyra_LinearOpWithSolveHelpers.hpp" 
   19 #include "Stratimikos_Config.h" 
   20 #ifdef HAVE_STRATIMIKOS_THYRATPETRAADAPTERS 
   22 #  include <MatrixMarket_Tpetra.hpp> 
   23 #  include <TpetraExt_MatrixMatrix.hpp> 
   73                           const std::string& residualScalingType)
 
  102     if (solverValidParams->
isParameter (
"Implicit Residual Scaling")) {
 
  103       solverParams->
set (
"Implicit Residual Scaling", residualScalingType);
 
  105     if (solverValidParams->
isParameter (
"Explicit Residual Scaling")) {
 
  106       solverParams->
set (
"Explicit Residual Scaling", residualScalingType);
 
  119 template<
class Scalar>
 
  121   :convergenceTestFrequency_(-1),
 
  122   isExternalPrec_(false),
 
  123   supportSolveUse_(SUPPORT_SOLVE_UNSPECIFIED),
 
  132 template<
class Scalar>
 
  134   const RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > &lp,
 
  136   const RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > &iterativeSolver,
 
  137   const RCP<
const LinearOpSourceBase<Scalar> > &fwdOpSrc,
 
  138   const RCP<
const PreconditionerBase<Scalar> > &prec,
 
  139   const bool isExternalPrec_in,
 
  140   const RCP<
const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
 
  141   const ESupportSolveUse &supportSolveUse_in,
 
  142   const int convergenceTestFrequency
 
  150   this->setLinePrefix(
"BELOS/T");
 
  153   solverPL_ = solverPL;
 
  154   iterativeSolver_ = iterativeSolver;
 
  155   fwdOpSrc_ = fwdOpSrc;
 
  157   isExternalPrec_ = isExternalPrec_in;
 
  158   approxFwdOpSrc_ = approxFwdOpSrc;
 
  159   supportSolveUse_ = supportSolveUse_in;
 
  160   convergenceTestFrequency_ = convergenceTestFrequency;
 
  164     if (solverPL_->isParameter(
"Convergence Tolerance")) {
 
  170       if (solverPL_->isType<
double> (
"Convergence Tolerance")) {
 
  172           as<magnitude_type> (solverPL_->get<
double> (
"Convergence Tolerance"));
 
  174       else if (std::is_same_v<double, magnitude_type>) {
 
  178           true, std::invalid_argument, 
"BelosLinearOpWithSolve::initialize: " 
  179           "The \"Convergence Tolerance\" parameter, which you provided, must " 
  180           "have type double (the type of the magnitude of Scalar = double).");
 
  182       else if (solverPL_->isType<magnitude_type> (
"Convergence Tolerance")) {
 
  183         defaultTol_ = solverPL_->get<magnitude_type> (
"Convergence Tolerance");
 
  192           true, InvalidParameterType, 
"BelosLinearOpWithSolve::initialize: " 
  193           "The \"Convergence Tolerance\" parameter, which you provided, must " 
  194           "have type double (preferred) or the type of the magnitude of Scalar " 
  195           "= " << TypeNameTraits<Scalar>::name () << 
", which is " <<
 
  196           TypeNameTraits<magnitude_type>::name () << 
" in this case.  You can " 
  197           "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
 
  201     if (solverPL_->isParameter(
"Timer Label") && solverPL_->isType<std::string>(
"Timer Label")) {
 
  202       label_ = solverPL_->get<std::string>(
"Timer Label");
 
  203       lp_->setLabel(label_);
 
  205     if (solverPL_->isParameter(
"Filename LHS") && solverPL_->isType<std::string>(
"Filename LHS")) {
 
  206       filenameLHS_ = solverPL_->get<std::string>(
"Filename LHS");
 
  208     if (solverPL_->isParameter(
"Filename RHS") && solverPL_->isType<std::string>(
"Filename RHS")) {
 
  209       filenameRHS_ = solverPL_->get<std::string>(
"Filename RHS");
 
  214       iterativeSolver->getValidParameters();
 
  220     if (defaultPL->
isType<
double> (
"Convergence Tolerance")) {
 
  222         as<magnitude_type> (defaultPL->
get<
double> (
"Convergence Tolerance"));
 
  224     else if (std::is_same_v<double, magnitude_type>) {
 
  228         true, std::invalid_argument, 
"BelosLinearOpWithSolve::initialize: " 
  229         "The \"Convergence Tolerance\" parameter, which you provided, must " 
  230         "have type double (the type of the magnitude of Scalar = double).");
 
  232     else if (defaultPL->
isType<magnitude_type> (
"Convergence Tolerance")) {
 
  233       defaultTol_ = defaultPL->
get<magnitude_type> (
"Convergence Tolerance");
 
  242         true, InvalidParameterType, 
"BelosLinearOpWithSolve::initialize: " 
  243         "The \"Convergence Tolerance\" parameter, which you provided, must " 
  244         "have type double (preferred) or the type of the magnitude of Scalar " 
  245         "= " << TypeNameTraits<Scalar>::name () << 
", which is " <<
 
  246         TypeNameTraits<magnitude_type>::name () << 
" in this case.  You can " 
  247         "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
 
  253 template<
class Scalar>
 
  258     _fwdOpSrc = fwdOpSrc_;
 
  264 template<
class Scalar>
 
  275 template<
class Scalar>
 
  278   return isExternalPrec_;
 
  282 template<
class Scalar>
 
  287     _approxFwdOpSrc = approxFwdOpSrc_;
 
  289   return _approxFwdOpSrc;
 
  293 template<
class Scalar>
 
  296   return supportSolveUse_;
 
  300 template<
class Scalar>
 
  302   RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > *lp,
 
  304   RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > *iterativeSolver,
 
  305   RCP<
const LinearOpSourceBase<Scalar> > *fwdOpSrc,
 
  306   RCP<
const PreconditionerBase<Scalar> > *prec,
 
  307   bool *isExternalPrec_in,
 
  308   RCP<
const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
 
  309   ESupportSolveUse *supportSolveUse_in
 
  313   if (solverPL) *solverPL = solverPL_;
 
  314   if (iterativeSolver) *iterativeSolver = iterativeSolver_;
 
  315   if (fwdOpSrc) *fwdOpSrc = fwdOpSrc_;
 
  316   if (prec) *prec = prec_;
 
  317   if (isExternalPrec_in) *isExternalPrec_in = isExternalPrec_;
 
  318   if (approxFwdOpSrc) *approxFwdOpSrc = approxFwdOpSrc_;
 
  319   if (supportSolveUse_in) *supportSolveUse_in = supportSolveUse_;
 
  326   isExternalPrec_ = 
false;
 
  328   supportSolveUse_ = SUPPORT_SOLVE_UNSPECIFIED;
 
  335 template<
class Scalar>
 
  340     return lp_->getOperator()->range();
 
  345 template<
class Scalar>
 
  350     return lp_->getOperator()->domain();
 
  355 template<
class Scalar>
 
  366 template<
class Scalar>
 
  369   std::ostringstream oss;
 
  373     oss << 
"iterativeSolver=\'"<<iterativeSolver_->description()<<
"\'";
 
  374     oss << 
",fwdOp=\'"<<lp_->getOperator()->description()<<
"\'";
 
  375     if (lp_->getLeftPrec().get())
 
  376       oss << 
",leftPrecOp=\'"<<lp_->getLeftPrec()->description()<<
"\'";
 
  377     if (lp_->getRightPrec().get())
 
  378       oss << 
",rightPrecOp=\'"<<lp_->getRightPrec()->description()<<
"\'";
 
  387 template<
class Scalar>
 
  395   using Teuchos::describe;
 
  402       *out << this->description() << std::endl;
 
  409         << 
"rangeDim=" << this->range()->dim()
 
  410         << 
",domainDim=" << this->domain()->dim() << 
"}\n";
 
  411       if (lp_->getOperator().get()) {
 
  414           << 
"iterativeSolver = "<<describe(*iterativeSolver_,verbLevel)
 
  415           << 
"fwdOp = " << describe(*lp_->getOperator(),verbLevel);
 
  416         if (lp_->getLeftPrec().get())
 
  417           *out << 
"leftPrecOp = "<<describe(*lp_->getLeftPrec(),verbLevel);
 
  418         if (lp_->getRightPrec().get())
 
  419           *out << 
"rightPrecOp = "<<describe(*lp_->getRightPrec(),verbLevel);
 
  435 template<
class Scalar>
 
  438   return ::Thyra::opSupported(*lp_->getOperator(),M_trans);
 
  442 template<
class Scalar>
 
  444   const EOpTransp M_trans,
 
  445   const MultiVectorBase<Scalar> &X,
 
  446   const Ptr<MultiVectorBase<Scalar> > &Y,
 
  451   ::Thyra::apply<Scalar>(*lp_->getOperator(), M_trans, X, Y, alpha, beta);
 
  458 template<
class Scalar>
 
  466 template<
class Scalar>
 
  469   const Ptr<
const SolveCriteria<Scalar> > )
 const 
  472   if (real_trans(transp)==
NOTRANS) 
return true;
 
  490 template<
class Scalar>
 
  493   EOpTransp M_trans, 
const SolveMeasureType& solveMeasureType)
 const 
  495   SolveCriteria<Scalar> solveCriteria(solveMeasureType, SolveCriteria<Scalar>::unspecifiedTolerance());
 
  496   return solveSupportsNewImpl(M_trans, Teuchos::constOptInArg(solveCriteria));
 
  500 template<
class Scalar>
 
  503   const EOpTransp M_trans,
 
  504   const MultiVectorBase<Scalar> &
B,
 
  505   const Ptr<MultiVectorBase<Scalar> > &X,
 
  506   const Ptr<
const SolveCriteria<Scalar> > solveCriteria
 
  510   THYRA_FUNC_TIME_MONITOR(
"Stratimikos: BelosLOWS");
 
  513   using Teuchos::rcpFromRef;
 
  514   using Teuchos::rcpFromPtr;
 
  518   using Teuchos::parameterList;
 
  519   using Teuchos::describe;
 
  521   typedef typename ST::magnitudeType ScalarMag;
 
  522   Teuchos::Time totalTimer(
"Stratimikos::BelosLinearOpWithSolve::totalTime");
 
  523   totalTimer.
start(
true);
 
  525   assertSolveSupports(*
this, M_trans, solveCriteria);
 
  531   OSTab tab = this->getOSTab();
 
  533     *out << 
"\nStarting iterations with Belos:\n";
 
  535     *out << 
"Using forward operator = " << describe(*fwdOpSrc_->getOp(),verbLevel);
 
  536     *out << 
"Using iterative solver = " << describe(*iterativeSolver_,verbLevel);
 
  537     *out << 
"With #Eqns="<<B.range()->dim()<<
", #RHSs="<<B.domain()->dim()<<
" ...\n";
 
  540 #ifdef HAVE_STRATIMIKOS_THYRATPETRAADAPTERS 
  544   if (filenameLHS_ != 
"") {
 
  547       Tpetra::MatrixMarket::Writer<::Tpetra::CrsMatrix<Scalar> >::writeDenseFile(filenameLHS_+ 
"." + label_ + 
"." + std::to_string(counter_), *tmv, 
"", 
"");
 
  548     } 
catch (
const std::logic_error&) {
 
  549       *out << 
"Warning: Cannot write LHS multivector to file.\n";
 
  552   if (filenameRHS_ != 
"") {
 
  555       Tpetra::MatrixMarket::Writer<::Tpetra::CrsMatrix<Scalar> >::writeDenseFile(filenameRHS_+ 
"." + label_ + 
"." + std::to_string(counter_), *tmv, 
"", 
"");
 
  556     } 
catch (
const std::logic_error&) {
 
  557       *out << 
"Warning: Cannot write RHS multivector to file.\n";
 
  567   bool ret = lp_->setProblem( rcpFromPtr(X), rcpFromRef(B) );
 
  569     ret == 
false, CatastrophicSolveFailure
 
  570     ,
"Error, the Belos::LinearProblem could not be set for the current solve!" 
  583   SolveMeasureType solveMeasureType;
 
  586     solveMeasureType = solveCriteria->solveMeasureType;
 
  587     const ScalarMag requestedTol = solveCriteria->requestedTol;
 
  588     if (solveMeasureType.useDefault()) {
 
  589       tmpPL->
set(
"Convergence Tolerance", defaultTol_);
 
  591     else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_RHS)) {
 
  592       if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
 
  593         tmpPL->
set(
"Convergence Tolerance", requestedTol);
 
  596         tmpPL->
set(
"Convergence Tolerance", defaultTol_);
 
  598       setResidualScalingType (tmpPL, validPL, 
"Norm of RHS");
 
  600     else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_INIT_RESIDUAL)) {
 
  601       if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
 
  602         tmpPL->
set(
"Convergence Tolerance", requestedTol);
 
  605         tmpPL->
set(
"Convergence Tolerance", defaultTol_);
 
  607       setResidualScalingType (tmpPL, validPL, 
"Norm of Initial Residual");
 
  611       generalSolveCriteriaBelosStatusTest = createGeneralSolveCriteriaBelosStatusTest(
 
  612         *solveCriteria, convergenceTestFrequency_);
 
  614       generalSolveCriteriaBelosStatusTest->setOStream(out);
 
  615       generalSolveCriteriaBelosStatusTest->setVerbLevel(
incrVerbLevel(verbLevel, -1));
 
  618       tmpPL->
set(
"Convergence Tolerance", 1.0);
 
  621     if (
nonnull(solveCriteria->extraParameters)) {
 
  622       if (Teuchos::isParameterType<int>(*solveCriteria->extraParameters,
"Maximum Iterations")) {
 
  623         tmpPL->
set(
"Maximum Iterations", Teuchos::get<int>(*solveCriteria->extraParameters,
"Maximum Iterations"));
 
  629         validPL->
isParameter (
"Implicit Residual Scaling"))
 
  630       tmpPL->
set(
"Implicit Residual Scaling",
 
  631                  "Norm of Preconditioned Initial Residual");
 
  635     tmpPL->
set(
"Convergence Tolerance", defaultTol_);
 
  642   Belos::ReturnType belosSolveStatus;
 
  652     tmpPL->
set(
"Output Stream", outUsed);
 
  654     if (
nonnull(generalSolveCriteriaBelosStatusTest)) {
 
  655       iterativeSolver_->setUserConvStatusTest(generalSolveCriteriaBelosStatusTest);
 
  658       belosSolveStatus = iterativeSolver_->solve();
 
  660     catch (Belos::BelosError& ex) {
 
  662                                   CatastrophicSolveFailure,
 
  673   SolveStatus<Scalar> solveStatus;
 
  675   switch (belosSolveStatus) {
 
  676     case Belos::Unconverged: {
 
  677       solveStatus.solveStatus = SOLVE_STATUS_UNCONVERGED;
 
  688         solveStatus.achievedTol = iterativeSolver_->achievedTol();
 
  689       } 
catch (std::runtime_error&) {
 
  694     case Belos::Converged: {
 
  695       solveStatus.solveStatus = SOLVE_STATUS_CONVERGED;
 
  696       if (
nonnull(generalSolveCriteriaBelosStatusTest)) {
 
  701         const ArrayView<const ScalarMag> achievedTol =
 
  702           generalSolveCriteriaBelosStatusTest->achievedTol();
 
  704         for (
Ordinal i = 0; i < achievedTol.size(); ++i) {
 
  705           solveStatus.achievedTol = std::max(solveStatus.achievedTol, achievedTol[i]);
 
  712           solveStatus.achievedTol = iterativeSolver_->achievedTol();
 
  713         } 
catch (std::runtime_error&) {
 
  716           solveStatus.achievedTol = tmpPL->
get(
"Convergence Tolerance", defaultTol_);
 
  724   std::ostringstream ossmessage;
 
  726     << 
"The Belos solver " << (label_ != 
"" ? (
"\"" + label_  + 
"\" ") : 
"")
 
  727     << 
"of type \""<<iterativeSolver_->description()
 
  728     <<
"\" returned a solve status of \""<< 
toString(solveStatus.solveStatus) << 
"\"" 
  729     << 
" in " << iterativeSolver_->getNumIters() << 
" iterations" 
  732     *out << 
"\n" << ossmessage.str() << 
"\n";
 
  734   solveStatus.message = ossmessage.str();
 
  739   if (solveStatus.extraParameters.is_null()) {
 
  740     solveStatus.extraParameters = parameterList ();
 
  742   solveStatus.extraParameters->set (
"Belos/Iteration Count",
 
  743                                     iterativeSolver_->getNumIters());\
 
  745   solveStatus.extraParameters->set (
"Iteration Count",
 
  746                                     iterativeSolver_->getNumIters());\
 
  753   solveStatus.extraParameters->set (
"Belos/Achieved Tolerance",
 
  754                                     solveStatus.achievedTol);
 
  769 #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. 
 
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 
 
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
 
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 
 
void start(bool reset=false)
 
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. 
 
double totalElapsedTime(bool readCurrentTime=false) const 
 
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)