DUNE PDELab (2.7)

newton.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_PDELAB_SOLVER_NEWTON_HH
4#define DUNE_PDELAB_SOLVER_NEWTON_HH
5
8
9#include <dune/pdelab/constraints/common/constraints.hh>
10#include <dune/pdelab/solver/newtonerrors.hh>
11#include <dune/pdelab/solver/linesearch.hh>
12#include <dune/pdelab/solver/terminate.hh>
13#include <dune/pdelab/solver/utility.hh>
14
15namespace Dune::PDELab
16{
17 namespace Impl
18 {
19 // Some SFinae magic to execute setReuse(bool) on a backend
20 template<typename T1, typename = void>
21 struct HasSetReuse
22 : std::false_type
23 {};
24
25 template<typename T>
26 struct HasSetReuse<T, decltype(std::declval<T>().setReuse(true), void())>
27 : std::true_type
28 {};
29
30 template<typename T>
31 inline void setLinearSystemReuse(T& solver_backend, bool reuse, std::true_type)
32 {
33 if (!solver_backend.getReuse() && reuse)
34 dwarn << "WARNING: Newton needed to override your choice to reuse the linear system in order to work!" << std::endl;
35 solver_backend.setReuse(reuse);
36 }
37
38 template<typename T>
39 inline void setLinearSystemReuse(T&, bool, std::false_type)
40 {}
41
42 template<typename T>
43 inline void setLinearSystemReuse(T& solver_backend, bool reuse)
44 {
45 setLinearSystemReuse(solver_backend, reuse, HasSetReuse<T>());
46 }
47 }
48
49
62 template <typename GridOperator_, typename LinearSolver_>
64 {
65 public:
67 using GridOperator = GridOperator_;
68
70 using LinearSolver = LinearSolver_;
71
74
77
80
82 using Real = typename Dune::FieldTraits<typename Domain::ElementType>::real_type;
83
85 using Result = PDESolverResult<Real>;
86
89
91 const Result& result() const
92 {
93 if (not _resultValid)
94 DUNE_THROW(NewtonError, "NewtonMethod::result() called before NewtonMethod::solve()");
95 return _result;
96 }
97
98 virtual void prepareStep(Domain& solution)
99 {
100 _reassembled = false;
101 if (_result.defect/_previousDefect > _reassembleThreshold){
102 if (_hangingNodeModifications){
103 auto dirichletValues = solution;
104 // Set all non dirichlet values to zero
105 Dune::PDELab::set_shifted_dofs(_gridOperator.localAssembler().trialConstraints(), 0.0, dirichletValues);
106 // Set all constrained DOFs to zero in solution
107 Dune::PDELab::set_constrained_dofs(_gridOperator.localAssembler().trialConstraints(), 0.0, solution);
108 // Copy correct Dirichlet values back into solution vector
109 Dune::PDELab::copy_constrained_dofs(_gridOperator.localAssembler().trialConstraints(), dirichletValues, solution);
110 // Interpolate periodic constraints / hanging nodes
111 _gridOperator.localAssembler().backtransform(solution);
112 }
113 if (_verbosity>=3)
114 std::cout << " Reassembling matrix..." << std::endl;
115 *_jacobian = Real(0.0);
116 _gridOperator.jacobian(solution, *_jacobian);
117 _reassembled = true;
118 }
119
120 _linearReduction = _minLinearReduction;
121 if (not _fixedLinearReduction){
122 // Determine maximum defect, where Newton is converged.
123 using std::min;
124 using std::max;
125 auto stop_defect = max(_result.first_defect*_reduction, _absoluteLimit);
126
127 // To achieve second order convergence of newton we need a linear
128 // reduction of at least current_defect^2/prev_defect^2. For the
129 // last newton step a linear reduction of
130 // 1/10*end_defect/current_defect is sufficient for convergence.
131 if (stop_defect/(10*_result.defect) > _result.defect*_result.defect/(_previousDefect*_previousDefect))
132 _linearReduction = stop_defect/(10*_result.defect);
133 else
134 _linearReduction = min(_minLinearReduction, _result.defect*_result.defect/(_previousDefect*_previousDefect));
135 }
136
137 if (_verbosity >= 3)
138 std::cout << " requested linear reduction: "
139 << std::setw(12) << std::setprecision(4) << std::scientific
140 << _linearReduction << std::endl;
141 }
142
143 virtual void linearSolve()
144 {
145 if (_verbosity >= 4)
146 std::cout << " Solving linear system..." << std::endl;
147
148 // If the jacobian was not reassembled we might save some work in the solver backend
149 Impl::setLinearSystemReuse(_linearSolver, not _reassembled);
150
151 // Solve the linear system
152 _correction = 0.0;
153 _linearSolver.apply(*_jacobian, _correction, _residual, _linearReduction);
154
155 if (not _linearSolver.result().converged)
156 DUNE_THROW(NewtonLinearSolverError,
157 "NewtonMethod::linearSolve(): Linear solver did not converge "
158 "in " << _linearSolver.result().iterations << " iterations");
159 if (_verbosity >= 4)
160 std::cout << " linear solver iterations: "
161 << std::setw(12) << _linearSolver.result().iterations << std::endl
162 << " linear defect reduction: "
163 << std::setw(12) << std::setprecision(4) << std::scientific
164 << _linearSolver.result().reduction << std::endl;
165 }
166
168 virtual void apply(Domain& solution)
169 {
170 // Reset solver statistics
171 _result.clear();
172 _resultValid = true;
173
174 // Store old ios flags (will be reset when this goes out of scope)
175 ios_base_all_saver restorer(std::cout);
176
177 // Prepare time measuring
178 using Clock = std::chrono::steady_clock;
179 using Duration = Clock::duration;
180 auto assembler_time = Duration::zero();
181 auto linear_solver_time = Duration::zero();
182 auto line_search_time = Duration::zero();
183 auto to_seconds =
184 [](Duration duration){
185 return std::chrono::duration<double>(duration).count();
186 };
187 auto start_solve = Clock::now();
188
189 //=========================
190 // Calculate initial defect
191 //=========================
192 updateDefect(solution);
193 _result.first_defect = _result.defect;
194 _previousDefect = _result.defect;
195
196 if (_verbosity >= 2)
197 std::cout << " Initial defect: "
198 << std::setw(12) << std::setprecision(4) << std::scientific
199 << _result.defect << std::endl;
200
201 //==========================
202 // Calculate Jacobian matrix
203 //==========================
204 if (not _jacobian)
205 _jacobian = std::make_shared<Jacobian>(_gridOperator);
206
207 //=========================
208 // Nonlinear iteration loop
209 //=========================
210 while (not _terminate->terminate()){
211 if(_verbosity >= 3)
212 std::cout << " Newton iteration " << _result.iterations
213 << " --------------------------------" << std::endl;
214
215 //=============
216 // Prepare step
217 //=============
218 auto start = Clock::now();
219 try{
220 prepareStep(solution);
221 }
222 catch (...)
223 {
224 // Keep track of statistics when the method fails. We record
225 // independently the time spent in non-converging attempts.
226 // Check OneStepMethod to see how these data are propagated.
227 auto end = Clock::now();
228 assembler_time += end-start;
229 _result.assembler_time = to_seconds(assembler_time);
230 throw;
231 }
232 auto end = Clock::now();
233 assembler_time += end -start;
234 _result.assembler_time = to_seconds(assembler_time);
235
236 // Store defect
237 _previousDefect = _result.defect;
238
239 //====================
240 // Solve linear system
241 //====================
242 start = Clock::now();
243 try{
244 linearSolve();
245 }
246 catch (...)
247 {
248 // Separately catch statistics for linear solver failures.
249 end = Clock::now();
250 linear_solver_time += end-start;
251 _result.linear_solver_time = to_seconds(linear_solver_time);
252 _result.linear_solver_iterations = _linearSolver.result().iterations;
253 throw;
254 }
255 end = Clock::now();
256 linear_solver_time += end -start;
257 _result.linear_solver_time = to_seconds(linear_solver_time);
258 _result.linear_solver_iterations = _linearSolver.result().iterations;
259
260 //===================================
261 // Do line search and update solution
262 //===================================
263 start = Clock::now();
264 _lineSearch->lineSearch(solution, _correction);
265 // lineSearch(solution);
266 end = Clock::now();
267 line_search_time += end -start;
268
269 //========================================
270 // Store statistics and create some output
271 //========================================
272 _result.reduction = _result.defect/_result.first_defect;
273 _result.iterations++;
274 _result.conv_rate = std::pow(_result.reduction, 1.0/_result.iterations);
275
276 if (_verbosity >= 3)
277 std::cout << " linear solver time: "
278 << std::setw(12) << std::setprecision(4) << std::scientific
279 << to_seconds(linear_solver_time) << std::endl
280 << " defect reduction (this iteration):"
281 << std::setw(12) << std::setprecision(4) << std::scientific
282 << _result.defect/_previousDefect << std::endl
283 << " defect reduction (total): "
284 << std::setw(12) << std::setprecision(4) << std::scientific
285 << _result.reduction << std::endl
286 << " new defect: "
287 << std::setw(12) << std::setprecision(4) << std::scientific
288 << _result.defect << std::endl;
289 if (_verbosity == 2)
290 std::cout << " Newton iteration "
291 << std::setw(2)
292 << _result.iterations
293 << ". New defect: "
294 << std::setw(12) << std::setprecision(4) << std::scientific
295 << _result.defect
296 << ". Reduction (this): "
297 << std::setw(12) << std::setprecision(4) << std::scientific
298 << _result.defect/_previousDefect
299 << ". Reduction (total): "
300 << std::setw(12) << std::setprecision(4) << std::scientific
301 << _result.reduction
302 << std::endl;
303 } // end while loop of nonlinear iterations
304
305 auto end_solve = Clock::now();
306 auto solve_time = end_solve - start_solve;
307 _result.elapsed = to_seconds(solve_time);
308
309 if (_verbosity == 1)
310 std::cout << " Newton converged after "
311 << std::setw(2)
312 << _result.iterations
313 << " iterations. Reduction: "
314 << std::setw(12) << std::setprecision(4) << std::scientific
315 << _result.reduction
316 << " (" << std::setprecision(4) << _result.elapsed << "s)"
317 << std::endl;
318
319 if (not _keepMatrix)
320 _jacobian.reset();
321 }
322
324 virtual void updateDefect(Domain& solution)
325 {
326 if (_hangingNodeModifications){
327 auto dirichletValues = solution;
328 // Set all non dirichlet values to zero
329 Dune::PDELab::set_shifted_dofs(_gridOperator.localAssembler().trialConstraints(), 0.0, dirichletValues);
330 // Set all constrained DOFs to zero in solution
331 Dune::PDELab::set_constrained_dofs(_gridOperator.localAssembler().trialConstraints(), 0.0, solution);
332 // Copy correct Dirichlet values back into solution vector
333 Dune::PDELab::copy_constrained_dofs(_gridOperator.localAssembler().trialConstraints(), dirichletValues, solution);
334 // Interpolate periodic constraints / hanging nodes
335 _gridOperator.localAssembler().backtransform(solution);
336 }
337
338 _residual = 0.0;
339 _gridOperator.residual(solution, _residual);
340
341 // Use the maximum norm as a stopping criterion. This helps loosen the tolerance
342 // when solving for stationary solutions of nonlinear time-dependent problems.
343 // The default is to use the Euclidean norm (in the else-block) as before
344 if (_useMaxNorm){
345 auto rankMax = Backend::native(_residual).infinity_norm();
346 _result.defect = _gridOperator.testGridFunctionSpace().gridView().comm().max(rankMax);
347 }
348 else
349 _result.defect = _linearSolver.norm(_residual);
350 }
351
353 void setVerbosityLevel(unsigned int verbosity)
354 {
355 if (_gridOperator.trialGridFunctionSpace().gridView().comm().rank()>0)
356 _verbosity = 0;
357 else
358 _verbosity = verbosity;
359 }
360
362 unsigned int getVerbosityLevel() const
363 {
364 return _verbosity;
365 }
366
368 void setReduction(Real reduction)
369 {
370 _reduction = reduction;
371 }
372
375 {
376 return _reduction;
377 }
378
380 void setAbsoluteLimit(Real absoluteLimit)
381 {
382 _absoluteLimit = absoluteLimit;
383 }
384
385 Real getAbsoluteLimit() const
386 {
387 return _absoluteLimit;
388 }
389
391 void setKeepMatrix(bool b)
392 {
393 _keepMatrix = b;
394 }
395
397 void setUseMaxNorm(bool b)
398 {
399 _useMaxNorm = b;
400 }
401
404 {
405 _hangingNodeModifications = b;
406 }
407
408
410 bool keepMatrix() const
411 {
412 return _keepMatrix;
413 }
414
417 {
418 if(_jacobian)
419 _jacobian.reset();
420 }
421
429 void setMinLinearReduction(Real minLinearReduction)
430 {
431 _minLinearReduction = minLinearReduction;
432 }
433
439 void setFixedLinearReduction(bool fixedLinearReduction)
440 {
441 _fixedLinearReduction = fixedLinearReduction;
442 }
443
450 void setReassembleThreshold(Real reassembleThreshold)
451 {
452 _reassembleThreshold = reassembleThreshold;
453 }
454
484 void setParameters(const ParameterTree& parameterTree){
485 _verbosity = parameterTree.get("VerbosityLevel", _verbosity);
486 _reduction = parameterTree.get("Reduction", _reduction);
487 _absoluteLimit = parameterTree.get("AbsoluteLimit", _absoluteLimit);
488 _keepMatrix = parameterTree.get("KeepMatrix", _keepMatrix);
489 _useMaxNorm = parameterTree.get("UseMaxNorm", _useMaxNorm);
490 _hangingNodeModifications = parameterTree.get("HangingNodeModifications", _hangingNodeModifications);
491 _minLinearReduction = parameterTree.get("MinLinearReduction", _minLinearReduction);
492 _fixedLinearReduction = parameterTree.get("FixedLinearReduction", _fixedLinearReduction);
493 _reassembleThreshold = parameterTree.get("ReassembleThreshold", _reassembleThreshold);
494
495 // first create the linesearch, depending on the parameter
496 std::string lineSearchStrategy = parameterTree.get("LineSearchStrategy","hackbuschReusken");
497 auto strategy = lineSearchStrategyFromString(lineSearchStrategy);
498 _lineSearch = createLineSearch(*this, strategy);
499
500 // now set parameters
501 if (parameterTree.hasSub("Terminate")){
502 _terminate->setParameters(parameterTree.sub("Terminate"));
503 }
504 else{
505 ParameterTree terminateTree;
506 terminateTree["MaxIterations"] = std::to_string(parameterTree.get("MaxIterations", 40));
507 terminateTree["ForceIteration"] = std::to_string(parameterTree.get("ForceIteration", false));
508 _terminate->setParameters(terminateTree);
509 }
510 if (parameterTree.hasSub("LineSearch")){
511 _lineSearch->setParameters(parameterTree.sub("LineSearch"));
512 }
513 else{
514 ParameterTree lineSearchTree;
515 lineSearchTree["MaxIterations"] = std::to_string(parameterTree.get("LineSearchMaxIterations", 10));
516 lineSearchTree["DampingFactor"] = std::to_string(parameterTree.get("LineSearchDampingFactor", 0.5));
517 lineSearchTree["AcceptBest"] = std::to_string(parameterTree.get("LineSearchAcceptBest", false));
518 _lineSearch->setParameters(lineSearchTree);
519 }
520 }
521
523 void setTerminate(std::shared_ptr<TerminateInterface> terminate)
524 {
525 _terminate = terminate;
526 }
527
529 std::shared_ptr<TerminateInterface> getTerminate() const
530 {
531 return _terminate;
532 }
533
538 void setLineSearch(std::shared_ptr<LineSearch> lineSearch)
539 {
540 _lineSearch = lineSearch;
541 }
542
544 std::shared_ptr<LineSearch> getLineSearch() const
545 {
546 return _lineSearch;
547 }
548
550
554 const GridOperator& gridOperator,
555 LinearSolver& linearSolver)
556 : _gridOperator(gridOperator)
557 , _linearSolver(linearSolver)
558 , _residual(gridOperator.testGridFunctionSpace())
559 , _correction(gridOperator.trialGridFunctionSpace())
560 {
561 _terminate = std::make_shared<DefaultTerminate<NewtonMethod>> (*this);
562 _lineSearch = createLineSearch(*this, LineSearchStrategy::hackbuschReusken);
563 }
564
567 const GridOperator& gridOperator,
568 LinearSolver& linearSolver,
569 const ParameterTree& parameterTree)
570 : _gridOperator(gridOperator)
571 , _linearSolver(linearSolver)
572 , _residual(gridOperator.testGridFunctionSpace())
573 , _correction(gridOperator.trialGridFunctionSpace())
574
575 {
576 _terminate = std::make_shared<DefaultTerminate<NewtonMethod>> (*this);
577 setParameters(parameterTree);
578 }
579
580 private:
581 const GridOperator& _gridOperator;
582 LinearSolver& _linearSolver;
583
584 // Vectors and Jacobi matrix we set up only once
585 Range _residual;
586 Domain _correction;
587 std::shared_ptr<Jacobian> _jacobian;
588 std::shared_ptr<Domain> _previousSolution;
589
590 std::shared_ptr<TerminateInterface> _terminate;
591 std::shared_ptr<LineSearch> _lineSearch;
592
593 // Class for storing results
594 Result _result;
595 bool _resultValid = false; // result class only valid after calling apply
596 Real _previousDefect = 0.0;
597
598 // Remember if jacobian was reassembled in prepareStep
599 bool _reassembled = true; // will be set in prepare step
600 Real _linearReduction = 0.0; // will be set in prepare step
601
602 // User parameters
603 unsigned int _verbosity = 0;
604 Real _reduction = 1e-8;
605 Real _absoluteLimit = 1e-12;
606 bool _keepMatrix = true;
607 bool _useMaxNorm = false;
608
609 // Special treatment if we have hanging nodes
610 bool _hangingNodeModifications = false;
611
612 // User parameters for prepareStep()
613 Real _minLinearReduction = 1e-3;
614 bool _fixedLinearReduction = false;
615 Real _reassembleThreshold = 0.0;
616 };
617
618} // namespace Dune::PDELab
619
620#endif
Standard grid operator implementation.
Definition: gridoperator.hh:36
Abstract base class describing the line search interface.
Definition: linesearch.hh:15
Newton solver for solving non-linear problems.
Definition: newton.hh:64
typename GridOperator::Traits::Jacobian Jacobian
Type of the Jacobian matrix.
Definition: newton.hh:79
const Result & result() const
Return results.
Definition: newton.hh:91
NewtonMethod(const GridOperator &gridOperator, LinearSolver &linearSolver, const ParameterTree &parameterTree)
Construct Newton passing a parameter tree.
Definition: newton.hh:566
void setMinLinearReduction(Real minLinearReduction)
Set the minimal reduction in the linear solver.
Definition: newton.hh:429
void setTerminate(std::shared_ptr< TerminateInterface > terminate)
Set the termination criterion.
Definition: newton.hh:523
bool keepMatrix() const
Return whether the jacobian matrix is kept across calls to apply().
Definition: newton.hh:410
unsigned int getVerbosityLevel() const
Get verbosity level.
Definition: newton.hh:362
std::shared_ptr< TerminateInterface > getTerminate() const
Return a pointer to the stored termination criterion.
Definition: newton.hh:529
typename GridOperator::Traits::Domain Domain
Type of the domain (solution)
Definition: newton.hh:73
void setParameters(const ParameterTree &parameterTree)
Interpret a parameter tree as a set of options for the newton solver.
Definition: newton.hh:484
void setHangingNodeModifications(bool b)
Does the problem have hanging nodes.
Definition: newton.hh:403
void setFixedLinearReduction(bool fixedLinearReduction)
Set wether to use a fixed reduction in the linear solver.
Definition: newton.hh:439
std::shared_ptr< LineSearch > getLineSearch() const
Return a pointer to the stored line search.
Definition: newton.hh:544
void discardMatrix()
Discard the stored Jacobian matrix.
Definition: newton.hh:416
GridOperator_ GridOperator
Type of the grid operator.
Definition: newton.hh:67
void setLineSearch(std::shared_ptr< LineSearch > lineSearch)
Set the line search.
Definition: newton.hh:538
void setReduction(Real reduction)
Set reduction Newton needs to achieve.
Definition: newton.hh:368
typename Dune::FieldTraits< typename Domain::ElementType >::real_type Real
Number type.
Definition: newton.hh:82
void setKeepMatrix(bool b)
Set whether the jacobian matrix should be kept across calls to apply().
Definition: newton.hh:391
LinearSolver_ LinearSolver
Type of the linear solver.
Definition: newton.hh:70
void setUseMaxNorm(bool b)
Set whether to use the maximum norm for stopping criteria.
Definition: newton.hh:397
typename GridOperator::Traits::Range Range
Type of the range (residual)
Definition: newton.hh:76
void setAbsoluteLimit(Real absoluteLimit)
Set absolute convergence limit.
Definition: newton.hh:380
virtual void updateDefect(Domain &solution)
Update _residual and defect in _result.
Definition: newton.hh:324
virtual void apply(Domain &solution)
Solve the nonlinear problem using solution as initial guess and for storing the result.
Definition: newton.hh:168
PDESolverResult< Real > Result
Type of results.
Definition: newton.hh:85
NewtonMethod(const GridOperator &gridOperator, LinearSolver &linearSolver)
Construct Newton using default parameters with default parameters.
Definition: newton.hh:553
void setReassembleThreshold(Real reassembleThreshold)
Set a threshold, when the linear operator is reassembled.
Definition: newton.hh:450
void setVerbosityLevel(unsigned int verbosity)
Set how much output you get.
Definition: newton.hh:353
Real getReduction() const
Get reduction.
Definition: newton.hh:374
Hierarchical structure of string parameters.
Definition: parametertree.hh:35
std::string get(const std::string &key, const std::string &defaultValue) const
get value as string
Definition: parametertree.cc:183
ParameterTree & sub(const std::string &sub)
get substructure by name
Definition: parametertree.cc:101
bool hasSub(const std::string &sub) const
test for substructure
Definition: parametertree.cc:74
Utility class for storing and resetting stream attributes.
Definition: ios_state.hh:32
A few common exception classes.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
void set_constrained_dofs(const CG &cg, typename XG::ElementType x, XG &xg)
construct constraints from given boundary condition function
Definition: constraints.hh:796
void set_shifted_dofs(const CG &cg, typename XG::ElementType x, XG &xg)
Definition: constraints.hh:1014
void copy_constrained_dofs(const CG &cg, const XG &xgin, XG &xgout)
Definition: constraints.hh:936
auto min(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::min()
Definition: defaults.hh:87
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
DWarnType dwarn(std::cerr)
Stream for warnings indicating problems.
Definition: stdstreams.hh:159
Utility class for storing and resetting stream attributes.
STL namespace.
Dune::PDELab::Backend::Matrix< MB, Domain, Range, JF > Jacobian
The type of the jacobian.
Definition: gridoperatorutilities.hh:72
Dune::PDELab::Backend::Vector< GFSV, RF > Range
The type of the range (residual).
Definition: gridoperatorutilities.hh:65
Dune::PDELab::Backend::Vector< GFSU, DF > Domain
The type of the domain (solution).
Definition: gridoperatorutilities.hh:58
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)