From: wolf Date: Wed, 4 Apr 2001 12:48:02 +0000 (+0000) Subject: Change name and style of exception thrown upon failure of linear solvers. X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=558223c77cdb067761f99fac9b66ca679e24a3f6;p=dealii-svn.git Change name and style of exception thrown upon failure of linear solvers. git-svn-id: https://svn.dealii.org/trunk@4370 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/lac/include/lac/eigen.h b/deal.II/lac/include/lac/eigen.h index c64c643b8e..c2886da3ba 100644 --- a/deal.II/lac/include/lac/eigen.h +++ b/deal.II/lac/include/lac/eigen.h @@ -361,10 +361,13 @@ EigenPower::solve (double &value, memory.free(Vr); deallog.pop(); - // Output - AssertThrow(control().last_check() == SolverControl::success, - typename Solver::ExcNoConvergence(control().last_step(), - control().last_value())); + + // in case of failure: throw + // exception + if (control().last_check() != SolverControl::success) + throw SolverControl::NoConvergence (control().last_step(), + control().last_value()); + // otherwise exit as normal } //----------------------------------------------------------------------// @@ -470,10 +473,12 @@ EigenInverse::solve (double &value, deallog.pop(); - // Output - AssertThrow(control().last_check() == SolverControl::success, - typename Solver::ExcNoConvergence(control().last_step(), - control().last_value())); + // in case of failure: throw + // exception + if (control().last_check() != SolverControl::success) + throw SolverControl::NoConvergence (control().last_step(), + control().last_value()); + // otherwise exit as normal } #endif diff --git a/deal.II/lac/include/lac/solver.h b/deal.II/lac/include/lac/solver.h index f6e646c598..f638719d1d 100644 --- a/deal.II/lac/include/lac/solver.h +++ b/deal.II/lac/include/lac/solver.h @@ -153,11 +153,6 @@ class Solver * Access to control object. */ SolverControl& control() const; - -//TODO: [WB] move this class to SolverControl to avoid forcing templatized class qualification when catching such exceptions - DeclException2(ExcNoConvergence, int, double, - << "Iteration did not converge after " << arg1 - << " steps. Final residual " << arg2); protected: diff --git a/deal.II/lac/include/lac/solver_bicgstab.h b/deal.II/lac/include/lac/solver_bicgstab.h index 9e6676879c..9ecaee3433 100644 --- a/deal.II/lac/include/lac/solver_bicgstab.h +++ b/deal.II/lac/include/lac/solver_bicgstab.h @@ -376,9 +376,12 @@ SolverBicgstab::solve(const MATRIX &A, deallog.pop(); - AssertThrow(control().last_check() == SolverControl::success, - typename Solver::ExcNoConvergence(control().last_step(), - control().last_value())); + // in case of failure: throw + // exception + if (control().last_check() != SolverControl::success) + throw SolverControl::NoConvergence (control().last_step(), + control().last_value()); + // otherwise exit as normal } diff --git a/deal.II/lac/include/lac/solver_cg.h b/deal.II/lac/include/lac/solver_cg.h index bf81acefc4..c4f4479640 100644 --- a/deal.II/lac/include/lac/solver_cg.h +++ b/deal.II/lac/include/lac/solver_cg.h @@ -282,9 +282,12 @@ SolverCG::solve (const MATRIX &A, // Output deallog.pop(); - AssertThrow(control().last_check() == SolverControl::success, - typename Solver::ExcNoConvergence(control().last_step(), - control().last_value())); + // in case of failure: throw + // exception + if (control().last_check() != SolverControl::success) + throw SolverControl::NoConvergence (control().last_step(), + control().last_value()); + // otherwise exit as normal }; diff --git a/deal.II/lac/include/lac/solver_control.h b/deal.II/lac/include/lac/solver_control.h index 839e11e2ff..a57a7103ad 100644 --- a/deal.II/lac/include/lac/solver_control.h +++ b/deal.II/lac/include/lac/solver_control.h @@ -63,6 +63,46 @@ class SolverControl : public Subscriptor enum State { iterate = 0, success, failure }; + + /** + * Class to be thrown upon + * failing convergence of an + * iterative solver, when either + * the number of iterations + * exceeds the limit or the + * residual fails to reach the + * desired limit, e.g. in the + * case of a break-down. + * + * The residual in the last + * iteration, as well as the + * iteration number of the last + * step are stored in this object + * and can be recovered upon + * catching an exception of this + * class. + */ + class NoConvergence : public std::exception + { + public: + /** + * Constructor. + */ + NoConvergence (const unsigned int last_step, + const double last_residual); + + /** + * Iteration number of the + * last step. + */ + const unsigned int last_step; + + /** + * Residual in the last step. + */ + const double last_residual; + }; + /** * Constructor. The parameters diff --git a/deal.II/lac/include/lac/solver_gmres.h b/deal.II/lac/include/lac/solver_gmres.h index cef56c0066..d8f31db828 100644 --- a/deal.II/lac/include/lac/solver_gmres.h +++ b/deal.II/lac/include/lac/solver_gmres.h @@ -398,9 +398,12 @@ unsigned int dim = 0; deallog.pop(); - AssertThrow(control().last_check() == SolverControl::success, - typename Solver::ExcNoConvergence(control().last_step(), - control().last_value())); + // in case of failure: throw + // exception + if (control().last_check() != SolverControl::success) + throw SolverControl::NoConvergence (control().last_step(), + control().last_value()); + // otherwise exit as normal }; diff --git a/deal.II/lac/include/lac/solver_minres.h b/deal.II/lac/include/lac/solver_minres.h index fa7f86742c..1103da70d6 100644 --- a/deal.II/lac/include/lac/solver_minres.h +++ b/deal.II/lac/include/lac/solver_minres.h @@ -351,9 +351,12 @@ SolverMinRes::solve (const MATRIX &A, // Output deallog.pop (); - AssertThrow(control().last_check() == SolverControl::success, - typename Solver::ExcNoConvergence(control().last_step(), - control().last_value())); + // in case of failure: throw + // exception + if (control().last_check() != SolverControl::success) + throw SolverControl::NoConvergence (control().last_step(), + control().last_value()); + // otherwise exit as normal }; diff --git a/deal.II/lac/include/lac/solver_qmrs.h b/deal.II/lac/include/lac/solver_qmrs.h index 4aa2aedf9b..dac847778b 100644 --- a/deal.II/lac/include/lac/solver_qmrs.h +++ b/deal.II/lac/include/lac/solver_qmrs.h @@ -254,21 +254,22 @@ SolverQMRS::solve (const MATRIX &A, } while (state); - // Deallocate Memory - + // Deallocate Memory memory.free(Vv); memory.free(Vp); memory.free(Vq); memory.free(Vt); memory.free(Vd); - // Output - + // Output deallog.pop(); - AssertThrow(control().last_check() == SolverControl::success, - typename Solver::ExcNoConvergence(control().last_step(), - control().last_value())); + // in case of failure: throw + // exception + if (control().last_check() != SolverControl::success) + throw SolverControl::NoConvergence (control().last_step(), + control().last_value()); + // otherwise exit as normal }; diff --git a/deal.II/lac/include/lac/solver_richardson.h b/deal.II/lac/include/lac/solver_richardson.h index 103550ec44..71b87816d6 100644 --- a/deal.II/lac/include/lac/solver_richardson.h +++ b/deal.II/lac/include/lac/solver_richardson.h @@ -205,10 +205,13 @@ SolverRichardson::solve (const MATRIX &A, memory.free(Vd); deallog.pop(); - // Output - AssertThrow(control().last_check() == SolverControl::success, - typename Solver::ExcNoConvergence(control().last_step(), - control().last_value())); + + // in case of failure: throw + // exception + if (control().last_check() != SolverControl::success) + throw SolverControl::NoConvergence (control().last_step(), + control().last_value()); + // otherwise exit as normal } @@ -251,10 +254,12 @@ SolverRichardson::Tsolve (const MATRIX &A, memory.free(Vd); deallog.pop(); - // Output - AssertThrow(control().last_check() == SolverControl::success, - typename Solver::ExcNoConvergence(control().last_step(), - control().last_value())); + // in case of failure: throw + // exception + if (control().last_check() != SolverControl::success) + throw SolverControl::NoConvergence (control().last_step(), + control().last_value()); + // otherwise exit as normal } diff --git a/deal.II/lac/source/solver_control.cc b/deal.II/lac/source/solver_control.cc index 27f671efb7..68b801fd3f 100644 --- a/deal.II/lac/source/solver_control.cc +++ b/deal.II/lac/source/solver_control.cc @@ -21,6 +21,15 @@ /*----------------------- SolverControl ---------------------------------*/ +SolverControl::NoConvergence::NoConvergence (const unsigned int last_step, + const double last_residual) + : + last_step (last_step), + last_residual (last_residual) +{}; + + + SolverControl::SolverControl (const unsigned int maxiter, const double tolerance, const bool _log_history,