DOpE
newtonsolver.h
Go to the documentation of this file.
1 
24 #ifndef _NEWTON_SOLVER_H_
25 #define _NEWTON_SOLVER_H_
26 
27 #include <lac/vector.h>
28 #include <lac/block_sparsity_pattern.h>
29 #include <lac/block_sparse_matrix.h>
30 
31 #include <numerics/vector_tools.h>
32 #include <vector>
33 #include <iostream>
34 #include <fstream>
35 #include <iomanip>
36 
37 #include "parameterreader.h"
38 
39 
40 
41 
42 namespace DOpE
43 {
53  template <typename INTEGRATOR, typename LINEARSOLVER, typename VECTOR>
54  class NewtonSolver : public LINEARSOLVER
55  {
56  public:
57  NewtonSolver(INTEGRATOR &integrator, ParameterReader &param_reader);
58  ~NewtonSolver();
59 
60  static void declare_params(ParameterReader &param_reader);
61 
66  template<typename PROBLEM>
67  void ReInit(PROBLEM& pde);
68 
96  template<typename PROBLEM>
97  bool NonlinearSolve(PROBLEM& pde, VECTOR &solution, bool apply_boundary_values=true,
98  bool force_matrix_build=false,
99  int priority = 5, std::string algo_level = "\t\t ");
100 
101  protected:
102 
103  inline INTEGRATOR& GetIntegrator();
104 
105  private:
106  INTEGRATOR &_integrator;
107 
108  bool _build_matrix;
109 
110  double _nonlinear_global_tol, _nonlinear_tol, _nonlinear_rho;
111  double _linesearch_rho;
112  int _nonlinear_maxiter, _line_maxiter;
113  };
114 
115  /**********************************Implementation*******************************************/
116 
117 template <typename INTEGRATOR, typename LINEARSOLVER, typename VECTOR>
120  {
121  param_reader.SetSubsection("newtonsolver parameters");
122  param_reader.declare_entry("nonlinear_global_tol", "1.e-12",Patterns::Double(0),"global tolerance for the newton iteration");
123  param_reader.declare_entry("nonlinear_tol", "1.e-10",Patterns::Double(0),"relative tolerance for the newton iteration");
124  param_reader.declare_entry("nonlinear_maxiter", "10",Patterns::Integer(0),"maximal number of newton iterations");
125  param_reader.declare_entry("nonlinear_rho", "0.1",Patterns::Double(0),"minimal newton reduction, if actual reduction is less, matrix is rebuild ");
126 
127  param_reader.declare_entry("line_maxiter", "4",Patterns::Integer(0),"maximal number of linesearch steps");
128  param_reader.declare_entry("linesearch_rho", "0.9",Patterns::Double(0),"reduction rate for the linesearch damping paramete");
129 
130  LINEARSOLVER::declare_params(param_reader);
131  }
132 
133  /*******************************************************************************************/
134 
135  template <typename INTEGRATOR, typename LINEARSOLVER, typename VECTOR>
137  ::NewtonSolver(INTEGRATOR &integrator, ParameterReader &param_reader)
138  : LINEARSOLVER(param_reader), _integrator(integrator)
139  {
140  param_reader.SetSubsection("newtonsolver parameters");
141  _nonlinear_global_tol = param_reader.get_double ("nonlinear_global_tol");
142  _nonlinear_tol = param_reader.get_double ("nonlinear_tol");
143  _nonlinear_maxiter = param_reader.get_integer ("nonlinear_maxiter");
144  _nonlinear_rho = param_reader.get_double ("nonlinear_rho");
145 
146  _line_maxiter = param_reader.get_integer ("line_maxiter");
147  _linesearch_rho = param_reader.get_double ("linesearch_rho");
148 
149  }
150 
151  /*******************************************************************************************/
152 
153  template <typename INTEGRATOR, typename LINEARSOLVER, typename VECTOR>
156  {
157  }
158 
159  /*******************************************************************************************/
160  template <typename INTEGRATOR, typename LINEARSOLVER, typename VECTOR>
161  template<typename PROBLEM>
163  ::ReInit(PROBLEM& pde)
164  {
165  LINEARSOLVER::ReInit(pde);
166  }
167 
168  /*******************************************************************************************/
169  template <typename INTEGRATOR, typename LINEARSOLVER, typename VECTOR>
170  template<typename PROBLEM>
172  ::NonlinearSolve(PROBLEM& pde,
173  VECTOR &solution,
174  bool apply_boundary_values,
175  bool force_matrix_build,
176  int priority,
177  std::string algo_level)
178  {
179  bool build_matrix = force_matrix_build;
180  VECTOR residual;
181  VECTOR du;
182  std::stringstream out;
183  pde.GetOutputHandler()->InitNewtonOut(out);
184 
185  du.reinit(solution);
186  residual.reinit(solution);
187 
188  if(apply_boundary_values)
189  {
190  GetIntegrator().ApplyInitialBoundaryValues(pde,solution);
191  }
192  pde.GetDoFConstraints().distribute(solution);
193 
194  GetIntegrator().AddDomainData("last_newton_solution",&solution);
195 
196  GetIntegrator().ComputeNonlinearResidual(pde,residual);
197  residual *= -1.;
198 
199  pde.GetOutputHandler()->SetIterationNumber(0,"PDENewton");
200  pde.GetOutputHandler()->Write(residual,"Residual"+pde.GetType(),pde.GetDoFType());
201 
202  double res = residual.linfty_norm();
203  double firstres = res;
204  double lastres = res;
205 
206 
207  out<< algo_level << "Newton step: " <<0<<"\t Residual (abs.): "
208  <<pde.GetOutputHandler()->ZeroTolerance(res, 1.0)
209  <<"\n";
210 
211  out<< algo_level << "Newton step: " <<0<<"\t Residual (rel.): " << std::scientific << firstres/firstres;
212 
213 
214  pde.GetOutputHandler()->Write(out,priority);
215 
216  int iter=0;
217  while(res > _nonlinear_global_tol && res > firstres * _nonlinear_tol)
218  {
219  iter++;
220 
221  if(iter > _nonlinear_maxiter)
222  {
223  GetIntegrator().DeleteDomainData("last_newton_solution");
224  throw DOpEIterationException("Iteration count exceeded bounds!","NewtonSolver::NonlinearSolve");
225  }
226 
227  pde.GetOutputHandler()->SetIterationNumber(iter,"PDENewton");
228 
229  LINEARSOLVER::Solve(pde,GetIntegrator(),residual,du,build_matrix);
230  build_matrix = false;
231 
232  //Linesearch
233  {
234  solution += du;
235  GetIntegrator().ComputeNonlinearResidual(pde,residual);
236  residual *= -1.;
237 
238  pde.GetOutputHandler()->Write(residual,"Residual"+pde.GetType(),pde.GetDoFType());
239  pde.GetOutputHandler()->Write(du,"Update"+pde.GetType(),pde.GetDoFType());
240 
241  double newres = residual.linfty_norm();
242  int lineiter=0;
243  double rho = _linesearch_rho;
244  double alpha=1;
245 
246  while(newres > res)
247  {
248  out<< algo_level << "Newton step: " <<iter<<"\t Residual (rel.): "
249  <<pde.GetOutputHandler()->ZeroTolerance(newres/firstres, 1.0)
250  << "\t LineSearch {"<<lineiter<<"} ";
251 
252  pde.GetOutputHandler()->Write(out,priority+1);
253 
254  lineiter++;
255  if(lineiter > _line_maxiter)
256  {
257  GetIntegrator().DeleteDomainData("last_newton_solution");
258  throw DOpEIterationException("Line-Iteration count exceeded bounds!","NewtonSolver::NonlinearSolve");
259  }
260  solution.add(alpha*(rho-1.),du);
261  alpha*= rho;
262 
263  GetIntegrator().ComputeNonlinearResidual(pde,residual);
264  residual *= -1.;
265  pde.GetOutputHandler()->Write(residual,"Residual"+pde.GetType(),pde.GetDoFType());
266 
267  newres = residual.linfty_norm();
268 
269  }
270  if(res/lastres > _nonlinear_rho)
271  {
272  build_matrix=true;
273  }
274  lastres=res;
275  res=newres;
276 
277  out << algo_level
278  << "Newton step: "
279  <<iter
280  <<"\t Residual (rel.): "
281  << pde.GetOutputHandler()->ZeroTolerance(res/firstres, 1.0)
282  << "\t LineSearch {"
283  <<lineiter
284  <<"} ";
285 
286 
287  pde.GetOutputHandler()->Write(out,priority);
288 
289  }//End of Linesearch
290  }
291  GetIntegrator().DeleteDomainData("last_newton_solution");
292 
293  return build_matrix;
294  }
295 
296  /*******************************************************************************************/
297  template <typename INTEGRATOR, typename LINEARSOLVER, typename VECTOR>
300  {
301  return _integrator;
302  }
303 
304  /*******************************************************************************************/
305 
306 }
307 #endif
308 
309 
310 
311 
312 
Definition: dopeexception.h:76
double get_double(const std::string &entry_name)
Definition: parameterreader.h:115
void ReInit(PROBLEM &pde)
Definition: newtonsolver.h:163
Definition: parameterreader.h:36
INTEGRATOR & GetIntegrator()
Definition: newtonsolver.h:299
~NewtonSolver()
Definition: newtonsolver.h:155
void declare_entry(const std::string &entry, const std::string &default_value, const Patterns::PatternBase &pattern=Patterns::Anything(), const std::string &documentation=std::string())
Definition: parameterreader.h:98
Definition: newtonsolver.h:54
int get_integer(const std::string &entry_name)
Definition: parameterreader.h:126
NewtonSolver(INTEGRATOR &integrator, ParameterReader &param_reader)
Definition: newtonsolver.h:137
static void declare_params(ParameterReader &param_reader)
Definition: newtonsolver.h:119
void SetSubsection(const std::string subsection)
Definition: parameterreader.h:93
bool NonlinearSolve(PROBLEM &pde, VECTOR &solution, bool apply_boundary_values=true, bool force_matrix_build=false, int priority=5, std::string algo_level="\t\t ")
Definition: newtonsolver.h:172