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 : #include "FEProblemSolve.h"
11 :
12 : #include "FEProblem.h"
13 : #include "NonlinearSystemBase.h"
14 : #include "LinearSystem.h"
15 : #include "Convergence.h"
16 : #include "Executioner.h"
17 :
18 : std::set<std::string> const FEProblemSolve::_moose_line_searches = {"contact", "project"};
19 :
20 : const std::set<std::string> &
21 286959 : FEProblemSolve::mooseLineSearches()
22 : {
23 286959 : return _moose_line_searches;
24 : }
25 :
26 : InputParameters
27 435294 : FEProblemSolve::feProblemDefaultConvergenceParams()
28 : {
29 435294 : InputParameters params = emptyInputParameters();
30 :
31 435294 : params.addParam<unsigned int>("nl_max_its", 50, "Max Nonlinear Iterations");
32 435294 : params.addParam<unsigned int>("nl_forced_its", 0, "The Number of Forced Nonlinear Iterations");
33 435294 : params.addParam<unsigned int>("nl_max_funcs", 10000, "Max Nonlinear solver function evaluations");
34 435294 : params.addParam<Real>("nl_abs_tol", 1.0e-50, "Nonlinear Absolute Tolerance");
35 435294 : params.addParam<Real>("nl_rel_tol", 1.0e-8, "Nonlinear Relative Tolerance");
36 1305882 : params.addParam<Real>(
37 : "nl_div_tol",
38 870588 : 1.0e10,
39 : "Nonlinear Relative Divergence Tolerance. A negative value disables this check.");
40 1305882 : params.addParam<Real>(
41 : "nl_abs_div_tol",
42 870588 : 1.0e50,
43 : "Nonlinear Absolute Divergence Tolerance. A negative value disables this check.");
44 435294 : params.addParam<Real>("nl_abs_step_tol", 0., "Nonlinear Absolute step Tolerance");
45 435294 : params.addParam<Real>("nl_rel_step_tol", 0., "Nonlinear Relative step Tolerance");
46 1305882 : params.addParam<unsigned int>("n_max_nonlinear_pingpong",
47 870588 : 100,
48 : "The maximum number of times the nonlinear residual can ping pong "
49 : "before requesting halting the current evaluation and requesting "
50 : "timestep cut for transient simulations");
51 :
52 435294 : params.addParamNamesToGroup(
53 : "nl_max_its nl_forced_its nl_max_funcs nl_abs_tol nl_rel_tol "
54 : "nl_rel_step_tol nl_abs_step_tol nl_div_tol nl_abs_div_tol n_max_nonlinear_pingpong",
55 : "Nonlinear Solver");
56 :
57 435294 : return params;
58 0 : }
59 :
60 : InputParameters
61 286959 : FEProblemSolve::validParams()
62 : {
63 286959 : InputParameters params = MultiSystemSolveObject::validParams();
64 286959 : params += FEProblemSolve::feProblemDefaultConvergenceParams();
65 :
66 286959 : params.addParam<std::vector<std::vector<std::string>>>(
67 : "splitting",
68 : {},
69 : "Top-level splitting defining a hierarchical decomposition into "
70 : "subsystems to help the solver. Outer-vector of this vector-of-vector parameter correspond "
71 : "to each nonlinear system.");
72 :
73 286959 : std::set<std::string> line_searches = mooseLineSearches();
74 :
75 1147836 : std::set<std::string> alias_line_searches = {"default", "none", "basic"};
76 286959 : line_searches.insert(alias_line_searches.begin(), alias_line_searches.end());
77 286959 : std::set<std::string> petsc_line_searches = Moose::PetscSupport::getPetscValidLineSearches();
78 286959 : line_searches.insert(petsc_line_searches.begin(), petsc_line_searches.end());
79 286959 : std::string line_search_string = Moose::stringify(line_searches, " ");
80 286959 : MooseEnum line_search(line_search_string, "default");
81 286959 : std::string addtl_doc_str(" (Note: none = basic)");
82 286959 : params.addParam<MooseEnum>(
83 573918 : "line_search", line_search, "Specifies the line search type" + addtl_doc_str);
84 286959 : MooseEnum line_search_package("petsc moose", "petsc");
85 286959 : params.addParam<MooseEnum>("line_search_package",
86 : line_search_package,
87 : "The solver package to use to conduct the line-search");
88 :
89 860877 : params.addParam<unsigned>("contact_line_search_allowed_lambda_cuts",
90 573918 : 2,
91 : "The number of times lambda is allowed to be cut in half in the "
92 : "contact line search. We recommend this number be roughly bounded by 0 "
93 : "<= allowed_lambda_cuts <= 3");
94 286959 : params.addParam<Real>("contact_line_search_ltol",
95 : "The linear relative tolerance to be used while the contact state is "
96 : "changing between non-linear iterations. We recommend that this tolerance "
97 : "be looser than the standard linear tolerance");
98 :
99 286959 : params += Moose::PetscSupport::getPetscValidParams();
100 286959 : params.addParam<Real>("l_tol", 1.0e-5, "Linear Relative Tolerance");
101 286959 : params.addParam<Real>("l_abs_tol", 1.0e-50, "Linear Absolute Tolerance");
102 286959 : params.addParam<unsigned int>("l_max_its", 10000, "Max Linear Iterations");
103 286959 : params.addParam<std::vector<ConvergenceName>>(
104 : "nonlinear_convergence",
105 : "Name of the Convergence object(s) to use to assess convergence of the "
106 : "nonlinear system(s) solve. If not provided, the default Convergence "
107 : "associated with the Problem will be constructed internally.");
108 286959 : params.addParam<std::vector<ConvergenceName>>(
109 : "linear_convergence",
110 : "Name of the Convergence object(s) to use to assess convergence of the "
111 : "linear system(s) solve. If not provided, the linear solver tolerance parameters are used");
112 860877 : params.addParam<bool>(
113 : "snesmf_reuse_base",
114 573918 : true,
115 : "Specifies whether or not to reuse the base vector for matrix-free calculation");
116 860877 : params.addParam<bool>(
117 573918 : "skip_exception_check", false, "Specifies whether or not to skip exception check");
118 860877 : params.addParam<bool>("compute_initial_residual_before_preset_bcs",
119 573918 : false,
120 : "Use the residual norm computed *before* solution modifying objects like "
121 : "preset BCs are imposed in relative convergence check.");
122 286959 : params.deprecateParam(
123 : "compute_initial_residual_before_preset_bcs", "use_pre_SMO_residual", "12/31/2024");
124 860877 : params.addParam<bool>(
125 : "use_pre_SMO_residual",
126 573918 : false,
127 : "Compute the pre-SMO residual norm and use it in the relative convergence check. The "
128 : "pre-SMO residual is computed at the begining of the time step before solution-modifying "
129 : "objects are executed. Solution-modifying objects include preset BCs, constraints, "
130 : "predictors, etc.");
131 286959 : params.addParam<bool>("automatic_scaling", "Whether to use automatic scaling for the variables.");
132 286959 : params.addParam<std::vector<bool>>(
133 : "compute_scaling_once",
134 : {true},
135 : "Whether the scaling factors should only be computed once at the beginning of the simulation "
136 : "through an extra Jacobian evaluation. If this is set to false, then the scaling factors "
137 : "will be computed during an extra Jacobian evaluation at the beginning of every time step. "
138 : "Vector entries correspond to each nonlinear system.");
139 286959 : params.addParam<std::vector<bool>>(
140 : "off_diagonals_in_auto_scaling",
141 : {false},
142 : "Whether to consider off-diagonals when determining automatic scaling factors. Vector "
143 : "entries correspond to each nonlinear system.");
144 286959 : params.addRangeCheckedParam<std::vector<Real>>(
145 : "resid_vs_jac_scaling_param",
146 : {0},
147 : "0<=resid_vs_jac_scaling_param<=1",
148 : "A parameter that indicates the weighting of the residual vs the Jacobian in determining "
149 : "variable scaling parameters. A value of 1 indicates pure residual-based scaling. A value of "
150 : "0 indicates pure Jacobian-based scaling. Vector entries correspond to each nonlinear "
151 : "system.");
152 286959 : params.addParam<std::vector<std::vector<std::vector<std::string>>>>(
153 : "scaling_group_variables",
154 : "Name of variables that are grouped together for determining scale factors. (Multiple "
155 : "groups can be provided, separated by semicolon). Vector entries correspond to each "
156 : "nonlinear system.");
157 286959 : params.addParam<std::vector<std::vector<std::string>>>(
158 : "ignore_variables_for_autoscaling",
159 : "List of variables that do not participate in autoscaling. Vector entries correspond to each "
160 : "nonlinear system.");
161 860877 : params.addRangeCheckedParam<unsigned int>(
162 : "num_grids",
163 573918 : 1,
164 : "num_grids>0",
165 : "The number of grids to use for a grid sequencing algorithm. This includes the final grid, "
166 : "so num_grids = 1 indicates just one solve in a time-step");
167 286959 : params.addParam<std::vector<bool>>("residual_and_jacobian_together",
168 : {false},
169 : "Whether to compute the residual and Jacobian together. "
170 : "Vector entries correspond to each nonlinear system.");
171 :
172 860877 : params.addParam<bool>("reuse_preconditioner",
173 573918 : false,
174 : "If true reuse the previously calculated "
175 : "preconditioner for the linearized "
176 : "system across multiple solves "
177 : "spanning nonlinear iterations and time steps. "
178 : "The preconditioner resets as controlled by "
179 : "reuse_preconditioner_max_linear_its");
180 860877 : params.addParam<unsigned int>("reuse_preconditioner_max_linear_its",
181 573918 : 25,
182 : "Reuse the previously calculated "
183 : "preconditioner for the linear system "
184 : "until the number of linear iterations "
185 : "exceeds this number");
186 :
187 : // Multi-system fixed point
188 : // Defaults to false because of the difficulty of defining a good multi-system convergence
189 : // criterion, unless we add a default one to the simulation?
190 860877 : params.addParam<bool>(
191 : "multi_system_fixed_point",
192 573918 : false,
193 : "Whether to perform fixed point (Picard) iterations between the nonlinear systems.");
194 286959 : params.addParam<ConvergenceName>(
195 : "multi_system_fixed_point_convergence",
196 : "Convergence object to determine the convergence of the multi-system fixed point iteration. "
197 : "If unspecified, defaults to checking that every system is converged (based on their own "
198 : "convergence criterion)");
199 :
200 286959 : params.addParamNamesToGroup("l_tol l_abs_tol l_max_its reuse_preconditioner "
201 : "reuse_preconditioner_max_linear_its",
202 : "Linear Solver");
203 286959 : params.addParamNamesToGroup(
204 : "solve_type snesmf_reuse_base use_pre_SMO_residual "
205 : "num_grids residual_and_jacobian_together splitting nonlinear_convergence linear_convergence",
206 : "Nonlinear Solver");
207 286959 : params.addParamNamesToGroup(
208 : "automatic_scaling compute_scaling_once off_diagonals_in_auto_scaling "
209 : "scaling_group_variables resid_vs_jac_scaling_param ignore_variables_for_autoscaling",
210 : "Solver variable scaling");
211 286959 : params.addParamNamesToGroup("line_search line_search_package contact_line_search_ltol "
212 : "contact_line_search_allowed_lambda_cuts",
213 : "Solver line search");
214 286959 : params.addParamNamesToGroup("multi_system_fixed_point multi_system_fixed_point_convergence",
215 : "Multiple solver system");
216 286959 : params.addParamNamesToGroup("skip_exception_check", "Advanced");
217 :
218 573918 : return params;
219 860877 : }
220 :
221 57491 : FEProblemSolve::FEProblemSolve(Executioner & ex)
222 : : MultiSystemSolveObject(ex),
223 57491 : _num_grid_steps(cast_int<unsigned int>(getParam<unsigned int>("num_grids") - 1)),
224 57491 : _using_multi_sys_fp_iterations(getParam<bool>("multi_system_fixed_point")),
225 57491 : _multi_sys_fp_convergence(nullptr) // has not been created yet
226 : {
227 57491 : if (_moose_line_searches.find(getParam<MooseEnum>("line_search").operator std::string()) !=
228 114982 : _moose_line_searches.end())
229 0 : _problem.addLineSearch(_pars);
230 :
231 230940 : auto set_solver_params = [this, &ex](const SolverSystem & sys)
232 : {
233 57735 : const auto prefix = sys.prefix();
234 57735 : Moose::PetscSupport::storePetscOptions(_problem, prefix, ex);
235 57735 : Moose::PetscSupport::setConvergedReasonFlags(_problem, prefix);
236 :
237 : // Set solver parameter prefix and system number
238 57735 : auto & solver_params = _problem.solverParams(sys.number());
239 57735 : solver_params._prefix = prefix;
240 57735 : solver_params._solver_sys_num = sys.number();
241 57735 : };
242 :
243 : // Extract and store PETSc related settings on FEProblemBase
244 115226 : for (const auto * const sys : _systems)
245 57735 : set_solver_params(*sys);
246 :
247 : // Set linear solve parameters in the equation system
248 : // Nonlinear solve parameters are added in the DefaultNonlinearConvergence
249 57491 : EquationSystems & es = _problem.es();
250 57491 : es.parameters.set<Real>("linear solver tolerance") = getParam<Real>("l_tol");
251 57491 : es.parameters.set<Real>("linear solver absolute tolerance") = getParam<Real>("l_abs_tol");
252 57491 : es.parameters.set<unsigned int>("linear solver maximum iterations") =
253 114982 : getParam<unsigned int>("l_max_its");
254 57491 : es.parameters.set<bool>("reuse preconditioner") = getParam<bool>("reuse_preconditioner");
255 57491 : es.parameters.set<unsigned int>("reuse preconditioner maximum linear iterations") =
256 114982 : getParam<unsigned int>("reuse_preconditioner_max_linear_its");
257 :
258 : // Transfer to the Problem misc nonlinear solve optimization parameters
259 57491 : _problem.setSNESMFReuseBase(getParam<bool>("snesmf_reuse_base"),
260 114982 : _pars.isParamSetByUser("snesmf_reuse_base"));
261 57491 : _problem.skipExceptionCheck(getParam<bool>("skip_exception_check"));
262 :
263 57491 : if (isParamValid("nonlinear_convergence"))
264 : {
265 292 : if (_problem.onlyAllowDefaultNonlinearConvergence())
266 0 : mooseError("The selected problem does not allow 'nonlinear_convergence' to be set.");
267 292 : _problem.setNonlinearConvergenceNames(
268 : getParam<std::vector<ConvergenceName>>("nonlinear_convergence"));
269 : }
270 : else
271 57199 : _problem.setNeedToAddDefaultNonlinearConvergence();
272 57491 : if (isParamValid("linear_convergence"))
273 : {
274 12 : if (_problem.numLinearSystems() == 0)
275 0 : paramError(
276 : "linear_convergence",
277 : "Setting 'linear_convergence' is currently only possible for solving linear systems");
278 12 : _problem.setLinearConvergenceNames(
279 : getParam<std::vector<ConvergenceName>>("linear_convergence"));
280 : }
281 :
282 : // Check whether the user has explicitly requested automatic scaling and is using a solve type
283 : // without a matrix. If so, then we warn them
284 57932 : if ((_pars.isParamSetByUser("automatic_scaling") && getParam<bool>("automatic_scaling")) &&
285 441 : std::all_of(_systems.begin(),
286 : _systems.end(),
287 441 : [this](const auto & solver_sys)
288 441 : { return _problem.solverParams(solver_sys->number())._type == Moose::ST_JFNK; }))
289 : {
290 0 : paramWarning("automatic_scaling",
291 : "Automatic scaling isn't implemented for the case where you do not have a "
292 : "preconditioning matrix. No scaling will be applied");
293 0 : _problem.automaticScaling(false);
294 : }
295 : else
296 : // Check to see whether automatic_scaling has been specified anywhere, including at the
297 : // application level. No matter what: if we don't have a matrix, we don't do scaling
298 57491 : _problem.automaticScaling(
299 114982 : isParamValid("automatic_scaling")
300 57984 : ? getParam<bool>("automatic_scaling")
301 56998 : : (getMooseApp().defaultAutomaticScaling() &&
302 0 : std::any_of(_systems.begin(),
303 : _systems.end(),
304 0 : [this](const auto & solver_sys) {
305 0 : return _problem.solverParams(solver_sys->number())._type !=
306 0 : Moose::ST_JFNK;
307 : })));
308 :
309 57491 : if (!_using_multi_sys_fp_iterations && isParamValid("multi_system_fixed_point_convergence"))
310 4 : paramError("multi_system_fixed_point_convergence",
311 : "Cannot set a convergence object for multi-system fixed point iterations if "
312 : "'multi_system_fixed_point' is set to false");
313 57487 : if (_using_multi_sys_fp_iterations && !isParamValid("multi_system_fixed_point_convergence"))
314 4 : paramError("multi_system_fixed_point_convergence",
315 : "Must set a convergence object for multi-system fixed point iterations if using "
316 : "multi-system fixed point iterations");
317 :
318 : // Set the same parameters to every nonlinear system by default
319 57483 : int i_nl_sys = -1;
320 115182 : for (const auto i_sys : index_range(_systems))
321 : {
322 57711 : auto nl_ptr = dynamic_cast<NonlinearSystemBase *>(_systems[i_sys]);
323 : // Linear systems have very different parameters at the moment
324 57711 : if (!nl_ptr)
325 1119 : continue;
326 56592 : auto & nl = *nl_ptr;
327 56592 : i_nl_sys++;
328 :
329 56592 : nl.setPreSMOResidual(getParam<bool>("use_pre_SMO_residual"));
330 :
331 56592 : const auto & all_splittings = getParam<std::vector<std::vector<std::string>>>("splitting");
332 56592 : if (all_splittings.size())
333 0 : nl.setDecomposition(
334 0 : getParamFromNonlinearSystemVectorParam<std::vector<std::string>>("splitting", i_nl_sys));
335 : else
336 56592 : nl.setDecomposition({});
337 :
338 : const auto res_and_jac =
339 56592 : getParamFromNonlinearSystemVectorParam<bool>("residual_and_jacobian_together", i_nl_sys);
340 56584 : if (res_and_jac)
341 319 : nl.residualAndJacobianTogether();
342 :
343 : // Automatic scaling parameters
344 56584 : nl.computeScalingOnce(
345 56584 : getParamFromNonlinearSystemVectorParam<bool>("compute_scaling_once", i_nl_sys));
346 56584 : nl.autoScalingParam(
347 : getParamFromNonlinearSystemVectorParam<Real>("resid_vs_jac_scaling_param", i_nl_sys));
348 56584 : nl.offDiagonalsInAutoScaling(
349 56584 : getParamFromNonlinearSystemVectorParam<bool>("off_diagonals_in_auto_scaling", i_nl_sys));
350 56584 : if (isParamValid("scaling_group_variables"))
351 16 : nl.scalingGroupVariables(
352 32 : getParamFromNonlinearSystemVectorParam<std::vector<std::vector<std::string>>>(
353 : "scaling_group_variables", i_nl_sys));
354 56584 : if (isParamValid("ignore_variables_for_autoscaling"))
355 : {
356 : // Before setting ignore_variables_for_autoscaling, check that they are not present in
357 : // scaling_group_variables
358 14 : if (isParamValid("scaling_group_variables"))
359 : {
360 : const auto & ignore_variables_for_autoscaling =
361 : getParamFromNonlinearSystemVectorParam<std::vector<std::string>>(
362 4 : "ignore_variables_for_autoscaling", i_nl_sys);
363 : const auto & scaling_group_variables =
364 : getParamFromNonlinearSystemVectorParam<std::vector<std::vector<std::string>>>(
365 4 : "scaling_group_variables", i_nl_sys);
366 4 : for (const auto & group : scaling_group_variables)
367 8 : for (const auto & var_name : group)
368 8 : if (std::find(ignore_variables_for_autoscaling.begin(),
369 : ignore_variables_for_autoscaling.end(),
370 16 : var_name) != ignore_variables_for_autoscaling.end())
371 4 : paramError("ignore_variables_for_autoscaling",
372 : "Variables cannot be in a scaling grouping and also be ignored");
373 0 : }
374 10 : nl.ignoreVariablesForAutoscaling(
375 20 : getParamFromNonlinearSystemVectorParam<std::vector<std::string>>(
376 : "ignore_variables_for_autoscaling", i_nl_sys));
377 : }
378 : }
379 :
380 : // Multi-grid options
381 57471 : _problem.numGridSteps(_num_grid_steps);
382 57471 : }
383 :
384 : template <typename T>
385 : T
386 226378 : FEProblemSolve::getParamFromNonlinearSystemVectorParam(const std::string & param_name,
387 : unsigned int index) const
388 : {
389 226378 : const auto & param_vec = getParam<std::vector<T>>(param_name);
390 226378 : if (index > _num_nl_systems)
391 0 : paramError(param_name,
392 : "Vector parameter is requested at index (" + std::to_string(index) +
393 : ") which is larger than number of nonlinear systems (" +
394 0 : std::to_string(_num_nl_systems) + ").");
395 226378 : if (param_vec.size() == 0)
396 4 : paramError(
397 : param_name,
398 : "This parameter was passed to a routine which cannot handle empty vector parameters");
399 226374 : if (param_vec.size() != 1 && param_vec.size() != _num_nl_systems)
400 4 : paramError(param_name,
401 : "Vector parameter size (" + std::to_string(param_vec.size()) +
402 : ") is different than the number of nonlinear systems (" +
403 4 : std::to_string(_num_nl_systems) + ").");
404 :
405 : // User passed only one parameter, assume it applies to all nonlinear systems
406 226370 : if (param_vec.size() == 1)
407 226370 : return param_vec[0];
408 : else
409 0 : return param_vec[index];
410 : }
411 :
412 : bool
413 304217 : FEProblemSolve::solve()
414 : {
415 : // This should be late enough to retrieve the convergence object.
416 : // TODO: Move this to a setup phase, which does not exist for SolveObjects
417 304217 : if (isParamValid("multi_system_fixed_point_convergence"))
418 915 : _multi_sys_fp_convergence =
419 915 : &_problem.getConvergence(getParam<ConvergenceName>("multi_system_fixed_point_convergence"));
420 :
421 : // Outer loop for multi-grid convergence
422 304217 : bool converged = false;
423 304217 : unsigned int num_fp_multisys_iters = 0;
424 604527 : for (MooseIndex(_num_grid_steps) grid_step = 0; grid_step <= _num_grid_steps; ++grid_step)
425 : {
426 : // Multi-system fixed point loop
427 : // Use a convergence object if provided, if not, use a reasonable default of every nested system
428 : // being converged
429 304268 : num_fp_multisys_iters = 0;
430 304268 : converged = false;
431 608507 : while (!converged)
432 : {
433 : // Loop over each system
434 616136 : for (const auto sys : _systems)
435 : {
436 311897 : const bool is_nonlinear = (dynamic_cast<NonlinearSystemBase *>(sys) != nullptr);
437 :
438 : // Call solve on the problem for that system
439 311897 : if (is_nonlinear)
440 306989 : _problem.solve(sys->number());
441 : else
442 : {
443 : const auto linear_sys_number =
444 4908 : cast_int<unsigned int>(sys->number() - _problem.numNonlinearSystems());
445 4908 : _problem.solveLinearSystem(linear_sys_number, &_problem.getPetscOptions());
446 :
447 : // This is for postprocessing purposes in case none of the objects
448 : // request the gradients.
449 : // TODO: Somehow collect information if the postprocessors
450 : // need gradients and if nothing needs this, just skip it
451 4908 : _problem.getLinearSystem(linear_sys_number).computeGradients();
452 : }
453 :
454 : // Check convergence
455 : const auto solve_name =
456 311819 : _systems.size() == 1 ? " Solve" : "System " + sys->name() + ": Solve";
457 311819 : if (_problem.shouldSolve())
458 : {
459 283618 : if (_problem.converged(sys->number()))
460 279738 : _console << COLOR_GREEN << solve_name << " Converged!" << COLOR_DEFAULT << std::endl;
461 : else
462 : {
463 3876 : _console << COLOR_RED << solve_name << " Did NOT Converge!" << COLOR_DEFAULT
464 3876 : << std::endl;
465 3876 : return false;
466 : }
467 : }
468 : else
469 28201 : _console << COLOR_GREEN << solve_name << " Skipped!" << COLOR_DEFAULT << std::endl;
470 311815 : }
471 :
472 : // Assess convergence of the multi-system fixed point iteration
473 304239 : if (!_using_multi_sys_fp_iterations)
474 299395 : converged = true;
475 : else
476 : {
477 4844 : converged = _multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
478 : Convergence::MooseConvergenceStatus::CONVERGED;
479 4844 : if (_multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
480 : Convergence::MooseConvergenceStatus::DIVERGED)
481 0 : break;
482 : }
483 304239 : num_fp_multisys_iters++;
484 : }
485 :
486 300310 : if (grid_step != _num_grid_steps)
487 51 : _problem.uniformRefine();
488 : }
489 :
490 300259 : if (_multi_sys_fp_convergence)
491 915 : return (_multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
492 915 : Convergence::MooseConvergenceStatus::CONVERGED);
493 : else
494 299344 : return converged;
495 : }
|