https://mooseframework.inl.gov
DefaultNonlinearConvergence.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 // MOOSE includes
12 #include "FEProblemBase.h"
13 #include "PetscSupport.h"
14 #include "NonlinearSystemBase.h"
16 
17 #include "libmesh/equation_systems.h"
18 
19 // PETSc includes
20 #include <petsc.h>
21 #include <petscmat.h>
22 #include <petscsnes.h>
23 
25 
28 {
31 
32  params.addPrivateParam<bool>("added_as_default", false);
33 
34  params.addClassDescription("Default convergence criteria for FEProblem.");
35 
36  return params;
37 }
38 
40  : DefaultConvergenceBase(parameters),
41  _fe_problem(*getCheckedPointerParam<FEProblemBase *>("_fe_problem_base")),
42  _nl_abs_div_tol(getSharedExecutionerParam<Real>("nl_abs_div_tol")),
43  _nl_rel_div_tol(getSharedExecutionerParam<Real>("nl_div_tol")),
44  _div_threshold(std::numeric_limits<Real>::max()),
45  _nl_forced_its(getSharedExecutionerParam<unsigned int>("nl_forced_its")),
46  _nl_max_pingpong(getSharedExecutionerParam<unsigned int>("n_max_nonlinear_pingpong")),
47  _nl_current_pingpong(0)
48 {
49  EquationSystems & es = _fe_problem.es();
50 
51  es.parameters.set<unsigned int>("nonlinear solver maximum iterations") =
52  getSharedExecutionerParam<unsigned int>("nl_max_its");
53  es.parameters.set<unsigned int>("nonlinear solver maximum function evaluations") =
54  getSharedExecutionerParam<unsigned int>("nl_max_funcs");
55  es.parameters.set<Real>("nonlinear solver absolute residual tolerance") =
56  getSharedExecutionerParam<Real>("nl_abs_tol");
57  es.parameters.set<Real>("nonlinear solver relative residual tolerance") =
58  getSharedExecutionerParam<Real>("nl_rel_tol");
59  es.parameters.set<Real>("nonlinear solver divergence tolerance") =
60  getSharedExecutionerParam<Real>("nl_div_tol");
61  es.parameters.set<Real>("nonlinear solver absolute step tolerance") =
62  getSharedExecutionerParam<Real>("nl_abs_step_tol");
63  es.parameters.set<Real>("nonlinear solver relative step tolerance") =
64  getSharedExecutionerParam<Real>("nl_rel_step_tol");
65 }
66 
67 void
69 {
71 
73  mooseError("DefaultNonlinearConvergence can only be used with nonlinear solves.");
74 }
75 
76 bool
78  const Real fnorm,
79  const Real ref_norm,
80  const Real rel_tol,
81  const Real /*abs_tol*/,
82  std::ostringstream & oss)
83 {
84  if (fnorm <= ref_norm * rel_tol)
85  {
86  oss << "Converged due to relative/normalized residual norm " << fnorm / ref_norm
87  << " < relative tolerance (" << rel_tol << ")\n";
88  return true;
89  }
90  else
91  return false;
92 }
93 
96 {
97  TIME_SECTION(_perfid_check_convergence);
98 
101 
102  // Needed by ResidualReferenceConvergence
104 
105  SNES snes = system.getSNES();
106 
107  // ||u||
108  PetscReal xnorm;
109  LibmeshPetscCallA(_fe_problem.comm().get(), SNESGetSolutionNorm(snes, &xnorm));
110 
111  // ||r||
112  PetscReal fnorm;
113  LibmeshPetscCallA(_fe_problem.comm().get(), SNESGetFunctionNorm(snes, &fnorm));
114 
115  // ||du||
116  PetscReal snorm;
117  LibmeshPetscCallA(_fe_problem.comm().get(), SNESGetUpdateNorm(snes, &snorm));
118 
119  // Get current number of function evaluations done by SNES
120  PetscInt nfuncs;
121  LibmeshPetscCallA(_fe_problem.comm().get(), SNESGetNumberFunctionEvals(snes, &nfuncs));
122 
123  // Get tolerances from SNES
124  PetscReal abs_tol, rel_tol, rel_step_tol;
125  PetscInt max_its, max_funcs;
126  LibmeshPetscCallA(
127  _fe_problem.comm().get(),
128  SNESGetTolerances(snes, &abs_tol, &rel_tol, &rel_step_tol, &max_its, &max_funcs));
129 
130 #if !PETSC_VERSION_LESS_THAN(3, 8, 4)
131  PetscBool force_iteration = PETSC_FALSE;
132  LibmeshPetscCallA(_fe_problem.comm().get(), SNESGetForceIteration(snes, &force_iteration));
133 
134  // if PETSc says to force iteration, then force at least one iteration
135  if (force_iteration && !(_nl_forced_its))
136  _nl_forced_its = 1;
137 
138  // if specified here to force iteration, but PETSc doesn't know, tell it
139  if (!force_iteration && (_nl_forced_its))
140  {
141  LibmeshPetscCallA(_fe_problem.comm().get(), SNESSetForceIteration(snes, PETSC_TRUE));
142  }
143 #endif
144 
145  // See if SNESSetFunctionDomainError() has been called. Note:
146  // SNESSetFunctionDomainError() and SNESGetFunctionDomainError()
147  // were added in different releases of PETSc.
148  PetscBool domainerror;
149  LibmeshPetscCallA(_fe_problem.comm().get(), SNESGetFunctionDomainError(snes, &domainerror));
150  if (domainerror)
152 
153  Real fnorm_old;
154  // This is the first residual before any iterations have been done, but after
155  // solution-modifying objects (if any) have been imposed on the solution vector.
156  // We save it, and use it to detect convergence if system.usePreSMOResidual() == false.
157  if (iter == 0)
158  {
159  system.setInitialResidual(fnorm);
160  fnorm_old = fnorm;
162  }
163  else
164  fnorm_old = system._last_nl_rnorm;
165 
166  // Check for nonlinear residual ping-pong.
167  // Ping-pong will always start from a residual increase
168  if ((_nl_current_pingpong % 2 == 1 && !(fnorm > fnorm_old)) ||
169  (_nl_current_pingpong % 2 == 0 && fnorm > fnorm_old))
171  else
173 
174  std::ostringstream oss;
175  if (fnorm != fnorm)
176  {
177  oss << "Failed to converge, residual norm is NaN\n";
179  }
180  else if ((iter >= _nl_forced_its) && fnorm < abs_tol)
181  {
182  oss << "Converged due to residual norm " << fnorm << " < " << abs_tol << '\n';
184  }
185  else if (nfuncs >= max_funcs)
186  {
187  oss << "Exceeded maximum number of residual evaluations: " << nfuncs << " > " << max_funcs
188  << '\n';
190  }
191  else if ((iter >= _nl_forced_its) && iter && fnorm > system._last_nl_rnorm &&
192  fnorm >= _div_threshold)
193  {
194  oss << "Nonlinear solve was blowing up!\n";
196  }
197  if ((iter >= _nl_forced_its) && iter && status == MooseConvergenceStatus::ITERATING)
198  {
199  const auto ref_residual = system.referenceResidual();
200  if (checkRelativeConvergence(iter, fnorm, ref_residual, rel_tol, abs_tol, oss))
202  else if (snorm < rel_step_tol * xnorm)
203  {
204  oss << "Converged due to small update length: " << snorm << " < " << rel_step_tol << " * "
205  << xnorm << '\n';
207  }
208  else if (_nl_rel_div_tol > 0 && fnorm > ref_residual * _nl_rel_div_tol)
209  {
210  oss << "Diverged due to relative residual " << ref_residual << " > divergence tolerance "
211  << _nl_rel_div_tol << " * relative residual " << ref_residual << '\n';
213  }
214  else if (_nl_abs_div_tol > 0 && fnorm > _nl_abs_div_tol)
215  {
216  oss << "Diverged due to residual " << fnorm << " > absolute divergence tolerance "
217  << _nl_abs_div_tol << '\n';
219  }
221  {
222  oss << "Diverged due to maximum nonlinear residual pingpong achieved" << '\n';
224  }
225  }
226 
227  system._last_nl_rnorm = fnorm;
228  system._current_nl_its = static_cast<unsigned int>(iter);
229 
230  std::string msg;
231  msg = oss.str();
232  if (_app.multiAppLevel() > 0)
234  if (msg.length() > 0)
235 #if !PETSC_VERSION_LESS_THAN(3, 17, 0)
236  LibmeshPetscCallA(_fe_problem.comm().get(), PetscInfo(snes, "%s", msg.c_str()));
237 #else
238  LibmeshPetscCallA(_fe_problem.comm().get(), PetscInfo(snes, msg.c_str()));
239 #endif
240 
241  verboseOutput(oss);
242 
243  return status;
244 }
registerMooseObject("MooseApp", DefaultNonlinearConvergence)
virtual void nonlinearConvergenceSetup()
Performs setup necessary for each call to checkConvergence.
void addPrivateParam(const std::string &name, const T &value)
These method add a parameter to the InputParameters object which can be retrieved like any other para...
virtual SNES getSNES()=0
const Real _nl_rel_div_tol
Nonlinear relative divergence tolerance.
Base class for default convergence criteria.
DefaultNonlinearConvergence(const InputParameters &parameters)
unsigned int multiAppLevel() const
The MultiApp Level.
Definition: MooseApp.h:832
void verboseOutput(std::ostringstream &oss)
Outputs the stream to the console if verbose output is enabled.
Definition: Convergence.C:51
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const Parallel::Communicator & comm() const
virtual bool checkRelativeConvergence(const unsigned int it, const Real fnorm, const Real ref_norm, const Real rel_tol, const Real abs_tol, std::ostringstream &oss)
Check the relative convergence of the nonlinear solution.
virtual MooseConvergenceStatus checkConvergence(unsigned int iter) override
Returns convergence status.
const Real _nl_abs_div_tol
Nonlinear absolute divergence tolerance.
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
virtual const std::string & name() const
Get the name of the class.
Definition: MooseBase.h:57
auto max(const L &left, const R &right)
const unsigned int _nl_max_pingpong
Maximum number of nonlinear ping-pong iterations for a solve.
Nonlinear system to be solved.
MPI_Status status
void indentMessage(const std::string &prefix, std::string &message, const char *color=COLOR_CYAN, bool dont_indent_first_line=true, const std::string &post_prefix=": ")
Indents the supplied message given the prefix and color.
Definition: MooseUtils.C:734
Real referenceResidual() const
The reference residual used in relative convergence check.
NonlinearSystemBase & currentNonlinearSystem()
virtual libMesh::EquationSystems & es() override
static InputParameters feProblemDefaultConvergenceParams()
virtual void checkIterationType(IterationType) const
Perform checks related to the iteration type.
Definition: Convergence.h:48
MooseApp & _app
The MOOSE application this is associated with.
Definition: MooseBase.h:84
MooseConvergenceStatus
Status returned by calls to checkConvergence.
Definition: Convergence.h:33
PerfID _perfid_check_convergence
Performance ID for checkConvergence.
Definition: Convergence.h:77
static InputParameters validParams()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
unsigned int _nl_current_pingpong
Current number of nonlinear ping-pong iterations for the current solve.
Class for containing MooseEnum item information.
Definition: MooseEnumItem.h:18
T & set(const std::string &)
void setInitialResidual(Real r)
Record the initial residual (for later relative convergence check)
unsigned int _nl_forced_its
Number of iterations to force.
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
const Real _div_threshold
Divergence threshold value.
Default nonlinear convergence criteria for FEProblem.
void ErrorVector unsigned int
static InputParameters validParams()
virtual void checkIterationType(IterationType it_type) const override
Perform checks related to the iteration type.