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 : #include "MooseUtils.h"
19 :
20 : std::set<std::string> const FEProblemSolve::_moose_line_searches = {"contact", "project"};
21 :
22 : const std::set<std::string> &
23 156854 : FEProblemSolve::mooseLineSearches()
24 : {
25 156854 : return _moose_line_searches;
26 : }
27 :
28 : InputParameters
29 287903 : FEProblemSolve::feProblemDefaultConvergenceParams()
30 : {
31 287903 : InputParameters params = emptyInputParameters();
32 :
33 1151612 : params.addParam<unsigned int>("nl_max_its", 50, "Max Nonlinear Iterations");
34 1151612 : params.addParam<unsigned int>("nl_forced_its", 0, "The Number of Forced Nonlinear Iterations");
35 1151612 : params.addParam<unsigned int>("nl_max_funcs", 10000, "Max Nonlinear solver function evaluations");
36 1151612 : params.addParam<Real>("nl_abs_tol", 1.0e-50, "Nonlinear Absolute Tolerance");
37 1151612 : params.addParam<Real>("nl_rel_tol", 1.0e-8, "Nonlinear Relative Tolerance");
38 863709 : params.addParam<Real>(
39 : "nl_div_tol",
40 575806 : 1.0e10,
41 : "Nonlinear Relative Divergence Tolerance. A negative value disables this check.");
42 863709 : params.addParam<Real>(
43 : "nl_abs_div_tol",
44 575806 : 1.0e50,
45 : "Nonlinear Absolute Divergence Tolerance. A negative value disables this check.");
46 1151612 : params.addParam<Real>("nl_abs_step_tol", 0., "Nonlinear Absolute step Tolerance");
47 1151612 : params.addParam<Real>("nl_rel_step_tol", 0., "Nonlinear Relative step Tolerance");
48 863709 : params.addParam<unsigned int>("n_max_nonlinear_pingpong",
49 575806 : 100,
50 : "The maximum number of times the nonlinear residual can ping pong "
51 : "before requesting halting the current evaluation and requesting "
52 : "timestep cut for transient simulations");
53 :
54 863709 : params.addParamNamesToGroup(
55 : "nl_max_its nl_forced_its nl_max_funcs nl_abs_tol nl_rel_tol "
56 : "nl_rel_step_tol nl_abs_step_tol nl_div_tol nl_abs_div_tol n_max_nonlinear_pingpong",
57 : "Nonlinear Solver");
58 :
59 287903 : return params;
60 0 : }
61 :
62 : InputParameters
63 156854 : FEProblemSolve::validParams()
64 : {
65 156854 : InputParameters params = MultiSystemSolveObject::validParams();
66 156854 : params += FEProblemSolve::feProblemDefaultConvergenceParams();
67 :
68 156854 : std::set<std::string> line_searches = mooseLineSearches();
69 :
70 156854 : std::set<std::string> alias_line_searches = {"default", "none", "basic"};
71 156854 : line_searches.insert(alias_line_searches.begin(), alias_line_searches.end());
72 156854 : std::set<std::string> petsc_line_searches = Moose::PetscSupport::getPetscValidLineSearches();
73 156854 : line_searches.insert(petsc_line_searches.begin(), petsc_line_searches.end());
74 627416 : std::string line_search_string = Moose::stringify(line_searches, " ");
75 313708 : MooseEnum line_search(line_search_string, "default");
76 156854 : std::string addtl_doc_str(" (Note: none = basic)");
77 156854 : params.addParam<MooseEnum>(
78 313708 : "line_search", line_search, "Specifies the line search type" + addtl_doc_str);
79 627416 : MooseEnum line_search_package("petsc moose", "petsc");
80 627416 : params.addParam<MooseEnum>("line_search_package",
81 : line_search_package,
82 : "The solver package to use to conduct the line-search");
83 :
84 470562 : params.addParam<unsigned>("contact_line_search_allowed_lambda_cuts",
85 313708 : 2,
86 : "The number of times lambda is allowed to be cut in half in the "
87 : "contact line search. We recommend this number be roughly bounded by 0 "
88 : "<= allowed_lambda_cuts <= 3");
89 470562 : params.addParam<Real>("contact_line_search_ltol",
90 : "The linear relative tolerance to be used while the contact state is "
91 : "changing between non-linear iterations. We recommend that this tolerance "
92 : "be looser than the standard linear tolerance");
93 :
94 156854 : params += Moose::PetscSupport::getPetscValidParams();
95 627416 : params.addParam<Real>("l_tol", 1.0e-5, "Linear Relative Tolerance");
96 627416 : params.addParam<Real>("l_abs_tol", 1.0e-50, "Linear Absolute Tolerance");
97 627416 : params.addParam<unsigned int>("l_max_its", 10000, "Max Linear Iterations");
98 627416 : params.addParam<std::vector<ConvergenceName>>(
99 : "nonlinear_convergence",
100 : "Name of the Convergence object(s) to use to assess convergence of the "
101 : "nonlinear system(s) solve. If not provided, the default Convergence "
102 : "associated with the Problem will be constructed internally.");
103 627416 : params.addParam<std::vector<ConvergenceName>>(
104 : "linear_convergence",
105 : "Name of the Convergence object(s) to use to assess convergence of the "
106 : "linear system(s) solve. If not provided, the linear solver tolerance parameters are used");
107 470562 : params.addParam<bool>(
108 : "snesmf_reuse_base",
109 313708 : true,
110 : "Specifies whether or not to reuse the base vector for matrix-free calculation");
111 470562 : params.addParam<bool>(
112 313708 : "skip_exception_check", false, "Specifies whether or not to skip exception check");
113 470562 : params.addParam<bool>(
114 : "use_pre_SMO_residual",
115 313708 : false,
116 : "Compute the pre-SMO residual norm and use it in the relative convergence check. The "
117 : "pre-SMO residual is computed at the begining of the time step before solution-modifying "
118 : "objects are executed. Solution-modifying objects include preset BCs, constraints, "
119 : "predictors, etc.");
120 627416 : params.addParam<bool>("automatic_scaling", "Whether to use automatic scaling for the variables.");
121 627416 : params.addParam<std::vector<bool>>(
122 : "compute_scaling_once",
123 : {true},
124 : "Whether the scaling factors should only be computed once at the beginning of the simulation "
125 : "through an extra Jacobian evaluation. If this is set to false, then the scaling factors "
126 : "will be computed during an extra Jacobian evaluation at the beginning of every time step. "
127 : "Vector entries correspond to each nonlinear system.");
128 627416 : params.addParam<std::vector<bool>>(
129 : "off_diagonals_in_auto_scaling",
130 : {false},
131 : "Whether to consider off-diagonals when determining automatic scaling factors. Vector "
132 : "entries correspond to each nonlinear system.");
133 1254832 : params.addRangeCheckedParam<std::vector<Real>>(
134 : "resid_vs_jac_scaling_param",
135 : {0},
136 : "0<=resid_vs_jac_scaling_param<=1",
137 : "A parameter that indicates the weighting of the residual vs the Jacobian in determining "
138 : "variable scaling parameters. A value of 1 indicates pure residual-based scaling. A value of "
139 : "0 indicates pure Jacobian-based scaling. Vector entries correspond to each nonlinear "
140 : "system.");
141 627416 : params.addParam<std::vector<std::vector<std::vector<std::string>>>>(
142 : "scaling_group_variables",
143 : "Name of variables that are grouped together for determining scale factors. (Multiple "
144 : "groups can be provided, separated by semicolon). Vector entries correspond to each "
145 : "nonlinear system.");
146 627416 : params.addParam<std::vector<std::vector<std::string>>>(
147 : "ignore_variables_for_autoscaling",
148 : "List of variables that do not participate in autoscaling. Vector entries correspond to each "
149 : "nonlinear system.");
150 784270 : params.addRangeCheckedParam<unsigned int>(
151 : "num_grids",
152 313708 : 1,
153 : "num_grids>0",
154 : "The number of grids to use for a grid sequencing algorithm. This includes the final grid, "
155 : "so num_grids = 1 indicates just one solve in a time-step");
156 627416 : params.addParam<std::vector<bool>>("residual_and_jacobian_together",
157 : {false},
158 : "Whether to compute the residual and Jacobian together. "
159 : "Vector entries correspond to each nonlinear system.");
160 :
161 470562 : params.addParam<bool>("reuse_preconditioner",
162 313708 : false,
163 : "If true reuse the previously calculated "
164 : "preconditioner for the linearized "
165 : "system across multiple solves "
166 : "spanning nonlinear iterations and time steps. "
167 : "The preconditioner resets as controlled by "
168 : "reuse_preconditioner_max_linear_its");
169 470562 : params.addParam<unsigned int>("reuse_preconditioner_max_linear_its",
170 313708 : 25,
171 : "Reuse the previously calculated "
172 : "preconditioner for the linear system "
173 : "until the number of linear iterations "
174 : "exceeds this number");
175 :
176 : // Multi-system fixed point
177 : // Defaults to false because of the difficulty of defining a good multi-system convergence
178 : // criterion, unless we add a default one to the simulation?
179 470562 : params.addParam<bool>(
180 : "multi_system_fixed_point",
181 313708 : false,
182 : "Whether to perform fixed point (Picard) iterations between the nonlinear systems.");
183 1254832 : params.addRangeCheckedParam<std::vector<Real>>(
184 : "multi_system_fixed_point_relaxation_factor",
185 : {1.0},
186 : "multi_system_fixed_point_relaxation_factor>0 & multi_system_fixed_point_relaxation_factor<2",
187 : "Relaxation factor(s) applied to system solution updates during multi-system fixed point "
188 : "iterations; 1 disables relaxation. If one value is provided it is applied to every system; "
189 : "otherwise the vector must match the number/order of systems being solved.");
190 627416 : params.addParam<ConvergenceName>(
191 : "multi_system_fixed_point_convergence",
192 : "Convergence object to determine the convergence of the multi-system fixed point iteration. "
193 : "If unspecified, defaults to checking that every system is converged (based on their own "
194 : "convergence criterion)");
195 :
196 627416 : params.addParamNamesToGroup("l_tol l_abs_tol l_max_its reuse_preconditioner "
197 : "reuse_preconditioner_max_linear_its",
198 : "Linear Solver");
199 627416 : params.addParamNamesToGroup(
200 : "solve_type snesmf_reuse_base use_pre_SMO_residual "
201 : "num_grids residual_and_jacobian_together nonlinear_convergence linear_convergence",
202 : "Nonlinear Solver");
203 627416 : params.addParamNamesToGroup(
204 : "automatic_scaling compute_scaling_once off_diagonals_in_auto_scaling "
205 : "scaling_group_variables resid_vs_jac_scaling_param ignore_variables_for_autoscaling",
206 : "Solver variable scaling");
207 627416 : params.addParamNamesToGroup("line_search line_search_package contact_line_search_ltol "
208 : "contact_line_search_allowed_lambda_cuts",
209 : "Solver line search");
210 627416 : params.addParamNamesToGroup("multi_system_fixed_point multi_system_fixed_point_convergence "
211 : "multi_system_fixed_point_relaxation_factor",
212 : "Multiple solver system");
213 470562 : params.addParamNamesToGroup("skip_exception_check", "Advanced");
214 :
215 313708 : return params;
216 156854 : }
217 :
218 59709 : FEProblemSolve::FEProblemSolve(Executioner & ex)
219 : : MultiSystemSolveObject(ex),
220 59709 : _num_grid_steps(cast_int<unsigned int>(getParam<unsigned int>("num_grids") - 1)),
221 119418 : _using_multi_sys_fp_iterations(getParam<bool>("multi_system_fixed_point")),
222 119418 : _multi_sys_fp_convergence(nullptr) // has not been created yet
223 : {
224 179412 : if (_pars.isParamSetByUser("multi_system_fixed_point_relaxation_factor") &&
225 285 : !_using_multi_sys_fp_iterations)
226 0 : paramError("Can't use relaxation factors because multisystem fixed point iteration hasn't been "
227 : "enabled!");
228 :
229 59709 : setupMultiSystemFixedPointRelaxationFactors();
230 :
231 119418 : if (_moose_line_searches.find(getParam<MooseEnum>("line_search").operator std::string()) !=
232 119418 : _moose_line_searches.end())
233 0 : _problem.addLineSearch(_pars);
234 :
235 59977 : auto set_solver_params = [this, &ex](const SolverSystem & sys)
236 : {
237 59977 : const auto prefix = sys.prefix();
238 59977 : if (dynamic_cast<const LinearSystem *>(&sys))
239 1032 : Moose::PetscSupport::dontAddCommonSNESOptions(_problem, prefix);
240 59977 : Moose::PetscSupport::storePetscOptions(_problem, prefix, ex);
241 59977 : Moose::PetscSupport::setConvergedReasonFlags(_problem, prefix);
242 :
243 : // Set solver parameter prefix and system number
244 59977 : auto & solver_params = _problem.solverParams(sys.number());
245 59977 : solver_params._prefix = prefix;
246 59977 : solver_params._solver_sys_num = sys.number();
247 59977 : };
248 :
249 : // Extract and store PETSc related settings on FEProblemBase
250 119686 : for (const auto * const sys : _systems)
251 59977 : set_solver_params(*sys);
252 :
253 : // Set linear solve parameters in the equation system
254 : // Nonlinear solve parameters are added in the DefaultNonlinearConvergence
255 59709 : EquationSystems & es = _problem.es();
256 238836 : es.parameters.set<Real>("linear solver tolerance") = getParam<Real>("l_tol");
257 238836 : es.parameters.set<Real>("linear solver absolute tolerance") = getParam<Real>("l_abs_tol");
258 59709 : es.parameters.set<unsigned int>("linear solver maximum iterations") =
259 179127 : getParam<unsigned int>("l_max_its");
260 238836 : es.parameters.set<bool>("reuse preconditioner") = getParam<bool>("reuse_preconditioner");
261 59709 : es.parameters.set<unsigned int>("reuse preconditioner maximum linear iterations") =
262 179127 : getParam<unsigned int>("reuse_preconditioner_max_linear_its");
263 :
264 : // Transfer to the Problem misc nonlinear solve optimization parameters
265 59709 : _problem.setSNESMFReuseBase(getParam<bool>("snesmf_reuse_base"),
266 179127 : _pars.isParamSetByUser("snesmf_reuse_base"));
267 119418 : _problem.skipExceptionCheck(getParam<bool>("skip_exception_check"));
268 :
269 179127 : if (isParamValid("nonlinear_convergence"))
270 : {
271 306 : if (_problem.onlyAllowDefaultNonlinearConvergence())
272 0 : mooseError("The selected problem does not allow 'nonlinear_convergence' to be set.");
273 918 : _problem.setNonlinearConvergenceNames(
274 : getParam<std::vector<ConvergenceName>>("nonlinear_convergence"));
275 : }
276 : else
277 59403 : _problem.setNeedToAddDefaultNonlinearConvergence();
278 179127 : if (isParamValid("linear_convergence"))
279 : {
280 134 : if (_problem.numLinearSystems() == 0)
281 0 : paramError(
282 : "linear_convergence",
283 : "Setting 'linear_convergence' is currently only possible for solving linear systems");
284 402 : _problem.setLinearConvergenceNames(
285 : getParam<std::vector<ConvergenceName>>("linear_convergence"));
286 : }
287 :
288 : // Check whether the user has explicitly requested automatic scaling and is using a solve type
289 : // without a matrix. If so, then we warn them
290 180432 : if ((_pars.isParamSetByUser("automatic_scaling") && getParam<bool>("automatic_scaling")) &&
291 403 : std::all_of(_systems.begin(),
292 : _systems.end(),
293 403 : [this](const auto & solver_sys)
294 403 : { return _problem.solverParams(solver_sys->number())._type == Moose::ST_JFNK; }))
295 : {
296 0 : paramWarning("automatic_scaling",
297 : "Automatic scaling isn't implemented for the case where you do not have a "
298 : "preconditioning matrix. No scaling will be applied");
299 0 : _problem.automaticScaling(false);
300 : }
301 : else
302 : // Check to see whether automatic_scaling has been specified anywhere, including at the
303 : // application level. No matter what: if we don't have a matrix, we don't do scaling
304 59709 : _problem.automaticScaling(
305 179127 : isParamValid("automatic_scaling")
306 61062 : ? getParam<bool>("automatic_scaling")
307 59258 : : (getMooseApp().defaultAutomaticScaling() &&
308 0 : std::any_of(_systems.begin(),
309 : _systems.end(),
310 0 : [this](const auto & solver_sys)
311 : {
312 0 : return _problem.solverParams(solver_sys->number())._type !=
313 0 : Moose::ST_JFNK;
314 : })));
315 :
316 177563 : if (!_using_multi_sys_fp_iterations && isParamValid("multi_system_fixed_point_convergence"))
317 6 : paramError("multi_system_fixed_point_convergence",
318 : "Cannot set a convergence object for multi-system fixed point iterations if "
319 : "'multi_system_fixed_point' is set to false");
320 61270 : if (_using_multi_sys_fp_iterations && !isParamValid("multi_system_fixed_point_convergence"))
321 6 : paramError("multi_system_fixed_point_convergence",
322 : "Must set a convergence object for multi-system fixed point iterations if using "
323 : "multi-system fixed point iterations");
324 :
325 : // Set the same parameters to every nonlinear system by default
326 59703 : int i_nl_sys = -1;
327 119653 : for (const auto i_sys : index_range(_systems))
328 : {
329 59959 : auto nl_ptr = dynamic_cast<NonlinearSystemBase *>(_systems[i_sys]);
330 : // Linear systems have very different parameters at the moment
331 59959 : if (!nl_ptr)
332 1032 : continue;
333 58927 : auto & nl = *nl_ptr;
334 58927 : i_nl_sys++;
335 :
336 117854 : nl.setPreSMOResidual(getParam<bool>("use_pre_SMO_residual"));
337 :
338 : const auto res_and_jac =
339 117854 : getParamFromNonlinearSystemVectorParam<bool>("residual_and_jacobian_together", i_nl_sys);
340 58921 : if (res_and_jac)
341 460 : nl.residualAndJacobianTogether();
342 :
343 : // Automatic scaling parameters
344 58921 : nl.computeScalingOnce(
345 117842 : getParamFromNonlinearSystemVectorParam<bool>("compute_scaling_once", i_nl_sys));
346 117842 : nl.autoScalingParam(
347 : getParamFromNonlinearSystemVectorParam<Real>("resid_vs_jac_scaling_param", i_nl_sys));
348 58921 : nl.offDiagonalsInAutoScaling(
349 117842 : getParamFromNonlinearSystemVectorParam<bool>("off_diagonals_in_auto_scaling", i_nl_sys));
350 176763 : if (isParamValid("scaling_group_variables"))
351 15 : nl.scalingGroupVariables(
352 60 : getParamFromNonlinearSystemVectorParam<std::vector<std::vector<std::string>>>(
353 : "scaling_group_variables", i_nl_sys));
354 176763 : 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 36 : if (isParamValid("scaling_group_variables"))
359 : {
360 : const auto & ignore_variables_for_autoscaling =
361 : getParamFromNonlinearSystemVectorParam<std::vector<std::string>>(
362 6 : "ignore_variables_for_autoscaling", i_nl_sys);
363 : const auto & scaling_group_variables =
364 : getParamFromNonlinearSystemVectorParam<std::vector<std::vector<std::string>>>(
365 6 : "scaling_group_variables", i_nl_sys);
366 3 : for (const auto & group : scaling_group_variables)
367 6 : for (const auto & var_name : group)
368 6 : if (std::find(ignore_variables_for_autoscaling.begin(),
369 : ignore_variables_for_autoscaling.end(),
370 12 : var_name) != ignore_variables_for_autoscaling.end())
371 6 : paramError("ignore_variables_for_autoscaling",
372 : "Variables cannot be in a scaling grouping and also be ignored");
373 0 : }
374 9 : nl.ignoreVariablesForAutoscaling(
375 36 : getParamFromNonlinearSystemVectorParam<std::vector<std::string>>(
376 : "ignore_variables_for_autoscaling", i_nl_sys));
377 : }
378 : }
379 :
380 : // Multi-grid options
381 59694 : _problem.numGridSteps(_num_grid_steps);
382 59694 : }
383 :
384 : template <typename T>
385 : T
386 235720 : FEProblemSolve::getParamFromNonlinearSystemVectorParam(const std::string & param_name,
387 : unsigned int index) const
388 : {
389 235720 : const auto & param_vec = getParam<std::vector<T>>(param_name);
390 235720 : 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 235720 : if (param_vec.size() == 0)
396 3 : paramError(
397 : param_name,
398 : "This parameter was passed to a routine which cannot handle empty vector parameters");
399 235717 : if (param_vec.size() != 1 && param_vec.size() != _num_nl_systems)
400 3 : paramError(param_name,
401 : "Vector parameter size (" + std::to_string(param_vec.size()) +
402 : ") is different than the number of nonlinear systems (" +
403 3 : std::to_string(_num_nl_systems) + ").");
404 :
405 : // User passed only one parameter, assume it applies to all nonlinear systems
406 235714 : if (param_vec.size() == 1)
407 235714 : return param_vec[0];
408 : else
409 0 : return param_vec[index];
410 : }
411 :
412 : void
413 56801 : FEProblemSolve::initialSetup()
414 : {
415 56801 : MultiSystemSolveObject::initialSetup();
416 56801 : convergenceSetup();
417 : // Keep track of the solution warnings from the setup
418 : // before a count reset at the beginning of the time step
419 56795 : if (!_app.isRecovering())
420 : {
421 52579 : _app.solutionInvalidity().syncIteration();
422 52579 : _app.solutionInvalidity().accumulateIterationIntoTimeStepOccurences();
423 52579 : _app.solutionInvalidity().accumulateTimeStepIntoTotalOccurences(0);
424 : }
425 56795 : }
426 :
427 : void
428 56801 : FEProblemSolve::convergenceSetup()
429 : {
430 : // nonlinear
431 56801 : const auto conv_names = _problem.getNonlinearConvergenceNames();
432 112702 : for (const auto & conv_name : conv_names)
433 : {
434 55907 : auto & conv = _problem.getConvergence(conv_name);
435 55907 : conv.checkIterationType(ConvergenceIterationTypes::NONLINEAR);
436 : }
437 :
438 : // linear
439 170385 : if (isParamValid("linear_convergence"))
440 : {
441 268 : const auto conv_names = getParam<std::vector<ConvergenceName>>("linear_convergence");
442 268 : for (const auto & conv_name : conv_names)
443 : {
444 134 : auto & conv = _problem.getConvergence(conv_name);
445 134 : conv.checkIterationType(ConvergenceIterationTypes::LINEAR);
446 : }
447 134 : }
448 :
449 : // multisystem fixed point
450 170385 : if (isParamValid("multi_system_fixed_point_convergence"))
451 : {
452 773 : _multi_sys_fp_convergence =
453 1546 : &_problem.getConvergence(getParam<ConvergenceName>("multi_system_fixed_point_convergence"));
454 773 : _multi_sys_fp_convergence->checkIterationType(
455 : ConvergenceIterationTypes::MULTISYSTEM_FIXED_POINT);
456 : }
457 56795 : }
458 :
459 : void
460 59709 : FEProblemSolve::setupMultiSystemFixedPointRelaxationFactors()
461 : {
462 : _multi_sys_fp_relax_factors =
463 119418 : getParam<std::vector<Real>>("multi_system_fixed_point_relaxation_factor");
464 59709 : if (_multi_sys_fp_relax_factors.size() == 1)
465 59709 : _multi_sys_fp_relax_factors.resize(_systems.size(), _multi_sys_fp_relax_factors[0]);
466 0 : else if (_multi_sys_fp_relax_factors.size() != _systems.size())
467 0 : paramError("multi_system_fixed_point_relaxation_factor",
468 0 : "Must provide either 1 value or " + Moose::stringify(_systems.size()) +
469 : " values (one per system in the solve order).");
470 59709 : }
471 :
472 : bool
473 308802 : FEProblemSolve::solve()
474 : {
475 : // Outer loop for multi-grid convergence
476 308802 : bool converged = false;
477 308802 : unsigned int num_fp_multisys_iters = 0;
478 :
479 615736 : for (MooseIndex(_num_grid_steps) grid_step = 0; grid_step <= _num_grid_steps; ++grid_step)
480 : {
481 : // Multi-system fixed point loop
482 : // Use a convergence object if provided, if not, use a reasonable default of every nested system
483 : // being converged
484 308852 : num_fp_multisys_iters = 0;
485 308852 : converged = false;
486 640946 : while (!converged)
487 : {
488 : // Loop over each system
489 670083 : for (const auto sys_i : index_range(_systems))
490 : {
491 337989 : auto * const sys = _systems[sys_i];
492 337989 : const bool is_nonlinear = (dynamic_cast<NonlinearSystemBase *>(sys) != nullptr);
493 : const Real fp_relax =
494 337989 : _using_multi_sys_fp_iterations ? _multi_sys_fp_relax_factors[sys_i] : 1.0;
495 : const bool apply_fp_relax =
496 337989 : _using_multi_sys_fp_iterations && !MooseUtils::absoluteFuzzyEqual(fp_relax, 1.0);
497 337989 : if (apply_fp_relax)
498 : {
499 21100 : sys->setFixedPointRelaxationFactor(fp_relax);
500 21100 : sys->saveOldSolutionForFixedPointRelaxation();
501 : }
502 :
503 : // Call solve on the problem for that system
504 337989 : if (is_nonlinear)
505 311903 : _problem.solve(sys->number());
506 : else
507 : {
508 : const auto linear_sys_number =
509 26086 : cast_int<unsigned int>(sys->number() - _problem.numNonlinearSystems());
510 26086 : _problem.solveLinearSystem(linear_sys_number, &_problem.getPetscOptions());
511 : }
512 :
513 : // Check convergence
514 : const auto solve_name =
515 667908 : _systems.size() == 1 ? " Solve" : "System " + sys->name() + ": Solve";
516 337931 : if (_problem.shouldSolve())
517 : {
518 307354 : if (_problem.converged(sys->number()))
519 : {
520 305494 : if (apply_fp_relax)
521 21100 : sys->applyFixedPointRelaxation();
522 305494 : _console << COLOR_GREEN << solve_name << " Converged!" << COLOR_DEFAULT << std::endl;
523 : }
524 : else
525 : {
526 1857 : _console << COLOR_RED << solve_name << " Did NOT Converge!" << COLOR_DEFAULT
527 1857 : << std::endl;
528 1857 : if (apply_fp_relax)
529 0 : sys->clearFixedPointRelaxation();
530 1857 : return false;
531 : }
532 : }
533 : else
534 30577 : _console << COLOR_GREEN << solve_name << " Skipped!" << COLOR_DEFAULT << std::endl;
535 :
536 336071 : if (!is_nonlinear)
537 : {
538 : const auto linear_sys_number =
539 26077 : cast_int<unsigned int>(sys->number() - _problem.numNonlinearSystems());
540 26077 : auto & linear_sys = _problem.getLinearSystem(linear_sys_number);
541 :
542 : // This is for postprocessing purposes in case none of the objects request the gradients.
543 : // TODO: Somehow collect information if the postprocessors need gradients and if nothing
544 : // needs this, just skip it
545 26077 : linear_sys.computeGradients();
546 : }
547 :
548 336071 : if (apply_fp_relax)
549 21100 : sys->clearFixedPointRelaxation();
550 337928 : }
551 :
552 : // Assess convergence of the multi-system fixed point iteration
553 332094 : if (!_using_multi_sys_fp_iterations)
554 306171 : converged = true;
555 : else
556 : {
557 25923 : converged = _multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
558 : Convergence::MooseConvergenceStatus::CONVERGED;
559 25923 : if (_multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
560 : Convergence::MooseConvergenceStatus::DIVERGED)
561 0 : break;
562 : }
563 332094 : num_fp_multisys_iters++;
564 : }
565 :
566 306934 : if (grid_step != _num_grid_steps)
567 50 : _problem.uniformRefine();
568 : }
569 :
570 306884 : if (_multi_sys_fp_convergence)
571 763 : return (_multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
572 763 : Convergence::MooseConvergenceStatus::CONVERGED);
573 : else
574 306121 : return converged;
575 : }
|