Reference documentation for deal.II version 8.5.1
solver_control.cc
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 #include <deal.II/base/logstream.h>
17 #include <deal.II/base/parameter_handler.h>
18 #include <deal.II/base/signaling_nan.h>
19 #include <deal.II/lac/solver_control.h>
20 
21 #include <cmath>
22 #include <sstream>
23 
24 DEAL_II_NAMESPACE_OPEN
25 
26 /*----------------------- SolverControl ---------------------------------*/
27 
28 
29 SolverControl::SolverControl (const unsigned int maxiter,
30  const double tolerance,
31  const bool m_log_history,
32  const bool m_log_result)
33  :
34  maxsteps(maxiter),
35  tol(tolerance),
36  lcheck(failure),
37  initial_val(numbers::signaling_nan<double>()),
38  lvalue(numbers::signaling_nan<double>()),
39  lstep(0),
40  check_failure(false),
41  relative_failure_residual(0),
42  failure_residual(0),
43  m_log_history(m_log_history),
44  m_log_frequency(1),
45  m_log_result(m_log_result),
46  history_data_enabled(false)
47 {}
48 
49 
50 
52 {}
53 
54 
55 
57 SolverControl::check (const unsigned int step,
58  const double check_value)
59 {
60  // if this is the first time we
61  // come here, then store the
62  // residual for later comparisons
63  if (step==0)
64  {
65  initial_val = check_value;
67  history_data.resize(maxsteps);
68  }
69 
70  if (m_log_history && ((step % m_log_frequency) == 0))
71  deallog << "Check " << step << "\t" << check_value << std::endl;
72 
73  lstep = step;
74  lvalue = check_value;
75 
76  if (step==0)
77  {
78  if (check_failure)
80 
81  if (m_log_result)
82  deallog << "Starting value " << check_value << std::endl;
83  }
84 
86  history_data[step] = check_value;
87 
88  if (check_value <= tol)
89  {
90  if (m_log_result)
91  deallog << "Convergence step " << step
92  << " value " << check_value << std::endl;
93  lcheck = success;
94  return success;
95  }
96 
97  if ((step >= maxsteps) ||
98  numbers::is_nan(check_value) ||
99  (check_failure && (check_value > failure_residual))
100  )
101  {
102  if (m_log_result)
103  deallog << "Failure step " << step
104  << " value " << check_value << std::endl;
105  lcheck = failure;
106  return failure;
107  }
108 
109  lcheck = iterate;
110  return iterate;
111 }
112 
113 
114 
117 {
118  return lcheck;
119 }
120 
121 
122 double
124 {
125  return initial_val;
126 }
127 
128 
129 double
131 {
132  return lvalue;
133 }
134 
135 
136 unsigned int
138 {
139  return lstep;
140 }
141 
142 
143 unsigned int
145 {
146  if (f==0)
147  f = 1;
148  unsigned int old = m_log_frequency;
149  m_log_frequency = f;
150  return old;
151 }
152 
153 
154 void
156 {
157  history_data_enabled = true;
158 }
159 
160 
161 double
163 {
164  if (lstep == 0)
165  return 0.;
166 
169  Assert (history_data[0] > 0., ExcInternalError());
171 
172  return std::pow(history_data[lstep]/history_data[0], 1./lstep);
173 }
174 
175 
176 
177 double
178 SolverControl::step_reduction(unsigned int step) const
179 {
182  Assert (step <=lstep, ExcIndexRange(step,1,lstep+1));
183  Assert (step>0, ExcIndexRange(step,1,lstep+1));
184 
185  return history_data[step]/history_data[@ref step_1 "step-1"];
186 }
187 
188 
189 double
191 {
192  return step_reduction(lstep);
193 }
194 
195 
196 void
198 {
199  param.declare_entry ("Max steps", "100", Patterns::Integer());
200  param.declare_entry ("Tolerance", "1.e-10", Patterns::Double());
201  param.declare_entry ("Log history", "false", Patterns::Bool());
202  param.declare_entry ("Log frequency", "1", Patterns::Integer());
203  param.declare_entry ("Log result", "true", Patterns::Bool());
204 }
205 
206 
208 {
209  set_max_steps (param.get_integer("Max steps"));
210  set_tolerance (param.get_double("Tolerance"));
211  log_history (param.get_bool("Log history"));
212  log_result (param.get_bool("Log result"));
213  log_frequency (param.get_integer("Log frequency"));
214 }
215 
216 /*----------------------- ReductionControl ---------------------------------*/
217 
218 
220  const double tol,
221  const double red,
222  const bool m_log_history,
223  const bool m_log_result)
224  :
225  SolverControl (n, tol, m_log_history, m_log_result),
226  reduce(red),
227  reduced_tol(numbers::signaling_nan<double>())
228 {}
229 
230 
232  :
233  SolverControl(c),
234  reduce(numbers::signaling_nan<double>()),
235  reduced_tol(numbers::signaling_nan<double>())
236 {
237  set_reduction(0.);
238 }
239 
240 
243 {
245  set_reduction(0.);
246  return *this;
247 }
248 
249 
251 {}
252 
253 
255 ReductionControl::check (const unsigned int step,
256  const double check_value)
257 {
258  // if this is the first time we
259  // come here, then store the
260  // residual for later comparisons
261  if (step==0)
262  {
263  initial_val = check_value;
264  reduced_tol = check_value * reduce;
265  }
266 
267  // check whether desired reduction
268  // has been achieved. also check
269  // for equality in case initial
270  // residual already was zero
271  if (check_value <= reduced_tol)
272  {
273  if (m_log_result)
274  deallog << "Convergence step " << step
275  << " value " << check_value << std::endl;
276  lstep = step;
277  lvalue = check_value;
278 
279  lcheck = success;
280  return success;
281  }
282  else
283  return SolverControl::check(step, check_value);
284 }
285 
286 
287 
288 void
290 {
292  param.declare_entry("Reduction", "1.e-2", Patterns::Double());
293 }
294 
295 
296 void
298 {
300  set_reduction (param.get_double("Reduction"));
301 }
302 
303 
304 /*---------------------- IterationNumberControl -----------------------------*/
305 
306 
308  const double tolerance,
309  const bool m_log_history,
310  const bool m_log_result)
311  :
312  SolverControl (n, tolerance, m_log_history, m_log_result) {}
313 
314 
316 {}
317 
318 
320 IterationNumberControl::check (const unsigned int step,
321  const double check_value)
322 {
323  // check whether the given number of iterations was reached, and return
324  // success in that case. Otherwise, go on to the check of the base class.
325  if (step >= this->maxsteps)
326  {
327  if (m_log_result)
328  deallog << "Convergence step " << step
329  << " value " << check_value << std::endl;
330  lstep = step;
331  lvalue = check_value;
332 
333  lcheck = success;
334  return success;
335  }
336  else
337  return SolverControl::check(step, check_value);
338 }
339 
340 DEAL_II_NAMESPACE_CLOSE
double initial_value() const
Stop iteration, goal not reached.
long int get_integer(const std::string &entry_string) const
bool log_history() const
virtual State check(const unsigned int step, const double check_value)
bool history_data_enabled
Continue iteration.
Subscriptor & operator=(const Subscriptor &)
Definition: subscriptor.cc:154
static void declare_parameters(ParameterHandler &param)
double step_reduction(unsigned int step) const
double final_reduction() const
void parse_parameters(ParameterHandler &param)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
double failure_residual
double get_double(const std::string &entry_name) const
static void declare_parameters(ParameterHandler &param)
unsigned int log_frequency(unsigned int)
Stop iteration, goal reached.
#define Assert(cond, exc)
Definition: exceptions.h:313
unsigned int m_log_frequency
double set_reduction(const double)
static ::ExceptionBase & ExcHistoryDataRequired()
virtual ~ReductionControl()
bool is_nan(const double x)
Definition: numbers.h:262
bool get_bool(const std::string &entry_name) const
double relative_failure_residual
double last_value() const
double average_reduction() const
unsigned int lstep
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
virtual ~SolverControl()
void declare_entry(const std::string &entry, const std::string &default_value, const Patterns::PatternBase &pattern=Patterns::Anything(), const std::string &documentation=std::string())
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)
static ::ExceptionBase & ExcInternalError()