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