Reference documentation for deal.II version 9.0.0
solver_control.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_solver_control_h
17 #define dealii_solver_control_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/subscriptor.h>
22 
23 #include <vector>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 class ParameterHandler;
28 
31 
64 class SolverControl : public Subscriptor
65 {
66 public:
67 
72  enum State
73  {
75  iterate = 0,
80  };
81 
82 
83 
95  {
96  public:
97  NoConvergence (const unsigned int last_step,
98  const double last_residual)
100  {}
101 
102  virtual ~NoConvergence () noexcept = default;
103 
104  virtual void print_info (std::ostream &out) const
105  {
106  out << "Iterative method reported convergence failure in step "
107  << last_step << ". The residual in the last step was " << last_residual
108  << ".\n\n"
109  << "This error message can indicate that you have simply not allowed "
110  << "a sufficiently large number of iterations for your iterative solver "
111  << "to converge. This often happens when you increase the size of your "
112  << "problem. In such cases, the last residual will likely still be very "
113  << "small, and you can make the error go away by increasing the allowed "
114  << "number of iterations when setting up the SolverControl object that "
115  << "determines the maximal number of iterations you allow."
116  << "\n\n"
117  << "The other situation where this error may occur is when your matrix "
118  << "is not invertible (e.g., your matrix has a null-space), or if you "
119  << "try to apply the wrong solver to a matrix (e.g., using CG for a "
120  << "matrix that is not symmetric or not positive definite). In these "
121  << "cases, the residual in the last iteration is likely going to be large."
122  << std::endl;
123  }
124 
128  const unsigned int last_step;
129 
133  const double last_residual;
134  };
135 
136 
137 
149  SolverControl (const unsigned int n = 100,
150  const double tol = 1.e-10,
151  const bool log_history = false,
152  const bool log_result = true);
153 
158  virtual ~SolverControl() = default;
159 
163  static void declare_parameters (ParameterHandler &param);
164 
168  void parse_parameters (ParameterHandler &param);
169 
187  virtual State check (const unsigned int step,
188  const double check_value);
189 
193  State last_check() const;
194 
198  double initial_value() const;
199 
204  double last_value() const;
205 
209  unsigned int last_step() const;
210 
214  unsigned int max_steps () const;
215 
219  unsigned int set_max_steps (const unsigned int);
220 
226  void set_failure_criterion (const double rel_failure_residual);
227 
232  void clear_failure_criterion ();
233 
237  double tolerance () const;
238 
242  double set_tolerance (const double);
243 
247  void enable_history_data();
248 
252  const std::vector<double> &get_history_data() const;
253 
259  double average_reduction() const;
266  double final_reduction() const;
267 
273  double step_reduction(unsigned int step) const;
274 
278  void log_history (const bool);
279 
283  bool log_history () const;
284 
288  unsigned int log_frequency (unsigned int);
289 
293  void log_result (const bool);
294 
298  bool log_result () const;
299 
306 
307 protected:
311  unsigned int maxsteps;
312 
316  double tol;
317 
322 
326  double initial_val;
327 
331  double lvalue;
332 
336  unsigned int lstep;
337 
343 
348 
356 
361 
365  unsigned int m_log_frequency;
366 
373 
378 
385  std::vector<double> history_data;
386 };
387 
388 
403 {
404 public:
410  ReductionControl (const unsigned int maxiter = 100,
411  const double tolerance = 1.e-10,
412  const double reduce = 1.e-2,
413  const bool log_history = false,
414  const bool log_result = true);
415 
420  ReductionControl (const SolverControl &c);
421 
427 
432  virtual ~ReductionControl() = default;
433 
437  static void declare_parameters (ParameterHandler &param);
438 
442  void parse_parameters (ParameterHandler &param);
443 
449  virtual State check (const unsigned int step,
450  const double check_value);
451 
455  double reduction () const;
456 
460  double set_reduction (const double);
461 
462 protected:
466  double reduce;
467 
472  double reduced_tol;
473 };
474 
488 {
489 public:
494  IterationNumberControl (const unsigned int maxiter = 100,
495  const double tolerance = 1e-12,
496  const bool log_history = false,
497  const bool log_result = true);
498 
504 
511 
516  virtual ~IterationNumberControl() = default;
517 
523  virtual State check (const unsigned int step,
524  const double check_value);
525 };
526 
527 
541 {
542 public:
549  ConsecutiveControl (const unsigned int maxiter = 100,
550  const double tolerance = 1.e-10,
551  const unsigned int n_consecutive_iterations = 2,
552  const bool log_history = false,
553  const bool log_result = false);
554 
560 
567 
572  virtual ~ConsecutiveControl() = default;
573 
578  virtual State check (const unsigned int step,
579  const double check_value);
580 
581 protected:
587 
592 
593 };
594 
596 //---------------------------------------------------------------------------
597 
598 #ifndef DOXYGEN
599 
600 inline unsigned int
602 {
603  return maxsteps;
604 }
605 
606 
607 
608 inline unsigned int
609 SolverControl::set_max_steps (const unsigned int newval)
610 {
611  unsigned int old = maxsteps;
612  maxsteps = newval;
613  return old;
614 }
615 
616 
617 
618 inline void
619 SolverControl::set_failure_criterion (const double rel_failure_residual)
620 {
621  relative_failure_residual=rel_failure_residual;
622  check_failure=true;
623 }
624 
625 
626 
627 inline void
629 {
632  check_failure=false;
633 }
634 
635 
636 
637 inline double
639 {
640  return tol;
641 }
642 
643 
644 
645 inline double
646 SolverControl::set_tolerance (const double t)
647 {
648  double old = tol;
649  tol = t;
650  return old;
651 }
652 
653 
654 inline void
655 SolverControl::log_history (const bool newval)
656 {
657  m_log_history = newval;
658 }
659 
660 
661 
662 inline bool
664 {
665  return m_log_history;
666 }
667 
668 
669 inline void
670 SolverControl::log_result (const bool newval)
671 {
672  m_log_result = newval;
673 }
674 
675 
676 inline bool
678 {
679  return m_log_result;
680 }
681 
682 
683 inline double
685 {
686  return reduce;
687 }
688 
689 
690 inline double
691 ReductionControl::set_reduction (const double t)
692 {
693  double old = reduce;
694  reduce = t;
695  return old;
696 }
697 
698 #endif // DOXYGEN
699 
700 DEAL_II_NAMESPACE_CLOSE
701 
702 #endif
double initial_value() const
Stop iteration, goal not reached.
bool log_history() const
virtual State check(const unsigned int step, const double check_value)
bool history_data_enabled
Continue iteration.
const std::vector< double > & get_history_data() const
virtual ~IterationNumberControl()=default
static void declare_parameters(ParameterHandler &param)
unsigned int n_consecutive_iterations
double step_reduction(unsigned int step) const
double final_reduction() const
void parse_parameters(ParameterHandler &param)
virtual ~ReductionControl()=default
ConsecutiveControl(const unsigned int maxiter=100, const double tolerance=1.e-10, const unsigned int n_consecutive_iterations=2, const bool log_history=false, const bool log_result=false)
void clear_failure_criterion()
ConsecutiveControl & operator=(const SolverControl &c)
IterationNumberControl & operator=(const SolverControl &c)
virtual ~SolverControl()=default
double reduction() const
double failure_residual
static void declare_parameters(ParameterHandler &param)
const unsigned int last_step
unsigned int log_frequency(unsigned int)
Stop iteration, goal reached.
unsigned int m_log_frequency
double set_reduction(const double)
#define DeclException0(Exception0)
Definition: exceptions.h:323
static ::ExceptionBase & ExcHistoryDataRequired()
virtual State check(const unsigned int step, const double check_value)
double relative_failure_residual
double last_value() const
double average_reduction() const
unsigned int lstep
double tolerance() const
unsigned int last_step() const
IterationNumberControl(const unsigned int maxiter=100, const double tolerance=1e-12, const bool log_history=false, const bool log_result=true)
unsigned int n_converged_iterations
State last_check() const
unsigned int max_steps() const
void enable_history_data()
double set_tolerance(const double)
ReductionControl(const unsigned int maxiter=100, const double tolerance=1.e-10, const double reduce=1.e-2, const bool log_history=false, const bool log_result=true)
unsigned int maxsteps
virtual ~ConsecutiveControl()=default
virtual State check(const unsigned int step, const double check_value)
void parse_parameters(ParameterHandler &param)
std::vector< double > history_data
SolverControl(const unsigned int n=100, const double tol=1.e-10, const bool log_history=false, const bool log_result=true)
bool log_result() const
unsigned int set_max_steps(const unsigned int)
virtual State check(const unsigned int step, const double check_value)
ReductionControl & operator=(const SolverControl &c)
virtual void print_info(std::ostream &out) const
void set_failure_criterion(const double rel_failure_residual)