11 #include "Thyra_LinearOpWithSolveHelpers.hpp" 
   47 class SetAztecSolveState {
 
   53     const Thyra::SolveMeasureType &solveMeasureType
 
   55   ~SetAztecSolveState();
 
   66 SetAztecSolveState::SetAztecSolveState(
 
   70   const Thyra::SolveMeasureType &solveMeasureType
 
   72   :aztecSolver_(aztecSolver.assert_not_null())
 
   77   verbLevel_ = verbLevel;
 
   83       fancyOStream_= Teuchos::tab(
 
   88       aztecSolver_->SetOutputStream(*fancyOStream_);
 
   89       aztecSolver_->SetErrorStream(*fancyOStream_);
 
   99     outputFrequency_ = aztecSolver_->GetAllAztecOptions()[AZ_output];
 
  100     aztecSolver_->SetAztecOption(AZ_output,0);
 
  104   convergenceTest_ = aztecSolver_->GetAztecOption(AZ_conv);
 
  105   if (solveMeasureType.useDefault())
 
  111       Thyra::SOLVE_MEASURE_NORM_RESIDUAL,
 
  112       Thyra::SOLVE_MEASURE_NORM_RHS
 
  116     aztecSolver_->SetAztecOption(AZ_conv,AZ_rhs);
 
  120       Thyra::SOLVE_MEASURE_NORM_RESIDUAL,
 
  121       Thyra::SOLVE_MEASURE_NORM_INIT_RESIDUAL
 
  125     aztecSolver_->SetAztecOption(AZ_conv,AZ_r0);
 
  134 SetAztecSolveState::~SetAztecSolveState()
 
  140       aztecSolver_->SetOutputStream(std::cout);
 
  141       aztecSolver_->SetErrorStream(std::cerr);
 
  142       *fancyOStream_ << 
"\n";
 
  146     aztecSolver_->SetAztecOption(AZ_output,outputFrequency_);
 
  150   aztecSolver_->SetAztecOption(AZ_conv,convergenceTest_);
 
  165   const int fwdDefaultMaxIterations_in
 
  166   ,
const double fwdDefaultTol_in
 
  167   ,
const int adjDefaultMaxIterations_in
 
  168   ,
const double adjDefaultTol_in
 
  169   ,
const bool outputEveryRhs_in
 
  171   :fwdDefaultMaxIterations_(fwdDefaultMaxIterations_in)
 
  172   ,fwdDefaultTol_(fwdDefaultTol_in)
 
  173   ,adjDefaultMaxIterations_(adjDefaultMaxIterations_in)
 
  174   ,adjDefaultTol_(adjDefaultTol_in)
 
  175   ,outputEveryRhs_(outputEveryRhs_in)
 
  176   ,isExternalPrec_(false)
 
  177   ,allowInexactFwdSolve_(false)
 
  178   ,allowInexactAdjSolve_(false)
 
  179   ,aztecSolverScalar_(0.0)
 
  184   const RCP<
const LinearOpBase<double> > &fwdOp
 
  185   ,
const RCP<
const LinearOpSourceBase<double> > &fwdOpSrc
 
  186   ,
const RCP<
const PreconditionerBase<double> > &prec
 
  187   ,
const bool isExternalPrec_in
 
  188   ,
const RCP<
const LinearOpSourceBase<double> > &approxFwdOpSrc
 
  190   ,
const bool allowInexactFwdSolve
 
  192   ,
const bool allowInexactAdjSolve
 
  193   ,
const double aztecSolverScalar
 
  211   const std::string fwdOpLabel = fwdOp_->getObjectLabel();
 
  212   if (fwdOpLabel.length())
 
  213     this->setObjectLabel( 
"lows("+fwdOpLabel+
")" );
 
  249   return _approxFwdOpSrc;
 
  254   RCP<
const LinearOpBase<double> > *fwdOp,
 
  255   RCP<
const LinearOpSourceBase<double> > *fwdOpSrc,
 
  256   RCP<
const PreconditionerBase<double> > *prec,
 
  257   bool *isExternalPrec_inout,
 
  258   RCP<
const LinearOpSourceBase<double> > *approxFwdOpSrc,
 
  260   bool *allowInexactFwdSolve,
 
  262   bool *allowInexactAdjSolve,
 
  263   double *aztecSolverScalar
 
  266   if (fwdOp) *fwdOp = 
fwdOp_;
 
  268   if (prec) *prec = 
prec_;
 
  319   std::ostringstream oss;
 
  323     oss << 
"fwdOp="<<
fwdOp_->description()<<
"";
 
  337   using Teuchos::describe;
 
  349         << 
"rangeDim=" << this->
range()->dim()
 
  350         << 
",domainDim="<< this->
domain()->dim() << 
"}\n";
 
  365             << 
"Aztec Fwd Mat = " 
  369             << 
"Aztec Fwd Prec Op = " 
  373             << 
"Aztec Fwd Prec Mat = " 
  383             << 
"Aztec Adj Mat = " 
  387             << 
"Aztec Adj Prec Op = " 
  391             << 
"Aztec Adj Prec Mat = " 
  410   return ::Thyra::opSupported(*
fwdOp_,M_trans);
 
  415   const EOpTransp M_trans,
 
  416   const MultiVectorBase<double> &X,
 
  417   const Ptr<MultiVectorBase<double> > &Y,
 
  422   Thyra::apply( *
fwdOp_, M_trans, X, Y, alpha, beta );
 
  432   if (real_trans(M_trans)==
NOTRANS) 
return true;
 
  439   EOpTransp M_trans, 
const SolveMeasureType& solveMeasureType
 
  442   if (real_trans(M_trans)==
NOTRANS) {
 
  443     if (solveMeasureType.useDefault())
 
  449         SOLVE_MEASURE_NORM_RESIDUAL,
 
  450         SOLVE_MEASURE_NORM_RHS
 
  460         SOLVE_MEASURE_NORM_RESIDUAL,
 
  461         SOLVE_MEASURE_NORM_INIT_RESIDUAL
 
  476     else if (solveMeasureType.useDefault())
 
  482         SOLVE_MEASURE_NORM_RESIDUAL,
 
  483         SOLVE_MEASURE_NORM_RHS
 
  493         SOLVE_MEASURE_NORM_RESIDUAL,
 
  494         SOLVE_MEASURE_NORM_INIT_RESIDUAL
 
  513   const EOpTransp M_trans,
 
  514   const MultiVectorBase<double> &
B,
 
  515   const Ptr<MultiVectorBase<double> > &X,
 
  516   const Ptr<
const SolveCriteria<double> > solveCriteria
 
  521   using Teuchos::rcpFromRef;
 
  522   using Teuchos::rcpFromPtr;
 
  524   typedef SolveCriteria<double> SC;
 
  525   typedef SolveStatus<double> SS;
 
  527   THYRA_FUNC_TIME_MONITOR(
"Stratimikos: AztecOOLOWS");
 
  529   totalTimer.start(
true);
 
  533   OSTab tab = this->getOSTab();
 
  535     *out << 
"\nSolving block system using AztecOO ...\n\n";
 
  541   SolveMeasureType solveMeasureType;
 
  543     solveMeasureType = solveCriteria->solveMeasureType;
 
  544     assertSupportsSolveMeasureType(*
this, M_trans, solveMeasureType);
 
  550   const EOpTransp aztecOpTransp = real_trans(M_trans);
 
  557   const Epetra_Operator
 
  558     *aztecOp = aztecSolver->GetUserOperator();
 
  564     &opRangeMap = aztecOp->OperatorRangeMap(),
 
  565     &opDomainMap = aztecOp->OperatorDomainMap();
 
  570   double tol = ( aztecOpTransp==
NOTRANS ? fwdDefaultTol() : adjDefaultTol() );
 
  571   int maxIterations = ( aztecOpTransp==
NOTRANS 
  572     ? fwdDefaultMaxIterations() : adjDefaultMaxIterations() );
 
  573   bool isDefaultSolveCriteria = 
true;
 
  575     if ( solveCriteria->requestedTol != SC::unspecifiedTolerance() ) {
 
  576       tol = solveCriteria->requestedTol;
 
  577       isDefaultSolveCriteria = 
false;
 
  579     if (
nonnull(solveCriteria->extraParameters)) {
 
  580       maxIterations = solveCriteria->extraParameters->
get(
"Maximum Iterations",maxIterations);
 
  594   if (opWrapper == 0) {
 
  603   int totalIterations = 0;
 
  604   SolveStatus<double> solveStatus;
 
  605   solveStatus.solveStatus = SOLVE_STATUS_CONVERGED;
 
  606   solveStatus.achievedTol = -1.0;
 
  612   const int m = B.domain()->dim();
 
  614   for( 
int j = 0; j < m; ++j ) {
 
  616     THYRA_FUNC_TIME_MONITOR_DIFF(
"Stratimikos: AztecOOLOWS:SingleSolve", SingleSolve);
 
  629     if (opWrapper == 0) {
 
  630       epetra_b_j = rcpFromRef(*const_cast<Epetra_Vector*>((*epetra_B)(j)));
 
  631       epetra_x_j = rcpFromRef(*(*epetra_X)(j));
 
  635         epetra_b_j = 
rcp(
new Epetra_Vector(opRangeMap));
 
  636         epetra_x_j = 
rcp(
new Epetra_Vector(opDomainMap));
 
  646     aztecSolver->SetRHS(&*epetra_b_j);
 
  647     aztecSolver->SetLHS(&*epetra_x_j);
 
  655         setAztecSolveState(aztecSolver,out,verbLevel,solveMeasureType);
 
  656       aztecSolver->Iterate( maxIterations, tol );
 
  673     if (opWrapper != 0) {
 
  681     const int iterations = aztecSolver->NumIters();
 
  682     const double achievedTol = aztecSolver->ScaledResidual();
 
  683     const double *AZ_status = aztecSolver->GetAztecStatus();
 
  684     std::ostringstream oss;
 
  685     bool converged = 
false;
 
  686     if (AZ_status[AZ_why]==AZ_normal) { oss << 
"Aztec returned AZ_normal."; converged = 
true; }
 
  687     else if (AZ_status[AZ_why]==AZ_param) oss << 
"Aztec returned AZ_param.";
 
  688     else if (AZ_status[AZ_why]==AZ_breakdown) oss << 
"Aztec returned AZ_breakdown.";
 
  689     else if (AZ_status[AZ_why]==AZ_loss) oss << 
"Aztec returned AZ_loss.";
 
  690     else if (AZ_status[AZ_why]==AZ_ill_cond) oss << 
"Aztec returned AZ_ill_cond.";
 
  691     else if (AZ_status[AZ_why]==AZ_maxits) oss << 
"Aztec returned AZ_maxits.";
 
  692     else oss << 
"Aztec returned an unknown status?";
 
  693     oss << 
"  Iterations = " << iterations << 
".";
 
  694     oss << 
"  Achieved Tolerance = " << achievedTol << 
".";
 
  696     if (out.
get() && 
static_cast<int>(verbLevel) > static_cast<int>(
Teuchos::VERB_NONE) && outputEveryRhs())
 
  699     solveStatus.achievedTol = 
TEUCHOS_MAX(solveStatus.achievedTol, achievedTol);
 
  702     totalIterations += iterations;
 
  704     solveStatus.message = oss.str();
 
  705     if ( isDefaultSolveCriteria ) {
 
  706       switch(solveStatus.solveStatus) {
 
  707         case SOLVE_STATUS_UNKNOWN:
 
  710         case SOLVE_STATUS_CONVERGED:
 
  711           solveStatus.solveStatus = ( converged ? SOLVE_STATUS_CONVERGED : SOLVE_STATUS_UNCONVERGED );
 
  713         case SOLVE_STATUS_UNCONVERGED:
 
  722   aztecSolver->UnsetLHSRHS();
 
  734   SolveStatus<double> overallSolveStatus;
 
  735   if (isDefaultSolveCriteria) {
 
  736     overallSolveStatus.solveStatus = SOLVE_STATUS_UNKNOWN;
 
  737     overallSolveStatus.achievedTol = SS::unknownTolerance();
 
  740     overallSolveStatus.solveStatus = solveStatus.solveStatus;
 
  741     overallSolveStatus.achievedTol = solveStatus.achievedTol;
 
  743   std::ostringstream oss;
 
  746     << ( overallSolveStatus.solveStatus==SOLVE_STATUS_CONVERGED ? 
"converged" : 
"unconverged" )
 
  747     << 
" on m = "<<m<<
" RHSs using " << totalIterations << 
" cumulative iterations" 
  748     << 
" for an average of " << (totalIterations/m) << 
" iterations/RHS and" 
  749     << 
" total CPU time of "<<totalTimer.totalElapsedTime()<<
" sec.";
 
  750   overallSolveStatus.message = oss.str();
 
  753   if (overallSolveStatus.extraParameters.is_null()) {
 
  754     overallSolveStatus.extraParameters = Teuchos::parameterList ();
 
  756   overallSolveStatus.extraParameters->set (
"AztecOO/Iteration Count",
 
  759   overallSolveStatus.extraParameters->set (
"Iteration Count",
 
  761   overallSolveStatus.extraParameters->set (
"AztecOO/Achieved Tolerance",
 
  762                                             overallSolveStatus.achievedTol);
 
  769       << 
"\nTotal solve time = "<<totalTimer.totalElapsedTime()<<
" sec\n";
 
  771   return overallSolveStatus;
 
RCP< Epetra_MultiVector > get_Epetra_MultiVector(const Epetra_Map &map, const RCP< MultiVectorBase< double > > &mv)
 
virtual bool opSupportedImpl(EOpTransp M_trans) const 
 
std::string typeName(const T &t)
 
bool is_null(const boost::shared_ptr< T > &p)
 
void copyEpetraIntoThyra(const Epetra_MultiVector &x, const Ptr< VectorBase< double > > &thyraVec) const 
 
RCP< const PreconditionerBase< double > > prec_
 
RCP< T > create_weak() const 
 
basic_FancyOStream< CharT, Traits > & o() const 
 
RCP< const LinearOpSourceBase< double > > extract_fwdOpSrc()
Extract the forward LinearOpBase<double> object so that it can be modified. 
 
RCP< AztecOO > aztecFwdSolver_
 
bool isExternalPrec() const 
Determine if the preconditioner was external or not. 
 
std::string description() const 
 
virtual bool solveSupportsSolveMeasureTypeImpl(EOpTransp M_trans, const SolveMeasureType &solveMeasureType) const 
 
virtual void applyImpl(const EOpTransp M_trans, const MultiVectorBase< double > &X, const Ptr< MultiVectorBase< double > > &Y, const double alpha, const double beta) const 
 
basic_OSTab< char > OSTab
 
virtual bool solveSupportsImpl(EOpTransp M_trans) const 
 
RCP< const VectorSpaceBase< double > > domain() const 
 
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const 
 
void start(bool reset=false)
 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
 
EpetraExt::ModelEvaluator::MPDerivative convert(const ModelEvaluatorBase::MPDerivative &derivative, const RCP< const Epetra_Map > &fnc_map, const RCP< const Epetra_Map > &var_map)
 
bool allowInexactAdjSolve_
 
virtual std::string description() const 
 
bool allowInexactFwdSolve_
 
#define TEUCHOS_MAX(x, y)
 
TypeTo implicit_cast(const TypeFrom &t)
 
double aztecSolverScalar_
 
RCP< const LinearOpBase< double > > fwdOp_
 
RCP< const LinearOpSourceBase< double > > fwdOpSrc_
 
void copyThyraIntoEpetra(const VectorBase< double > &thyraVec, Epetra_MultiVector &x) const 
 
bool nonnull(const boost::shared_ptr< T > &p)
 
void uninitialize(RCP< const LinearOpBase< double > > *fwdOp=NULL, RCP< const LinearOpSourceBase< double > > *fwdOpSrc=NULL, RCP< const PreconditionerBase< double > > *prec=NULL, bool *isExternalPrec=NULL, RCP< const LinearOpSourceBase< double > > *approxFwdOpSrc=NULL, RCP< AztecOO > *aztecFwdSolver=NULL, bool *allowInexactFwdSolve=NULL, RCP< AztecOO > *aztecAdjSolver=NULL, bool *allowInexactAdjSolve=NULL, double *aztecSolverScalar=NULL)
Uninitialize. 
 
RCP< const PreconditionerBase< double > > extract_prec()
Extract the preconditioner. 
 
double totalElapsedTime(bool readCurrentTime=false) const 
 
#define TEUCHOS_ASSERT(assertion_test)
 
void initialize(const RCP< const LinearOpBase< double > > &fwdOp, const RCP< const LinearOpSourceBase< double > > &fwdOpSrc, const RCP< const PreconditionerBase< double > > &prec, const bool isExternalPrec, const RCP< const LinearOpSourceBase< double > > &approxFwdOpSrc, const RCP< AztecOO > &aztecFwdSolver, const bool allowInexactFwdSolve=false, const RCP< AztecOO > &aztecAdjSolver=Teuchos::null, const bool allowInexactAdjSolve=false, const double aztecSolverScalar=1.0)
Sets up this object. 
 
RCP< const LinearOpBase< double > > clone() const 
 
RCP< const LinearOpSourceBase< double > > approxFwdOpSrc_
 
RCP< const VectorSpaceBase< double > > range() const 
 
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
 
SolveStatus< double > solveImpl(const EOpTransp M_trans, const MultiVectorBase< double > &B, const Ptr< MultiVectorBase< double > > &X, const Ptr< const SolveCriteria< double > > solveCriteria) const 
 
RCP< AztecOO > aztecAdjSolver_
 
RCP< const LinearOpSourceBase< double > > extract_approxFwdOpSrc()
Extract the approximate forward LinearOpBase<double> object used to build the preconditioner. 
 
AztecOOLinearOpWithSolve(const int fwdDefaultMaxIterations=400, const double fwdDefaultTol=1e-6, const int adjDefaultMaxIterations=400, const double adjDefaultTol=1e-6, const bool outputEveryRhs=false)