Reference documentation for deal.II version 8.5.1
solver_control.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2016 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 () DEAL_II_NOEXCEPT {}
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();
159 
163  static void declare_parameters (ParameterHandler &param);
164 
168  void parse_parameters (ParameterHandler &param);
169 
191  virtual State check (const unsigned int step,
192  const double check_value);
193 
197  State last_check() const;
198 
202  double initial_value() const;
203 
208  double last_value() const;
209 
213  unsigned int last_step() const;
214 
218  unsigned int max_steps () const;
219 
223  unsigned int set_max_steps (const unsigned int);
224 
230  void set_failure_criterion (const double rel_failure_residual);
231 
236  void clear_failure_criterion ();
237 
241  double tolerance () const;
242 
246  double set_tolerance (const double);
247 
251  void enable_history_data();
252 
258  double average_reduction() const;
265  double final_reduction() const;
266 
272  double step_reduction(unsigned int step) const;
273 
277  void log_history (const bool);
278 
282  bool log_history () const;
283 
287  unsigned int log_frequency (unsigned int);
288 
292  void log_result (const bool);
293 
297  bool log_result () const;
298 
305 
306 protected:
310  unsigned int maxsteps;
311 
315  double tol;
316 
321 
325  double initial_val;
326 
330  double lvalue;
331 
335  unsigned int lstep;
336 
342 
347 
355 
360 
364  unsigned int m_log_frequency;
365 
372 
377 
384  std::vector<double> history_data;
385 };
386 
387 
402 {
403 public:
409  ReductionControl (const unsigned int maxiter = 100,
410  const double tolerance = 1.e-10,
411  const double reduce = 1.e-2,
412  const bool log_history = false,
413  const bool log_result = true);
414 
419  ReductionControl (const SolverControl &c);
420 
426 
431  virtual ~ReductionControl();
432 
436  static void declare_parameters (ParameterHandler &param);
437 
441  void parse_parameters (ParameterHandler &param);
442 
448  virtual State check (const unsigned int step,
449  const double check_value);
450 
454  double reduction () const;
455 
459  double set_reduction (const double);
460 
461 protected:
465  double reduce;
466 
471  double reduced_tol;
472 };
473 
487 {
488 public:
493  IterationNumberControl (const unsigned int maxiter = 100,
494  const double tolerance = 1e-12,
495  const bool log_history = false,
496  const bool log_result = true);
497 
503 
510 
515  virtual ~IterationNumberControl();
516 
522  virtual State check (const unsigned int step,
523  const double check_value);
524 };
525 
527 //---------------------------------------------------------------------------
528 
529 #ifndef DOXYGEN
530 
531 inline unsigned int
533 {
534  return maxsteps;
535 }
536 
537 
538 
539 inline unsigned int
540 SolverControl::set_max_steps (const unsigned int newval)
541 {
542  unsigned int old = maxsteps;
543  maxsteps = newval;
544  return old;
545 }
546 
547 
548 
549 inline void
550 SolverControl::set_failure_criterion (const double rel_failure_residual)
551 {
552  relative_failure_residual=rel_failure_residual;
553  check_failure=true;
554 }
555 
556 
557 
558 inline void
560 {
563  check_failure=false;
564 }
565 
566 
567 
568 inline double
570 {
571  return tol;
572 }
573 
574 
575 
576 inline double
577 SolverControl::set_tolerance (const double t)
578 {
579  double old = tol;
580  tol = t;
581  return old;
582 }
583 
584 
585 inline void
586 SolverControl::log_history (const bool newval)
587 {
588  m_log_history = newval;
589 }
590 
591 
592 
593 inline bool
595 {
596  return m_log_history;
597 }
598 
599 
600 inline void
601 SolverControl::log_result (const bool newval)
602 {
603  m_log_result = newval;
604 }
605 
606 
607 inline bool
609 {
610  return m_log_result;
611 }
612 
613 
614 inline double
616 {
617  return reduce;
618 }
619 
620 
621 inline double
622 ReductionControl::set_reduction (const double t)
623 {
624  double old = reduce;
625  reduce = t;
626  return old;
627 }
628 
629 #endif // DOXYGEN
630 
631 DEAL_II_NAMESPACE_CLOSE
632 
633 #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.
static void declare_parameters(ParameterHandler &param)
double step_reduction(unsigned int step) const
double final_reduction() const
void parse_parameters(ParameterHandler &param)
void clear_failure_criterion()
IterationNumberControl & operator=(const SolverControl &c)
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:541
static ::ExceptionBase & ExcHistoryDataRequired()
virtual ~ReductionControl()
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)
State last_check() const
unsigned int max_steps() const
virtual ~SolverControl()
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 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)