Stratimikos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Thyra_AztecOOLinearOpWithSolve.cpp
1 // @HEADER
2 // *****************************************************************************
3 // Stratimikos: Thyra-based strategies for linear solvers
4 //
5 // Copyright 2006 NTESS and the Stratimikos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #include "Thyra_AztecOOLinearOpWithSolve.hpp"
11 #include "Thyra_LinearOpWithSolveHelpers.hpp"
13 #include "Thyra_EpetraOperatorWrapper.hpp"
14 #include "Teuchos_BLAS_types.hpp"
15 #include "Teuchos_TimeMonitor.hpp"
16 #include "Teuchos_Time.hpp"
17 #include "Teuchos_implicit_cast.hpp"
18 
19 
20 namespace {
21 
22 #if 0
23 // unused implementation detail internal to this file
24 inline
25 Teuchos::ETransp convert( Thyra::EOpTransp trans_in )
26 {
27  Teuchos::ETransp trans_out;
28  switch(trans_in) {
29  case Thyra::NOTRANS:
30  trans_out = Teuchos::NO_TRANS;
31  break;
32  case Thyra::TRANS:
33  trans_out = Teuchos::TRANS;
34  break;
35  default:
36  TEUCHOS_TEST_FOR_EXCEPT(true); // Should never get here!
37  }
38  return trans_out;
39 }
40 #endif // 0
41 
42 
43 // This class sets some solve instance specific state and then sets it back to
44 // the default state on destruction. But using the destructor to unset the
45 // state we can be sure that the state is rest correctly even if an exception
46 // is thrown.
47 class SetAztecSolveState {
48 public:
49  SetAztecSolveState(
50  const Teuchos::RCP<AztecOO> &aztecSolver,
51  const Teuchos::RCP<Teuchos::FancyOStream> &fancyOStream,
52  const Teuchos::EVerbosityLevel verbLevel,
53  const Thyra::SolveMeasureType &solveMeasureType
54  );
55  ~SetAztecSolveState();
56 private:
57  Teuchos::RCP<AztecOO> aztecSolver_;
59  Teuchos::EVerbosityLevel verbLevel_;
60  int outputFrequency_;
61  int convergenceTest_;
62  SetAztecSolveState(); // Not defined and not to be called!
63 };
64 
65 
66 SetAztecSolveState::SetAztecSolveState(
67  const Teuchos::RCP<AztecOO> &aztecSolver,
68  const Teuchos::RCP<Teuchos::FancyOStream> &fancyOStream,
69  const Teuchos::EVerbosityLevel verbLevel,
70  const Thyra::SolveMeasureType &solveMeasureType
71  )
72  :aztecSolver_(aztecSolver.assert_not_null())
73  ,outputFrequency_(0)
74 {
75 
76  // Output state
77  verbLevel_ = verbLevel;
78  if ( Teuchos::VERB_NONE != verbLevel_ ) {
79  if (!is_null(fancyOStream)) {
80  // AztecOO puts in two tabs before it prints anything. Therefore,
81  // there is not much that we can do to improve the layout of the
82  // indentation so just leave it!
83  fancyOStream_= Teuchos::tab(
84  fancyOStream.create_weak(),
85  0, // Don't indent since AztecOO already puts in two tabs (not spaces!)
86  Teuchos::implicit_cast<std::string>("AZTECOO")
87  );
88  aztecSolver_->SetOutputStream(*fancyOStream_);
89  aztecSolver_->SetErrorStream(*fancyOStream_);
90  // Note, above we can not save the current output and error streams
91  // since AztecOO does not define functions to get them. In the
92  // future, AztecOO should define these functions if we are to avoid
93  // treading on each others print statements. However, since the
94  // AztecOO object is most likely owned by these Thyra wrappers, this
95  // should not be a problem.
96  }
97  }
98  else {
99  outputFrequency_ = aztecSolver_->GetAllAztecOptions()[AZ_output];
100  aztecSolver_->SetAztecOption(AZ_output,0);
101  }
102 
103  // Convergence test
104  convergenceTest_ = aztecSolver_->GetAztecOption(AZ_conv);
105  if (solveMeasureType.useDefault())
106  {
107  // Just use the default solve measure type already set!
108  }
109  else if (
110  solveMeasureType(
113  )
114  )
115  {
116  aztecSolver_->SetAztecOption(AZ_conv,AZ_rhs);
117  }
118  else if (
119  solveMeasureType(
122  )
123  )
124  {
125  aztecSolver_->SetAztecOption(AZ_conv,AZ_r0);
126  }
127  else {
128  TEUCHOS_TEST_FOR_EXCEPT("Invalid solve measure type, you should never get here!");
129  }
130 
131 }
132 
133 
134 SetAztecSolveState::~SetAztecSolveState()
135 {
136 
137  // Output state
138  if ( Teuchos::VERB_NONE != verbLevel_ ) {
139  if (!is_null(fancyOStream_)) {
140  aztecSolver_->SetOutputStream(std::cout);
141  aztecSolver_->SetErrorStream(std::cerr);
142  *fancyOStream_ << "\n";
143  }
144  }
145  else {
146  aztecSolver_->SetAztecOption(AZ_output,outputFrequency_);
147  }
148 
149  // Convergence test
150  aztecSolver_->SetAztecOption(AZ_conv,convergenceTest_);
151 
152 }
153 
154 
155 } // namespace
156 
157 
158 namespace Thyra {
159 
160 
161 // Constructors/initializers/accessors
162 
163 
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
170  )
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)
180 {}
181 
182 
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
189  ,const RCP<AztecOO> &aztecFwdSolver
190  ,const bool allowInexactFwdSolve
191  ,const RCP<AztecOO> &aztecAdjSolver
192  ,const bool allowInexactAdjSolve
193  ,const double aztecSolverScalar
194  )
195 {
196 #ifdef TEUCHOS_DEBUG
197  TEUCHOS_TEST_FOR_EXCEPT(fwdOp.get()==NULL);
198  TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc.get()==NULL);
199  TEUCHOS_TEST_FOR_EXCEPT(aztecFwdSolver.get()==NULL);
200 #endif
201  fwdOp_ = fwdOp;
202  fwdOpSrc_ = fwdOpSrc;
203  isExternalPrec_ = isExternalPrec_in;
204  prec_ = prec;
205  approxFwdOpSrc_ = approxFwdOpSrc;
206  aztecFwdSolver_ = aztecFwdSolver;
207  allowInexactFwdSolve_ = allowInexactFwdSolve;
208  aztecAdjSolver_ = aztecAdjSolver;
209  allowInexactAdjSolve_ = allowInexactAdjSolve;
210  aztecSolverScalar_ = aztecSolverScalar;
211  const std::string fwdOpLabel = fwdOp_->getObjectLabel();
212  if (fwdOpLabel.length())
213  this->setObjectLabel( "lows("+fwdOpLabel+")" );
214 }
215 
216 
219 {
221  _fwdOpSrc = fwdOpSrc_;
222  fwdOpSrc_ = Teuchos::null;
223  return _fwdOpSrc;
224 }
225 
226 
229 {
231  _prec = prec_;
232  prec_ = Teuchos::null;
233  return _prec;
234 }
235 
236 
238 {
239  return isExternalPrec_;
240 }
241 
242 
245 {
247  _approxFwdOpSrc = approxFwdOpSrc_;
248  approxFwdOpSrc_ = Teuchos::null;
249  return _approxFwdOpSrc;
250 }
251 
252 
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,
259  RCP<AztecOO> *aztecFwdSolver,
260  bool *allowInexactFwdSolve,
261  RCP<AztecOO> *aztecAdjSolver,
262  bool *allowInexactAdjSolve,
263  double *aztecSolverScalar
264  )
265 {
266  if (fwdOp) *fwdOp = fwdOp_;
267  if (fwdOpSrc) *fwdOpSrc = fwdOpSrc_;
268  if (prec) *prec = prec_;
269  if (isExternalPrec_inout) *isExternalPrec_inout = isExternalPrec_;
270  if (approxFwdOpSrc) *approxFwdOpSrc = approxFwdOpSrc_;
271  if (aztecFwdSolver) *aztecFwdSolver = aztecFwdSolver_;
272  if (allowInexactFwdSolve) *allowInexactFwdSolve = allowInexactFwdSolve_;
273  if (aztecAdjSolver) *aztecAdjSolver = aztecAdjSolver_;
274  if (allowInexactAdjSolve) *allowInexactAdjSolve = allowInexactAdjSolve_;
275  if (aztecSolverScalar) *aztecSolverScalar = aztecSolverScalar_;
276 
277  fwdOp_ = Teuchos::null;
278  fwdOpSrc_ = Teuchos::null;
279  prec_ = Teuchos::null;
280  isExternalPrec_ = false; // Just to make unique
281  approxFwdOpSrc_ = Teuchos::null;
282  aztecFwdSolver_ = Teuchos::null;
283  allowInexactFwdSolve_ = false;
284  aztecAdjSolver_ = Teuchos::null;
285  allowInexactAdjSolve_ = false;
286  aztecSolverScalar_ = 0.0;
287 }
288 
289 
290 // Overridden from LinearOpBase
291 
292 
295 {
296  return ( fwdOp_.get() ? fwdOp_->range() : Teuchos::null );
297 }
298 
299 
302 {
303  return ( fwdOp_.get() ? fwdOp_->domain() : Teuchos::null );
304 }
305 
306 
309 {
310  return Teuchos::null; // Not supported yet but could be
311 }
312 
313 
314 // Overridden from Teuchos::Describable
315 
316 
318 {
319  std::ostringstream oss;
321  if (fwdOp_.get()) {
322  oss << "{";
323  oss << "fwdOp="<<fwdOp_->description()<<"";
324  oss << "}";
325  }
326  return oss.str();
327 }
328 
329 
332  const Teuchos::EVerbosityLevel verbLevel
333  ) const
334 {
335  using Teuchos::OSTab;
336  using Teuchos::typeName;
337  using Teuchos::describe;
338  switch(verbLevel) {
340  case Teuchos::VERB_LOW:
341  out << this->description() << std::endl;
342  break;
344  case Teuchos::VERB_HIGH:
346  {
347  out
349  << "rangeDim=" << this->range()->dim()
350  << ",domainDim="<< this->domain()->dim() << "}\n";
351  OSTab tab(out);
352  if (!is_null(fwdOp_)) {
353  out << "fwdOp = " << describe(*fwdOp_,verbLevel);
354  }
355  if (!is_null(prec_)) {
356  out << "prec = " << describe(*prec_,verbLevel);
357  }
358  if (!is_null(aztecFwdSolver_)) {
359  if (aztecFwdSolver_->GetUserOperator())
360  out
361  << "Aztec Fwd Op = "
362  << typeName(*aztecFwdSolver_->GetUserOperator()) << "\n";
363  if (aztecFwdSolver_->GetUserMatrix())
364  out
365  << "Aztec Fwd Mat = "
366  << typeName(*aztecFwdSolver_->GetUserMatrix()) << "\n";
367  if (aztecFwdSolver_->GetPrecOperator())
368  out
369  << "Aztec Fwd Prec Op = "
370  << typeName(*aztecFwdSolver_->GetPrecOperator()) << "\n";
371  if (aztecFwdSolver_->GetPrecMatrix())
372  out
373  << "Aztec Fwd Prec Mat = "
374  << typeName(*aztecFwdSolver_->GetPrecMatrix()) << "\n";
375  }
376  if (!is_null(aztecAdjSolver_)) {
377  if (aztecAdjSolver_->GetUserOperator())
378  out
379  << "Aztec Adj Op = "
380  << typeName(*aztecAdjSolver_->GetUserOperator()) << "\n";
381  if (aztecAdjSolver_->GetUserMatrix())
382  out
383  << "Aztec Adj Mat = "
384  << typeName(*aztecAdjSolver_->GetUserMatrix()) << "\n";
385  if (aztecAdjSolver_->GetPrecOperator())
386  out
387  << "Aztec Adj Prec Op = "
388  << typeName(*aztecAdjSolver_->GetPrecOperator()) << "\n";
389  if (aztecAdjSolver_->GetPrecMatrix())
390  out
391  << "Aztec Adj Prec Mat = "
392  << typeName(*aztecAdjSolver_->GetPrecMatrix()) << "\n";
393  }
394  break;
395  }
396  default:
397  TEUCHOS_TEST_FOR_EXCEPT(true); // Should never get here!
398  }
399 }
400 
401 
402 // protected
403 
404 
405 // Overridden from LinearOpBase
406 
407 
409 {
410  return ::Thyra::opSupported(*fwdOp_,M_trans);
411 }
412 
413 
415  const EOpTransp M_trans,
416  const MultiVectorBase<double> &X,
417  const Ptr<MultiVectorBase<double> > &Y,
418  const double alpha,
419  const double beta
420  ) const
421 {
422  Thyra::apply( *fwdOp_, M_trans, X, Y, alpha, beta );
423 }
424 
425 
426 // Overridden from LinearOpWithSolveBase
427 
428 
429 bool
431 {
432  if (real_trans(M_trans)==NOTRANS) return true;
433  return (nonnull(aztecAdjSolver_));
434 }
435 
436 
437 bool
439  EOpTransp M_trans, const SolveMeasureType& solveMeasureType
440  ) const
441 {
442  if (real_trans(M_trans)==NOTRANS) {
443  if (solveMeasureType.useDefault())
444  {
445  return true;
446  }
447  else if (
448  solveMeasureType(
451  )
452  &&
453  allowInexactFwdSolve_
454  )
455  {
456  return true;
457  }
458  else if (
459  solveMeasureType(
462  )
463  &&
464  allowInexactFwdSolve_
465  )
466  {
467  return true;
468  }
469  }
470  else {
471  // TRANS
472  if (aztecAdjSolver_.get()==NULL)
473  {
474  return false;
475  }
476  else if (solveMeasureType.useDefault())
477  {
478  return true;
479  }
480  else if (
481  solveMeasureType(
484  )
485  &&
486  allowInexactFwdSolve_
487  )
488  {
489  return true;
490  }
491  else if (
492  solveMeasureType(
495  )
496  &&
497  allowInexactFwdSolve_
498  )
499  {
500  return true;
501  }
502  }
503  // If you get here then we don't support the solve measure type!
504  return false;
505 }
506 
507 
508 // Overridden from LinearOpWithSolveBase
509 
510 
513  const EOpTransp M_trans,
514  const MultiVectorBase<double> &B,
515  const Ptr<MultiVectorBase<double> > &X,
516  const Ptr<const SolveCriteria<double> > solveCriteria
517  ) const
518 {
519 
520  using Teuchos::rcp;
521  using Teuchos::rcpFromRef;
522  using Teuchos::rcpFromPtr;
523  using Teuchos::OSTab;
524  typedef SolveCriteria<double> SC;
525  typedef SolveStatus<double> SS;
526 
527  THYRA_FUNC_TIME_MONITOR("Stratimikos: AztecOOLOWS");
528  Teuchos::Time totalTimer(""), timer("");
529  totalTimer.start(true);
530 
531  RCP<Teuchos::FancyOStream> out = this->getOStream();
532  Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
533  OSTab tab = this->getOSTab();
534  if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_NONE))
535  *out << "\nSolving block system using AztecOO ...\n\n";
536 
537  //
538  // Validate input
539  //
540  TEUCHOS_ASSERT(this->solveSupportsImpl(M_trans));
541  SolveMeasureType solveMeasureType;
542  if (nonnull(solveCriteria)) {
543  solveMeasureType = solveCriteria->solveMeasureType;
544  assertSupportsSolveMeasureType(*this, M_trans, solveMeasureType);
545  }
546 
547  //
548  // Get the transpose argument
549  //
550  const EOpTransp aztecOpTransp = real_trans(M_trans);
551 
552  //
553  // Get the solver, operator, and preconditioner that we will use
554  //
556  aztecSolver = ( aztecOpTransp == NOTRANS ? aztecFwdSolver_ : aztecAdjSolver_ );
557  const Epetra_Operator
558  *aztecOp = aztecSolver->GetUserOperator();
559 
560  //
561  // Get the op(...) range and domain maps
562  //
563  const Epetra_Map
564  &opRangeMap = aztecOp->OperatorRangeMap(),
565  &opDomainMap = aztecOp->OperatorDomainMap();
566 
567  //
568  // Get the convergence criteria
569  //
570  double tol = ( aztecOpTransp==NOTRANS ? fwdDefaultTol() : adjDefaultTol() );
571  int maxIterations = ( aztecOpTransp==NOTRANS
572  ? fwdDefaultMaxIterations() : adjDefaultMaxIterations() );
573  bool isDefaultSolveCriteria = true;
574  if (nonnull(solveCriteria)) {
575  if ( solveCriteria->requestedTol != SC::unspecifiedTolerance() ) {
576  tol = solveCriteria->requestedTol;
577  isDefaultSolveCriteria = false;
578  }
579  if (nonnull(solveCriteria->extraParameters)) {
580  maxIterations = solveCriteria->extraParameters->get("Maximum Iterations",maxIterations);
581  }
582  }
583 
584  //
585  // Get Epetra_MultiVector views of B and X
586  //
587 
589  RCP<Epetra_MultiVector> epetra_X;
590 
591  const EpetraOperatorWrapper* opWrapper =
592  dynamic_cast<const EpetraOperatorWrapper*>(aztecOp);
593 
594  if (opWrapper == 0) {
595  epetra_B = get_Epetra_MultiVector(opRangeMap, rcpFromRef(B));
596  epetra_X = get_Epetra_MultiVector(opDomainMap, rcpFromPtr(X));
597  }
598 
599  //
600  // Use AztecOO to solve each RHS one at a time (which is all that I can do anyway)
601  //
602 
603  int totalIterations = 0;
604  SolveStatus<double> solveStatus;
605  solveStatus.solveStatus = SOLVE_STATUS_CONVERGED;
606  solveStatus.achievedTol = -1.0;
607 
608  /* Get the number of columns in the multivector. We use Thyra
609  * functions rather than Epetra functions to do this, as we
610  * might not yet have created an Epetra multivector. - KL */
611  //const int m = epetra_B->NumVectors();
612  const int m = B.domain()->dim();
613 
614  for( int j = 0; j < m; ++j ) {
615 
616  THYRA_FUNC_TIME_MONITOR_DIFF("Stratimikos: AztecOOLOWS:SingleSolve", SingleSolve);
617 
618  //
619  // Get Epetra_Vector views of B(:,j) and X(:,j)
620  // How this is done will depend on whether we have a true Epetra operator
621  // or we are wrapping a general Thyra operator in an Epetra operator.
622  //
623 
624  // We need to declare epetra_x_j as non-const because when we have a phony
625  // Epetra operator we'll have to copy a thyra vector into it.
626  RCP<Epetra_Vector> epetra_b_j;
627  RCP<Epetra_Vector> epetra_x_j;
628 
629  if (opWrapper == 0) {
630  epetra_b_j = rcpFromRef(*const_cast<Epetra_Vector*>((*epetra_B)(j)));
631  epetra_x_j = rcpFromRef(*(*epetra_X)(j));
632  }
633  else {
634  if (is_null(epetra_b_j)) {
635  epetra_b_j = rcp(new Epetra_Vector(opRangeMap));
636  epetra_x_j = rcp(new Epetra_Vector(opDomainMap));
637  }
638  opWrapper->copyThyraIntoEpetra(*B.col(j), *epetra_b_j);
639  opWrapper->copyThyraIntoEpetra(*X->col(j), *epetra_x_j);
640  }
641 
642  //
643  // Set the RHS and LHS
644  //
645 
646  aztecSolver->SetRHS(&*epetra_b_j);
647  aztecSolver->SetLHS(&*epetra_x_j);
648 
649  //
650  // Solve the linear system
651  //
652  timer.start(true);
653  {
654  SetAztecSolveState
655  setAztecSolveState(aztecSolver,out,verbLevel,solveMeasureType);
656  aztecSolver->Iterate( maxIterations, tol );
657  // NOTE: We ignore the returned status but get it below
658  }
659  timer.stop();
660 
661  //
662  // Scale the solution
663  // (Originally, this was at the end of the loop after all columns had been
664  // processed. It's moved here because we need to do it before copying the
665  // solution back into a Thyra vector. - KL
666  //
667  if (aztecSolverScalar_ != 1.0)
668  epetra_x_j->Scale(1.0/aztecSolverScalar_);
669 
670  //
671  // If necessary, convert the solution back to a non-epetra vector
672  //
673  if (opWrapper != 0) {
674  opWrapper->copyEpetraIntoThyra(*epetra_x_j, X->col(j).ptr());
675  }
676 
677  //
678  // Set the return solve status
679  //
680 
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 << ".";
695  oss << " Total time = " << timer.totalElapsedTime() << " sec.";
696  if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_NONE) && outputEveryRhs())
697  Teuchos::OSTab(out).o() << "j="<<j<<": " << oss.str() << "\n";
698 
699  solveStatus.achievedTol = TEUCHOS_MAX(solveStatus.achievedTol, achievedTol);
700  // Note, achieveTol may actually be greater than tol due to ill conditioning and roundoff!
701 
702  totalIterations += iterations;
703 
704  solveStatus.message = oss.str();
705  if ( isDefaultSolveCriteria ) {
706  switch(solveStatus.solveStatus) {
708  // Leave overall unknown!
709  break;
711  solveStatus.solveStatus = ( converged ? SOLVE_STATUS_CONVERGED : SOLVE_STATUS_UNCONVERGED );
712  break;
714  // Leave overall unconverged!
715  break;
716  default:
717  TEUCHOS_TEST_FOR_EXCEPT(true); // Should never get here!
718  }
719  }
720  }
721 
722  aztecSolver->UnsetLHSRHS();
723 
724  //
725  // Release the Epetra_MultiVector views of X and B
726  //
727  epetra_X = Teuchos::null;
728  epetra_B = Teuchos::null;
729 
730  //
731  // Update the overall solve criteria
732  //
733  totalTimer.stop();
734  SolveStatus<double> overallSolveStatus;
735  if (isDefaultSolveCriteria) {
736  overallSolveStatus.solveStatus = SOLVE_STATUS_UNKNOWN;
737  overallSolveStatus.achievedTol = SS::unknownTolerance();
738  }
739  else {
740  overallSolveStatus.solveStatus = solveStatus.solveStatus;
741  overallSolveStatus.achievedTol = solveStatus.achievedTol;
742  }
743  std::ostringstream oss;
744  oss
745  << "AztecOO solver "
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();
751 
752  // Added these statistics following what was done for Belos
753  if (overallSolveStatus.extraParameters.is_null()) {
754  overallSolveStatus.extraParameters = Teuchos::parameterList ();
755  }
756  overallSolveStatus.extraParameters->set ("AztecOO/Iteration Count",
757  totalIterations);
758  // package independent version of the same
759  overallSolveStatus.extraParameters->set ("Iteration Count",
760  totalIterations);
761  overallSolveStatus.extraParameters->set ("AztecOO/Achieved Tolerance",
762  overallSolveStatus.achievedTol);
763 
764  //
765  // Report the overall time
766  //
767  if (out.get() && static_cast<int>(verbLevel) >= static_cast<int>(Teuchos::VERB_LOW))
768  *out
769  << "\nTotal solve time = "<<totalTimer.totalElapsedTime()<<" sec\n";
770 
771  return overallSolveStatus;
772 
773 }
774 
775 
776 } // end namespace Thyra
RCP< Epetra_MultiVector > get_Epetra_MultiVector(const Epetra_Map &map, const RCP< MultiVectorBase< double > > &mv)
virtual bool opSupportedImpl(EOpTransp M_trans) const
bool is_null(const boost::shared_ptr< T > &p)
EOpTransp
void copyEpetraIntoThyra(const Epetra_MultiVector &x, const Ptr< VectorBase< double > > &thyraVec) const
bool useDefault() const
int UnsetLHSRHS()
basic_OSTab< char > OSTab
std::string message
RCP< T > create_weak() const
const double * GetAztecStatus() const
basic_FancyOStream< CharT, Traits > & o() const
RCP< const LinearOpSourceBase< double > > extract_fwdOpSrc()
Extract the forward LinearOpBase&lt;double&gt; object so that it can be modified.
bool isExternalPrec() const
Determine if the preconditioner was external or not.
double ScaledResidual() const
EOpTransp real_trans(EOpTransp transp)
T * get() const
virtual bool solveSupportsSolveMeasureTypeImpl(EOpTransp M_trans, const SolveMeasureType &solveMeasureType) const
virtual const Epetra_Map & OperatorDomainMap() const =0
virtual void applyImpl(const EOpTransp M_trans, const MultiVectorBase< double > &X, const Ptr< MultiVectorBase< double > > &Y, const double alpha, const double beta) const
virtual bool solveSupportsImpl(EOpTransp M_trans) const
RCP< const VectorSpaceBase< double > > domain() const
virtual const Epetra_Map & OperatorRangeMap() const =0
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)
int SetLHS(Epetra_MultiVector *X)
double stop()
RCP< const VectorBase< Scalar > > col(Ordinal j) const
virtual std::string description() const
TypeTo implicit_cast(const TypeFrom &t)
Epetra_Operator * GetPrecOperator() const
virtual void setObjectLabel(const std::string &objectLabel)
int SetRHS(Epetra_MultiVector *B)
ESolveStatus solveStatus
int NumIters() const
Epetra_RowMatrix * GetPrecMatrix() const
Epetra_Operator * GetUserOperator() const
void copyThyraIntoEpetra(const VectorBase< double > &thyraVec, Epetra_MultiVector &x) const
virtual RCP< const VectorSpaceBase< Scalar > > domain() const =0
bool nonnull(const boost::shared_ptr< T > &p)
ScalarMag achievedTol
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.
Epetra_RowMatrix * GetUserMatrix() const
double totalElapsedTime(bool readCurrentTime=false) const
RCP< ParameterList > extraParameters
#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
int Iterate(long long MaxIters, double Tolerance)
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
std::string typeName(const T &t)
RCP< const LinearOpSourceBase< double > > extract_approxFwdOpSrc()
Extract the approximate forward LinearOpBase&lt;double&gt; 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)

Generated on Thu Nov 21 2024 09:22:16 for Stratimikos by doxygen 1.8.5