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