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