newton_solver.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 #include "libmesh/diff_system.h"
20 #include "libmesh/dof_map.h"
22 #include "libmesh/linear_solver.h"
23 #include "libmesh/newton_solver.h"
24 #include "libmesh/numeric_vector.h"
25 #include "libmesh/sparse_matrix.h"
26 
27 namespace libMesh
28 {
29 
30 // SIGN from Numerical Recipes
31 template <typename T>
32 inline
33 T SIGN(T a, T b)
34 {
35  return b >= 0 ? std::abs(a) : -std::abs(a);
36 }
37 
39  Real last_residual,
40  Real & current_residual,
41  NumericVector<Number> & newton_iterate,
42  const NumericVector<Number> & linear_solution)
43 {
44  // Take a full step if we got a residual reduction or if we
45  // aren't substepping
46  if ((current_residual < last_residual) ||
48  (!require_finite_residual || !libmesh_isnan(current_residual))))
49  return 1.;
50 
51  // The residual vector
53 
54  Real ax = 0.; // First abscissa, don't take negative steps
55  Real cx = 1.; // Second abscissa, don't extrapolate steps
56 
57  // Find bx, a step length that gives lower residual than ax or cx
58  Real bx = 1.;
59 
60  while (libmesh_isnan(current_residual) ||
61  (current_residual > last_residual &&
63  {
64  // Reduce step size to 1/2, 1/4, etc.
65  Real substepdivision;
66  if (brent_line_search && !libmesh_isnan(current_residual))
67  {
68  substepdivision = std::min(0.5, last_residual/current_residual);
69  substepdivision = std::max(substepdivision, tol*2.);
70  }
71  else
72  substepdivision = 0.5;
73 
74  newton_iterate.add (bx * (1.-substepdivision),
75  linear_solution);
76  newton_iterate.close();
77  bx *= substepdivision;
78  if (verbose)
79  libMesh::out << " Shrinking Newton step to "
80  << bx << std::endl;
81 
82  // We may need to localize a parallel solution
83  _system.update();
84 
85  // Check residual with fractional Newton step
86  _system.assembly (true, false);
87 
88  rhs.close();
89  current_residual = rhs.l2_norm();
90  if (verbose)
91  libMesh::out << " Current Residual: "
92  << current_residual << std::endl;
93 
94  if (bx/2. < minsteplength &&
95  (libmesh_isnan(current_residual) ||
96  (current_residual > last_residual)))
97  {
98  libMesh::out << "Inexact Newton step FAILED at step "
99  << _outer_iterations << std::endl;
100 
102  {
103  libmesh_convergence_failure();
104  }
105  else
106  {
107  libMesh::out << "Continuing anyway ..." << std::endl;
109  return bx;
110  }
111  }
112  } // end while (current_residual > last_residual)
113 
114  // Now return that reduced-residual step, or use Brent's method to
115  // find a more optimal step.
116 
117  if (!brent_line_search)
118  return bx;
119 
120  // Brent's method adapted from Numerical Recipes in C, ch. 10.2
121  Real e = 0.;
122 
123  Real x = bx, w = bx, v = bx;
124 
125  // Residuals at bx
126  Real fx = current_residual,
127  fw = current_residual,
128  fv = current_residual;
129 
130  // Max iterations for Brent's method loop
131  const unsigned int max_i = 20;
132 
133  // for golden ratio steps
134  const Real golden_ratio = 1.-(std::sqrt(5.)-1.)/2.;
135 
136  for (unsigned int i=1; i <= max_i; i++)
137  {
138  Real xm = (ax+cx)*0.5;
139  Real tol1 = tol * std::abs(x) + tol*tol;
140  Real tol2 = 2.0 * tol1;
141 
142  // Test if we're done
143  if (std::abs(x-xm) <= (tol2 - 0.5 * (cx - ax)))
144  return x;
145 
146  Real d;
147 
148  // Construct a parabolic fit
149  if (std::abs(e) > tol1)
150  {
151  Real r = (x-w)*(fx-fv);
152  Real q = (x-v)*(fx-fw);
153  Real p = (x-v)*q-(x-w)*r;
154  q = 2. * (q-r);
155  if (q > 0.)
156  p = -p;
157  else
158  q = std::abs(q);
159  if (std::abs(p) >= std::abs(0.5*q*e) ||
160  p <= q * (ax-x) ||
161  p >= q * (cx-x))
162  {
163  // Take a golden section step
164  e = x >= xm ? ax-x : cx-x;
165  d = golden_ratio * e;
166  }
167  else
168  {
169  // Take a parabolic fit step
170  d = p/q;
171  if (x+d-ax < tol2 || cx-(x+d) < tol2)
172  d = SIGN(tol1, xm - x);
173  }
174  }
175  else
176  {
177  // Take a golden section step
178  e = x >= xm ? ax-x : cx-x;
179  d = golden_ratio * e;
180  }
181 
182  Real u = std::abs(d) >= tol1 ? x+d : x + SIGN(tol1,d);
183 
184  // Assemble the residual at the new steplength u
185  newton_iterate.add (bx - u, linear_solution);
186  newton_iterate.close();
187  bx = u;
188  if (verbose)
189  libMesh::out << " Shrinking Newton step to "
190  << bx << std::endl;
191 
192  // We may need to localize a parallel solution
193  _system.update();
194  _system.assembly (true, false);
195 
196  rhs.close();
197  Real fu = current_residual = rhs.l2_norm();
198  if (verbose)
199  libMesh::out << " Current Residual: "
200  << fu << std::endl;
201 
202  if (fu <= fx)
203  {
204  if (u >= x)
205  ax = x;
206  else
207  cx = x;
208  v = w; w = x; x = u;
209  fv = fw; fw = fx; fx = fu;
210  }
211  else
212  {
213  if (u < x)
214  ax = u;
215  else
216  cx = u;
217  if (fu <= fw || w == x)
218  {
219  v = w; w = u;
220  fv = fw; fw = fu;
221  }
222  else if (fu <= fv || v == x || v == w)
223  {
224  v = u;
225  fv = fu;
226  }
227  }
228  }
229 
230  if (!quiet)
231  libMesh::out << "Warning! Too many iterations used in Brent line search!"
232  << std::endl;
233  return bx;
234 }
235 
236 
238  : Parent(s),
239  require_residual_reduction(true),
240  require_finite_residual(true),
241  brent_line_search(true),
242  track_linear_convergence(false),
243  minsteplength(1e-5),
244  linear_tolerance_multiplier(1e-3),
245  _linear_solver(LinearSolver<Number>::build(s.comm()))
246 {
247 }
248 
249 
250 
252 {
253 }
254 
255 
256 
258 {
259  Parent::init();
260 
261  if (libMesh::on_command_line("--solver-system-names"))
262  _linear_solver->init((_system.name()+"_").c_str());
263  else
264  _linear_solver->init();
265 
266  _linear_solver->init_names(_system);
267 }
268 
269 
270 
272 {
273  Parent::reinit();
274 
275  _linear_solver->clear();
276 
277  _linear_solver->init_names(_system);
278 }
279 
280 
281 
282 unsigned int NewtonSolver::solve()
283 {
284  LOG_SCOPE("solve()", "NewtonSolver");
285 
286  // Reset any prior solve result
288 
289  NumericVector<Number> & newton_iterate = *(_system.solution);
290 
291  std::unique_ptr<NumericVector<Number>> linear_solution_ptr = newton_iterate.zero_clone();
292  NumericVector<Number> & linear_solution = *linear_solution_ptr;
293  NumericVector<Number> & rhs = *(_system.rhs);
294 
295  newton_iterate.close();
296  linear_solution.close();
297  rhs.close();
298 
299 #ifdef LIBMESH_ENABLE_CONSTRAINTS
301 #endif
302 
303  SparseMatrix<Number> & matrix = *(_system.matrix);
304 
305  // Set starting linear tolerance
306  Real current_linear_tolerance = initial_linear_tolerance;
307 
308  // Start counting our linear solver steps
309  _inner_iterations = 0;
310 
311  // Now we begin the nonlinear loop
314  {
315  // We may need to localize a parallel solution
316  _system.update();
317 
318  if (verbose)
319  libMesh::out << "Assembling the System" << std::endl;
320 
321  _system.assembly(true, true);
322  rhs.close();
323  Real current_residual = rhs.l2_norm();
324 
325  if (libmesh_isnan(current_residual))
326  {
327  libMesh::out << " Nonlinear solver DIVERGED at step "
329  << " with residual Not-a-Number"
330  << std::endl;
331  libmesh_convergence_failure();
332  continue;
333  }
334 
335  if (current_residual <= absolute_residual_tolerance)
336  {
337  if (verbose)
338  libMesh::out << "Linear solve unnecessary; residual "
339  << current_residual
340  << " meets absolute tolerance "
342  << std::endl;
343 
344  // We're not doing a solve, but other code may reuse this
345  // matrix.
346  matrix.close();
347 
349  if (current_residual == 0)
350  {
353  if (absolute_step_tolerance > 0)
355  if (relative_step_tolerance > 0)
357  }
358 
359  break;
360  }
361 
362  // Prepare to take incomplete steps
363  Real last_residual = current_residual;
364 
365  max_residual_norm = std::max (current_residual,
367 
368  // Compute the l2 norm of the whole solution
369  Real norm_total = newton_iterate.l2_norm();
370 
372 
373  if (verbose)
374  libMesh::out << "Nonlinear Residual: "
375  << current_residual << std::endl;
376 
377  // Make sure our linear tolerance is low enough
378  current_linear_tolerance = std::min (current_linear_tolerance,
379  current_residual * linear_tolerance_multiplier);
380 
381  // But don't let it be too small
382  if (current_linear_tolerance < minimum_linear_tolerance)
383  {
384  current_linear_tolerance = minimum_linear_tolerance;
385  }
386 
387  // If starting the nonlinear solve with a really good initial guess, we dont want to set an absurd linear tolerance
388  current_linear_tolerance = std::max(current_linear_tolerance, absolute_residual_tolerance / current_residual / 10.0);
389 
390  // At this point newton_iterate is the current guess, and
391  // linear_solution is now about to become the NEGATIVE of the next
392  // Newton step.
393 
394  // Our best initial guess for the linear_solution is zero!
395  linear_solution.zero();
396 
397  if (verbose)
398  libMesh::out << "Linear solve starting, tolerance "
399  << current_linear_tolerance << std::endl;
400 
401  // Solve the linear system.
402  const std::pair<unsigned int, Real> rval =
403  _linear_solver->solve (matrix, _system.request_matrix("Preconditioner"),
404  linear_solution, rhs, current_linear_tolerance,
406 
408  {
409  LinearConvergenceReason linear_c_reason = _linear_solver->get_converged_reason();
410 
411  // Check if something went wrong during the linear solve
412  if (linear_c_reason < 0)
413  {
414  // The linear solver failed somehow
416  // Print a message
417  libMesh::out << "Linear solver failed during Newton step, dropping out."
418  << std::endl;
419  break;
420  }
421  }
422 
423  // We may need to localize a parallel solution
424  _system.update ();
425  // The linear solver may not have fit our constraints exactly
426 #ifdef LIBMESH_ENABLE_CONSTRAINTS
428  (_system, &linear_solution, /* homogeneous = */ true);
429 #endif
430 
431  const unsigned int linear_steps = rval.first;
432  libmesh_assert_less_equal (linear_steps, max_linear_iterations);
433  _inner_iterations += linear_steps;
434 
435  const bool linear_solve_finished =
436  !(linear_steps == max_linear_iterations);
437 
438  if (verbose)
439  libMesh::out << "Linear solve finished, step " << linear_steps
440  << ", residual " << rval.second
441  << std::endl;
442 
443  // Compute the l2 norm of the nonlinear update
444  Real norm_delta = linear_solution.l2_norm();
445 
446  if (verbose)
447  libMesh::out << "Trying full Newton step" << std::endl;
448  // Take a full Newton step
449  newton_iterate.add (-1., linear_solution);
450  newton_iterate.close();
451 
452  if (this->linear_solution_monitor.get())
453  {
454  // Compute the l2 norm of the whole solution
455  norm_total = newton_iterate.l2_norm();
456  rhs.close();
457  (*this->linear_solution_monitor)(linear_solution, norm_delta,
458  newton_iterate, norm_total,
459  rhs, rhs.l2_norm(), _outer_iterations);
460  }
461 
462  // Check residual with full Newton step, if that's useful for determining
463  // whether to line search, whether to quit early, or whether to die after
464  // hitting our max iteration count
465  if (this->require_residual_reduction ||
466  this->require_finite_residual ||
467  _outer_iterations+1 < max_nonlinear_iterations ||
469  {
470  _system.update ();
471  _system.assembly(true, false);
472 
473  rhs.close();
474  current_residual = rhs.l2_norm();
475  if (verbose)
476  libMesh::out << " Current Residual: "
477  << current_residual << std::endl;
478 
479  // don't fiddle around if we've already converged
480  if (test_convergence(current_residual, norm_delta,
481  linear_solve_finished &&
482  current_residual <= last_residual))
483  {
484  if (!quiet)
485  print_convergence(_outer_iterations, current_residual,
486  norm_delta, linear_solve_finished &&
487  current_residual <= last_residual);
489  break; // out of _outer_iterations for loop
490  }
491  }
492 
493  // since we're not converged, backtrack if necessary
494  Real steplength =
495  this->line_search(std::sqrt(TOLERANCE),
496  last_residual, current_residual,
497  newton_iterate, linear_solution);
498  norm_delta *= steplength;
499 
500  // Check to see if backtracking failed,
501  // and break out of the nonlinear loop if so...
503  {
505  break; // out of _outer_iterations for loop
506  }
507 
509  {
510  libMesh::out << " Nonlinear solver reached maximum step "
511  << max_nonlinear_iterations << ", latest evaluated residual "
512  << current_residual << std::endl;
514  {
516  libMesh::out << " Continuing..." << std::endl;
517  }
518  else
519  {
520  libmesh_convergence_failure();
521  }
522  continue;
523  }
524 
525  // Compute the l2 norm of the whole solution
526  norm_total = newton_iterate.l2_norm();
527 
529 
530  // Print out information for the
531  // nonlinear iterations.
532  if (verbose)
533  libMesh::out << " Nonlinear step: |du|/|u| = "
534  << norm_delta / norm_total
535  << ", |du| = " << norm_delta
536  << std::endl;
537 
538  // Terminate the solution iteration if the difference between
539  // this iteration and the last is sufficiently small.
540  if (test_convergence(current_residual, norm_delta / steplength,
541  linear_solve_finished))
542  {
543  if (!quiet)
544  print_convergence(_outer_iterations, current_residual,
545  norm_delta / steplength,
546  linear_solve_finished);
548  break; // out of _outer_iterations for loop
549  }
550  } // end nonlinear loop
551 
552  // The linear solver may not have fit our constraints exactly
553 #ifdef LIBMESH_ENABLE_CONSTRAINTS
555 #endif
556 
557  // We may need to localize a parallel solution
558  _system.update ();
559 
560  // Make sure we are returning something sensible as the
561  // _solve_result, except in the edge case where we weren't really asked to
562  // solve.
563  libmesh_assert (_solve_result != DiffSolver::INVALID_SOLVE_RESULT ||
565 
566  return _solve_result;
567 }
568 
569 
570 
571 bool NewtonSolver::test_convergence(Real current_residual,
572  Real step_norm,
573  bool linear_solve_finished)
574 {
575  // We haven't converged unless we pass a convergence test
576  bool has_converged = false;
577 
578  // Is our absolute residual low enough?
579  if (current_residual < absolute_residual_tolerance)
580  {
582  has_converged = true;
583  }
584 
585  // Is our relative residual low enough?
586  if ((current_residual / max_residual_norm) <
588  {
590  has_converged = true;
591  }
592 
593  // For incomplete linear solves, it's not safe to test step sizes
594  if (!linear_solve_finished)
595  {
596  return has_converged;
597  }
598 
599  // Is our absolute Newton step size small enough?
600  if (step_norm < absolute_step_tolerance)
601  {
603  has_converged = true;
604  }
605 
606  // Is our relative Newton step size small enough?
607  if (step_norm / max_solution_norm <
609  {
611  has_converged = true;
612  }
613 
614  return has_converged;
615 }
616 
617 
618 void NewtonSolver::print_convergence(unsigned int step_num,
619  Real current_residual,
620  Real step_norm,
621  bool linear_solve_finished)
622 {
623  // Is our absolute residual low enough?
624  if (current_residual < absolute_residual_tolerance)
625  {
626  libMesh::out << " Nonlinear solver converged, step " << step_num
627  << ", residual " << current_residual
628  << std::endl;
629  }
631  {
632  if (verbose)
633  libMesh::out << " Nonlinear solver current_residual "
634  << current_residual << " > "
635  << (absolute_residual_tolerance) << std::endl;
636  }
637 
638  // Is our relative residual low enough?
639  if ((current_residual / max_residual_norm) <
641  {
642  libMesh::out << " Nonlinear solver converged, step " << step_num
643  << ", residual reduction "
644  << current_residual / max_residual_norm
645  << " < " << relative_residual_tolerance
646  << std::endl;
647  }
649  {
650  if (verbose)
651  libMesh::out << " Nonlinear solver relative residual "
652  << (current_residual / max_residual_norm)
653  << " > " << relative_residual_tolerance
654  << std::endl;
655  }
656 
657  // For incomplete linear solves, it's not safe to test step sizes
658  if (!linear_solve_finished)
659  return;
660 
661  // Is our absolute Newton step size small enough?
662  if (step_norm < absolute_step_tolerance)
663  {
664  libMesh::out << " Nonlinear solver converged, step " << step_num
665  << ", absolute step size "
666  << step_norm
667  << " < " << absolute_step_tolerance
668  << std::endl;
669  }
670  else if (absolute_step_tolerance)
671  {
672  if (verbose)
673  libMesh::out << " Nonlinear solver absolute step size "
674  << step_norm
675  << " > " << absolute_step_tolerance
676  << std::endl;
677  }
678 
679  // Is our relative Newton step size small enough?
680  if (step_norm / max_solution_norm <
682  {
683  libMesh::out << " Nonlinear solver converged, step " << step_num
684  << ", relative step size "
685  << (step_norm / max_solution_norm)
686  << " < " << relative_step_tolerance
687  << std::endl;
688  }
689  else if (relative_step_tolerance)
690  {
691  if (verbose)
692  libMesh::out << " Nonlinear solver relative step size "
693  << (step_norm / max_solution_norm)
694  << " > " << relative_step_tolerance
695  << std::endl;
696  }
697 }
698 
699 } // namespace libMesh
bool continue_after_max_iterations
Definition: diff_solver.h:175
virtual void init() override
NewtonSolver(sys_type &system)
double abs(double a)
std::unique_ptr< LinearSolver< Number > > _linear_solver
virtual std::unique_ptr< NumericVector< T > > zero_clone() const =0
unsigned int max_nonlinear_iterations
Definition: diff_solver.h:157
Real absolute_residual_tolerance
Definition: diff_solver.h:192
virtual unsigned int solve() override
NumericVector< Number > * rhs
virtual void reinit()
Definition: diff_solver.C:60
static const Real TOLERANCE
long double max(long double a, double b)
unsigned int _solve_result
Definition: diff_solver.h:327
unsigned int _inner_iterations
Definition: diff_solver.h:314
virtual void zero()=0
virtual void init()
Definition: diff_solver.C:69
T SIGN(T a, T b)
Definition: newton_solver.C:33
virtual void assembly(bool, bool, bool=false, bool=false)
sys_type & _system
Definition: diff_solver.h:319
virtual Real l2_norm() const =0
unsigned int max_linear_iterations
Definition: diff_solver.h:149
std::unique_ptr< NumericVector< Number > > solution
Definition: system.h:1523
virtual void close()=0
const SparseMatrix< Number > * request_matrix(const std::string &mat_name) const
virtual void update()
Definition: system.C:408
virtual void close()=0
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool continue_after_backtrack_failure
Definition: diff_solver.h:181
bool libmesh_isnan(float a)
virtual void reinit() override
unsigned int _outer_iterations
Definition: diff_solver.h:309
SparseMatrix< Number > * matrix
bool test_convergence(Real current_residual, Real step_norm, bool linear_solve_finished)
bool on_command_line(std::string arg)
Definition: libmesh.C:876
const std::string & name() const
Definition: system.h:2017
virtual void add(const numeric_index_type i, const T value)=0
Real line_search(Real tol, Real last_residual, Real &current_residual, NumericVector< Number > &newton_iterate, const NumericVector< Number > &linear_solution)
Definition: newton_solver.C:38
void print_convergence(unsigned int step_num, Real current_residual, Real step_norm, bool linear_solve_finished)
OStreamProxy out(std::cout)
const DofMap & get_dof_map() const
Definition: system.h:2049
long double min(long double a, double b)
Real relative_residual_tolerance
Definition: diff_solver.h:193
std::unique_ptr< LinearSolutionMonitor > linear_solution_monitor
Definition: diff_solver.h:289
void enforce_constraints_exactly(const System &system, NumericVector< Number > *v=nullptr, bool homogeneous=false) const
Manages consistently variables, degrees of freedom, coefficient vectors, and matrices for implicit sy...