Reference documentation for deal.II version 9.1.0-pre
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.md at
12 // the top level directory of deal.II.
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 
22 #include <deal.II/base/subscriptor.h>
23 
24 #include <vector>
25 
26 DEAL_II_NAMESPACE_OPEN
27 
28 class ParameterHandler;
29 
32 
65 class SolverControl : public Subscriptor
66 {
67 public:
72  enum State
73  {
75  iterate = 0,
80  };
81 
82 
83 
95  {
96  public:
97  NoConvergence(const unsigned int last_step, const double last_residual)
98  : last_step(last_step)
99  , last_residual(last_residual)
100  {}
101 
102  virtual ~NoConvergence() noexcept override = default;
103 
104  virtual void
105  print_info(std::ostream &out) const override
106  {
107  out
108  << "Iterative method reported convergence failure in step " << last_step
109  << ". The residual in the last step was " << last_residual << ".\n\n"
110  << "This error message can indicate that you have simply not allowed "
111  << "a sufficiently large number of iterations for your iterative solver "
112  << "to converge. This often happens when you increase the size of your "
113  << "problem. In such cases, the last residual will likely still be very "
114  << "small, and you can make the error go away by increasing the allowed "
115  << "number of iterations when setting up the SolverControl object that "
116  << "determines the maximal number of iterations you allow."
117  << "\n\n"
118  << "The other situation where this error may occur is when your matrix "
119  << "is not invertible (e.g., your matrix has a null-space), or if you "
120  << "try to apply the wrong solver to a matrix (e.g., using CG for a "
121  << "matrix that is not symmetric or not positive definite). In these "
122  << "cases, the residual in the last iteration is likely going to be large."
123  << std::endl;
124  }
125 
129  const unsigned int last_step;
130 
134  const double last_residual;
135  };
136 
137 
138 
150  SolverControl(const unsigned int n = 100,
151  const double tol = 1.e-10,
152  const bool log_history = false,
153  const bool log_result = true);
154 
159  virtual ~SolverControl() override = default;
160 
164  static void
166 
170  void
172 
190  virtual State
191  check(const unsigned int step, const double check_value);
192 
196  State
197  last_check() const;
198 
202  double
203  initial_value() const;
204 
209  double
210  last_value() const;
211 
215  unsigned int
216  last_step() const;
217 
221  unsigned int
222  max_steps() const;
223 
227  unsigned int
228  set_max_steps(const unsigned int);
229 
235  void
236  set_failure_criterion(const double rel_failure_residual);
237 
242  void
244 
248  double
249  tolerance() const;
250 
254  double
255  set_tolerance(const double);
256 
260  void
262 
266  const std::vector<double> &
267  get_history_data() const;
268 
274  double
275  average_reduction() const;
282  double
283  final_reduction() const;
284 
290  double
291  step_reduction(unsigned int step) const;
292 
296  void
297  log_history(const bool);
298 
302  bool
303  log_history() const;
304 
308  unsigned int
309  log_frequency(unsigned int);
310 
314  void
315  log_result(const bool);
316 
320  bool
321  log_result() const;
322 
329 
330 protected:
334  unsigned int maxsteps;
335 
339  double tol;
340 
345 
349  double initial_val;
350 
354  double lvalue;
355 
359  unsigned int lstep;
360 
366 
371 
379 
384 
388  unsigned int m_log_frequency;
389 
396 
401 
408  std::vector<double> history_data;
409 };
410 
411 
426 {
427 public:
433  ReductionControl(const unsigned int maxiter = 100,
434  const double tolerance = 1.e-10,
435  const double reduce = 1.e-2,
436  const bool log_history = false,
437  const bool log_result = true);
438 
444 
450  operator=(const SolverControl &c);
451 
456  virtual ~ReductionControl() override = default;
457 
461  static void
463 
467  void
469 
475  virtual State
476  check(const unsigned int step, const double check_value) override;
477 
481  double
482  reduction() const;
483 
487  double
488  set_reduction(const double);
489 
490 protected:
494  double reduce;
495 
500  double reduced_tol;
501 };
502 
516 {
517 public:
522  IterationNumberControl(const unsigned int maxiter = 100,
523  const double tolerance = 1e-12,
524  const bool log_history = false,
525  const bool log_result = true);
526 
532 
539  operator=(const SolverControl &c);
540 
545  virtual ~IterationNumberControl() override = default;
546 
552  virtual State
553  check(const unsigned int step, const double check_value) override;
554 };
555 
556 
571 {
572 public:
579  ConsecutiveControl(const unsigned int maxiter = 100,
580  const double tolerance = 1.e-10,
581  const unsigned int n_consecutive_iterations = 2,
582  const bool log_history = false,
583  const bool log_result = false);
584 
590 
597  operator=(const SolverControl &c);
598 
603  virtual ~ConsecutiveControl() override = default;
604 
609  virtual State
610  check(const unsigned int step, const double check_value) override;
611 
612 protected:
618 
623 };
624 
626 //---------------------------------------------------------------------------
627 
628 #ifndef DOXYGEN
629 
630 inline unsigned int
632 {
633  return maxsteps;
634 }
635 
636 
637 
638 inline unsigned int
639 SolverControl::set_max_steps(const unsigned int newval)
640 {
641  unsigned int old = maxsteps;
642  maxsteps = newval;
643  return old;
644 }
645 
646 
647 
648 inline void
649 SolverControl::set_failure_criterion(const double rel_failure_residual)
650 {
651  relative_failure_residual = rel_failure_residual;
652  check_failure = true;
653 }
654 
655 
656 
657 inline void
659 {
661  failure_residual = 0;
662  check_failure = false;
663 }
664 
665 
666 
667 inline double
669 {
670  return tol;
671 }
672 
673 
674 
675 inline double
676 SolverControl::set_tolerance(const double t)
677 {
678  double old = tol;
679  tol = t;
680  return old;
681 }
682 
683 
684 inline void
685 SolverControl::log_history(const bool newval)
686 {
687  m_log_history = newval;
688 }
689 
690 
691 
692 inline bool
694 {
695  return m_log_history;
696 }
697 
698 
699 inline void
700 SolverControl::log_result(const bool newval)
701 {
702  m_log_result = newval;
703 }
704 
705 
706 inline bool
708 {
709  return m_log_result;
710 }
711 
712 
713 inline double
715 {
716  return reduce;
717 }
718 
719 
720 inline double
721 ReductionControl::set_reduction(const double t)
722 {
723  double old = reduce;
724  reduce = t;
725  return old;
726 }
727 
728 #endif // DOXYGEN
729 
730 DEAL_II_NAMESPACE_CLOSE
731 
732 #endif
Stop iteration, goal not reached.
double step_reduction(unsigned int step) const
virtual State check(const unsigned int step, const double check_value)
bool history_data_enabled
bool log_history() const
Continue iteration.
static void declare_parameters(ParameterHandler &param)
unsigned int n_consecutive_iterations
double tolerance() const
void clear_failure_criterion()
double average_reduction() const
double initial_value() const
double failure_residual
ExceptionBase operator=(const ExceptionBase &)=delete
double last_value() const
const unsigned int last_step
unsigned int log_frequency(unsigned int)
Stop iteration, goal reached.
unsigned int max_steps() const
unsigned int m_log_frequency
double set_reduction(const double)
#define DeclException0(Exception0)
Definition: exceptions.h:385
static::ExceptionBase & ExcHistoryDataRequired()
double relative_failure_residual
unsigned int lstep
virtual ~SolverControl() override=default
double reduction() const
State last_check() const
unsigned int n_converged_iterations
double final_reduction() const
void enable_history_data()
double set_tolerance(const double)
unsigned int maxsteps
const std::vector< double > & get_history_data() const
void parse_parameters(ParameterHandler &param)
std::vector< double > history_data
virtual void print_info(std::ostream &out) const override
bool log_result() const
SolverControl(const unsigned int n=100, const double tol=1.e-10, const bool log_history=false, const bool log_result=true)
unsigned int set_max_steps(const unsigned int)
void set_failure_criterion(const double rel_failure_residual)