Line data Source code
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
11 : #include "DefaultNonlinearConvergence.h"
12 : #include "FEProblemBase.h"
13 : #include "PetscSupport.h"
14 : #include "NonlinearSystemBase.h"
15 : #include "ConvergenceIterationTypes.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 :
24 : registerMooseObject("MooseApp", DefaultNonlinearConvergence);
25 :
26 : InputParameters
27 156908 : DefaultNonlinearConvergence::validParams()
28 : {
29 156908 : InputParameters params = DefaultConvergenceBase::validParams();
30 156908 : params += FEProblemSolve::feProblemDefaultConvergenceParams();
31 :
32 156908 : params.addPrivateParam<bool>("added_as_default", false);
33 :
34 156908 : params.addClassDescription("Default convergence criteria for FEProblem.");
35 :
36 156908 : return params;
37 0 : }
38 :
39 66673 : DefaultNonlinearConvergence::DefaultNonlinearConvergence(const InputParameters & parameters)
40 : : DefaultConvergenceBase(parameters),
41 66673 : _fe_problem(*getCheckedPointerParam<FEProblemBase *>("_fe_problem_base")),
42 66673 : _nl_abs_div_tol(getSharedExecutionerParam<Real>("nl_abs_div_tol")),
43 66673 : _nl_rel_div_tol(getSharedExecutionerParam<Real>("nl_div_tol")),
44 66673 : _div_threshold(std::numeric_limits<Real>::max()),
45 66673 : _nl_forced_its(getSharedExecutionerParam<unsigned int>("nl_forced_its")),
46 66673 : _nl_max_pingpong(getSharedExecutionerParam<unsigned int>("n_max_nonlinear_pingpong")),
47 66673 : _nl_current_pingpong(0)
48 : {
49 66673 : EquationSystems & es = _fe_problem.es();
50 :
51 66673 : es.parameters.set<unsigned int>("nonlinear solver maximum iterations") =
52 133346 : getSharedExecutionerParam<unsigned int>("nl_max_its");
53 66673 : es.parameters.set<unsigned int>("nonlinear solver maximum function evaluations") =
54 133346 : getSharedExecutionerParam<unsigned int>("nl_max_funcs");
55 66673 : es.parameters.set<Real>("nonlinear solver absolute residual tolerance") =
56 133346 : getSharedExecutionerParam<Real>("nl_abs_tol");
57 66673 : es.parameters.set<Real>("nonlinear solver relative residual tolerance") =
58 133346 : getSharedExecutionerParam<Real>("nl_rel_tol");
59 66673 : es.parameters.set<Real>("nonlinear solver divergence tolerance") =
60 133346 : getSharedExecutionerParam<Real>("nl_div_tol");
61 66673 : es.parameters.set<Real>("nonlinear solver absolute step tolerance") =
62 133346 : getSharedExecutionerParam<Real>("nl_abs_step_tol");
63 66673 : es.parameters.set<Real>("nonlinear solver relative step tolerance") =
64 133346 : getSharedExecutionerParam<Real>("nl_rel_step_tol");
65 66673 : }
66 :
67 : void
68 57138 : DefaultNonlinearConvergence::checkIterationType(IterationType it_type) const
69 : {
70 57138 : DefaultConvergenceBase::checkIterationType(it_type);
71 :
72 57138 : if (it_type != ConvergenceIterationTypes::NONLINEAR)
73 4 : mooseError("DefaultNonlinearConvergence can only be used with nonlinear solves.");
74 57134 : }
75 :
76 : bool
77 364950 : DefaultNonlinearConvergence::checkRelativeConvergence(const unsigned int /*it*/,
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 364950 : if (fnorm <= ref_norm * rel_tol)
85 : {
86 152217 : oss << "Converged due to relative/normalized residual norm " << fnorm / ref_norm
87 152217 : << " < relative tolerance (" << rel_tol << ")\n";
88 152217 : return true;
89 : }
90 : else
91 212733 : return false;
92 : }
93 :
94 : Convergence::MooseConvergenceStatus
95 812228 : DefaultNonlinearConvergence::checkConvergence(unsigned int iter)
96 : {
97 812228 : TIME_SECTION(_perfid_check_convergence);
98 :
99 812228 : NonlinearSystemBase & system = _fe_problem.currentNonlinearSystem();
100 812228 : MooseConvergenceStatus status = MooseConvergenceStatus::ITERATING;
101 :
102 : // Needed by ResidualReferenceConvergence
103 812228 : nonlinearConvergenceSetup();
104 :
105 812228 : SNES snes = system.getSNES();
106 :
107 : // ||u||
108 : PetscReal xnorm;
109 812228 : LibmeshPetscCallA(_fe_problem.comm().get(), SNESGetSolutionNorm(snes, &xnorm));
110 :
111 : // ||r||
112 : PetscReal fnorm;
113 812228 : LibmeshPetscCallA(_fe_problem.comm().get(), SNESGetFunctionNorm(snes, &fnorm));
114 :
115 : // ||du||
116 : PetscReal snorm;
117 812228 : LibmeshPetscCallA(_fe_problem.comm().get(), SNESGetUpdateNorm(snes, &snorm));
118 :
119 : // Get current number of function evaluations done by SNES
120 : PetscInt nfuncs;
121 812228 : 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 812228 : 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 812228 : PetscBool force_iteration = PETSC_FALSE;
132 812228 : LibmeshPetscCallA(_fe_problem.comm().get(), SNESGetForceIteration(snes, &force_iteration));
133 :
134 : // if PETSc says to force iteration, then force at least one iteration
135 812228 : if (force_iteration && !(_nl_forced_its))
136 672 : _nl_forced_its = 1;
137 :
138 : // if specified here to force iteration, but PETSc doesn't know, tell it
139 812228 : if (!force_iteration && (_nl_forced_its))
140 : {
141 116 : 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 812228 : LibmeshPetscCallA(_fe_problem.comm().get(), SNESGetFunctionDomainError(snes, &domainerror));
150 812228 : if (domainerror)
151 0 : status = MooseConvergenceStatus::DIVERGED;
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 812228 : if (iter == 0)
158 : {
159 308757 : system.setInitialResidual(fnorm);
160 308757 : fnorm_old = fnorm;
161 308757 : _nl_current_pingpong = 0;
162 : }
163 : else
164 503471 : 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 812228 : if ((_nl_current_pingpong % 2 == 1 && !(fnorm > fnorm_old)) ||
169 807542 : (_nl_current_pingpong % 2 == 0 && fnorm > fnorm_old))
170 9998 : _nl_current_pingpong += 1;
171 : else
172 802230 : _nl_current_pingpong = 0;
173 :
174 812228 : std::ostringstream oss;
175 812228 : if (fnorm != fnorm)
176 : {
177 0 : oss << "Failed to converge, residual norm is NaN\n";
178 0 : status = MooseConvergenceStatus::DIVERGED;
179 : }
180 812228 : else if ((iter >= _nl_forced_its) && fnorm < abs_tol)
181 : {
182 152092 : oss << "Converged due to residual norm " << fnorm << " < " << abs_tol << '\n';
183 152092 : status = MooseConvergenceStatus::CONVERGED;
184 : }
185 660136 : else if (nfuncs >= max_funcs)
186 : {
187 0 : oss << "Exceeded maximum number of residual evaluations: " << nfuncs << " > " << max_funcs
188 0 : << '\n';
189 0 : status = MooseConvergenceStatus::DIVERGED;
190 : }
191 660136 : else if ((iter >= _nl_forced_its) && iter && fnorm > system._last_nl_rnorm &&
192 5778 : fnorm >= _div_threshold)
193 : {
194 0 : oss << "Nonlinear solve was blowing up!\n";
195 0 : status = MooseConvergenceStatus::DIVERGED;
196 : }
197 812228 : if ((iter >= _nl_forced_its) && iter && status == MooseConvergenceStatus::ITERATING)
198 : {
199 385133 : const auto ref_residual = system.referenceResidual();
200 385133 : if (checkRelativeConvergence(iter, fnorm, ref_residual, rel_tol, abs_tol, oss))
201 154173 : status = MooseConvergenceStatus::CONVERGED;
202 230960 : else if (snorm < rel_step_tol * xnorm)
203 : {
204 0 : oss << "Converged due to small update length: " << snorm << " < " << rel_step_tol << " * "
205 0 : << xnorm << '\n';
206 0 : status = MooseConvergenceStatus::CONVERGED;
207 : }
208 230960 : else if (_nl_rel_div_tol > 0 && fnorm > ref_residual * _nl_rel_div_tol)
209 : {
210 6 : oss << "Diverged due to relative residual " << ref_residual << " > divergence tolerance "
211 6 : << _nl_rel_div_tol << " * relative residual " << ref_residual << '\n';
212 6 : status = MooseConvergenceStatus::DIVERGED;
213 : }
214 230954 : else if (_nl_abs_div_tol > 0 && fnorm > _nl_abs_div_tol)
215 : {
216 10 : oss << "Diverged due to residual " << fnorm << " > absolute divergence tolerance "
217 10 : << _nl_abs_div_tol << '\n';
218 10 : status = MooseConvergenceStatus::DIVERGED;
219 : }
220 230944 : else if (_nl_current_pingpong > _nl_max_pingpong)
221 : {
222 6 : oss << "Diverged due to maximum nonlinear residual pingpong achieved" << '\n';
223 6 : status = MooseConvergenceStatus::DIVERGED;
224 : }
225 : }
226 :
227 812228 : system._last_nl_rnorm = fnorm;
228 812228 : system._current_nl_its = static_cast<unsigned int>(iter);
229 :
230 812228 : std::string msg;
231 812228 : msg = oss.str();
232 812228 : if (_app.multiAppLevel() > 0)
233 206943 : MooseUtils::indentMessage(_app.name(), msg);
234 812228 : if (msg.length() > 0)
235 : #if !PETSC_VERSION_LESS_THAN(3, 17, 0)
236 306287 : 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 812228 : verboseOutput(oss);
242 :
243 812228 : return status;
244 812228 : }
|