https://mooseframework.inl.gov
FEProblemSolve.C
Go to the documentation of this file.
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"
18 
19 std::set<std::string> const FEProblemSolve::_moose_line_searches = {"contact", "project"};
20 
21 const std::set<std::string> &
23 {
24  return _moose_line_searches;
25 }
26 
29 {
31 
32  params.addParam<unsigned int>("nl_max_its", 50, "Max Nonlinear Iterations");
33  params.addParam<unsigned int>("nl_forced_its", 0, "The Number of Forced Nonlinear Iterations");
34  params.addParam<unsigned int>("nl_max_funcs", 10000, "Max Nonlinear solver function evaluations");
35  params.addParam<Real>("nl_abs_tol", 1.0e-50, "Nonlinear Absolute Tolerance");
36  params.addParam<Real>("nl_rel_tol", 1.0e-8, "Nonlinear Relative Tolerance");
37  params.addParam<Real>(
38  "nl_div_tol",
39  1.0e10,
40  "Nonlinear Relative Divergence Tolerance. A negative value disables this check.");
41  params.addParam<Real>(
42  "nl_abs_div_tol",
43  1.0e50,
44  "Nonlinear Absolute Divergence Tolerance. A negative value disables this check.");
45  params.addParam<Real>("nl_abs_step_tol", 0., "Nonlinear Absolute step Tolerance");
46  params.addParam<Real>("nl_rel_step_tol", 0., "Nonlinear Relative step Tolerance");
47  params.addParam<unsigned int>("n_max_nonlinear_pingpong",
48  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  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  return params;
59 }
60 
63 {
66 
67  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  std::set<std::string> line_searches = mooseLineSearches();
75 
76  std::set<std::string> alias_line_searches = {"default", "none", "basic"};
77  line_searches.insert(alias_line_searches.begin(), alias_line_searches.end());
78  std::set<std::string> petsc_line_searches = Moose::PetscSupport::getPetscValidLineSearches();
79  line_searches.insert(petsc_line_searches.begin(), petsc_line_searches.end());
80  std::string line_search_string = Moose::stringify(line_searches, " ");
81  MooseEnum line_search(line_search_string, "default");
82  std::string addtl_doc_str(" (Note: none = basic)");
83  params.addParam<MooseEnum>(
84  "line_search", line_search, "Specifies the line search type" + addtl_doc_str);
85  MooseEnum line_search_package("petsc moose", "petsc");
86  params.addParam<MooseEnum>("line_search_package",
87  line_search_package,
88  "The solver package to use to conduct the line-search");
89 
90  params.addParam<unsigned>("contact_line_search_allowed_lambda_cuts",
91  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  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 
101  params.addParam<Real>("l_tol", 1.0e-5, "Linear Relative Tolerance");
102  params.addParam<Real>("l_abs_tol", 1.0e-50, "Linear Absolute Tolerance");
103  params.addParam<unsigned int>("l_max_its", 10000, "Max Linear Iterations");
104  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  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  params.addParam<bool>(
114  "snesmf_reuse_base",
115  true,
116  "Specifies whether or not to reuse the base vector for matrix-free calculation");
117  params.addParam<bool>(
118  "skip_exception_check", false, "Specifies whether or not to skip exception check");
119  params.addParam<bool>("compute_initial_residual_before_preset_bcs",
120  false,
121  "Use the residual norm computed *before* solution modifying objects like "
122  "preset BCs are imposed in relative convergence check.");
123  params.deprecateParam(
124  "compute_initial_residual_before_preset_bcs", "use_pre_SMO_residual", "12/31/2024");
125  params.addParam<bool>(
126  "use_pre_SMO_residual",
127  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  params.addParam<bool>("automatic_scaling", "Whether to use automatic scaling for the variables.");
133  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  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  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  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  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  params.addRangeCheckedParam<unsigned int>(
163  "num_grids",
164  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  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  params.addParam<bool>("reuse_preconditioner",
174  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  params.addParam<unsigned int>("reuse_preconditioner_max_linear_its",
182  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  params.addParam<bool>(
192  "multi_system_fixed_point",
193  false,
194  "Whether to perform fixed point (Picard) iterations between the nonlinear systems.");
195  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  params.addParamNamesToGroup("l_tol l_abs_tol l_max_its reuse_preconditioner "
202  "reuse_preconditioner_max_linear_its",
203  "Linear Solver");
204  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  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  params.addParamNamesToGroup("line_search line_search_package contact_line_search_ltol "
213  "contact_line_search_allowed_lambda_cuts",
214  "Solver line search");
215  params.addParamNamesToGroup("multi_system_fixed_point multi_system_fixed_point_convergence",
216  "Multiple solver system");
217  params.addParamNamesToGroup("skip_exception_check", "Advanced");
218 
219  return params;
220 }
221 
224  _num_grid_steps(cast_int<unsigned int>(getParam<unsigned int>("num_grids") - 1)),
225  _using_multi_sys_fp_iterations(getParam<bool>("multi_system_fixed_point")),
226  _multi_sys_fp_convergence(nullptr) // has not been created yet
227 {
228  if (_moose_line_searches.find(getParam<MooseEnum>("line_search").operator std::string()) !=
229  _moose_line_searches.end())
231 
232  auto set_solver_params = [this, &ex](const SolverSystem & sys)
233  {
234  const auto prefix = sys.prefix();
237 
238  // Set solver parameter prefix and system number
239  auto & solver_params = _problem.solverParams(sys.number());
240  solver_params._prefix = prefix;
241  solver_params._solver_sys_num = sys.number();
242  };
243 
244  // Extract and store PETSc related settings on FEProblemBase
245  for (const auto * const sys : _systems)
246  set_solver_params(*sys);
247 
248  // Set linear solve parameters in the equation system
249  // Nonlinear solve parameters are added in the DefaultNonlinearConvergence
250  EquationSystems & es = _problem.es();
251  es.parameters.set<Real>("linear solver tolerance") = getParam<Real>("l_tol");
252  es.parameters.set<Real>("linear solver absolute tolerance") = getParam<Real>("l_abs_tol");
253  es.parameters.set<unsigned int>("linear solver maximum iterations") =
254  getParam<unsigned int>("l_max_its");
255  es.parameters.set<bool>("reuse preconditioner") = getParam<bool>("reuse_preconditioner");
256  es.parameters.set<unsigned int>("reuse preconditioner maximum linear iterations") =
257  getParam<unsigned int>("reuse_preconditioner_max_linear_its");
258 
259  // Transfer to the Problem misc nonlinear solve optimization parameters
260  _problem.setSNESMFReuseBase(getParam<bool>("snesmf_reuse_base"),
261  _pars.isParamSetByUser("snesmf_reuse_base"));
262  _problem.skipExceptionCheck(getParam<bool>("skip_exception_check"));
263 
264  if (isParamValid("nonlinear_convergence"))
265  {
267  mooseError("The selected problem does not allow 'nonlinear_convergence' to be set.");
269  getParam<std::vector<ConvergenceName>>("nonlinear_convergence"));
270  }
271  else
273  if (isParamValid("linear_convergence"))
274  {
275  if (_problem.numLinearSystems() == 0)
276  paramError(
277  "linear_convergence",
278  "Setting 'linear_convergence' is currently only possible for solving linear systems");
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  if ((_pars.isParamSetByUser("automatic_scaling") && getParam<bool>("automatic_scaling")) &&
286  std::all_of(_systems.begin(),
287  _systems.end(),
288  [this](const auto & solver_sys)
289  { return _problem.solverParams(solver_sys->number())._type == Moose::ST_JFNK; }))
290  {
291  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  _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
300  isParamValid("automatic_scaling")
301  ? getParam<bool>("automatic_scaling")
302  : (getMooseApp().defaultAutomaticScaling() &&
303  std::any_of(_systems.begin(),
304  _systems.end(),
305  [this](const auto & solver_sys) {
306  return _problem.solverParams(solver_sys->number())._type !=
308  })));
309 
310  if (!_using_multi_sys_fp_iterations && isParamValid("multi_system_fixed_point_convergence"))
311  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  if (_using_multi_sys_fp_iterations && !isParamValid("multi_system_fixed_point_convergence"))
315  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  int i_nl_sys = -1;
321  for (const auto i_sys : index_range(_systems))
322  {
323  auto nl_ptr = dynamic_cast<NonlinearSystemBase *>(_systems[i_sys]);
324  // Linear systems have very different parameters at the moment
325  if (!nl_ptr)
326  continue;
327  auto & nl = *nl_ptr;
328  i_nl_sys++;
329 
330  nl.setPreSMOResidual(getParam<bool>("use_pre_SMO_residual"));
331 
332  const auto & all_splittings = getParam<std::vector<std::vector<std::string>>>("splitting");
333  if (all_splittings.size())
334  nl.setDecomposition(
335  getParamFromNonlinearSystemVectorParam<std::vector<std::string>>("splitting", i_nl_sys));
336  else
337  nl.setDecomposition({});
338 
339  const auto res_and_jac =
340  getParamFromNonlinearSystemVectorParam<bool>("residual_and_jacobian_together", i_nl_sys);
341  if (res_and_jac)
342  nl.residualAndJacobianTogether();
343 
344  // Automatic scaling parameters
345  nl.computeScalingOnce(
346  getParamFromNonlinearSystemVectorParam<bool>("compute_scaling_once", i_nl_sys));
347  nl.autoScalingParam(
348  getParamFromNonlinearSystemVectorParam<Real>("resid_vs_jac_scaling_param", i_nl_sys));
349  nl.offDiagonalsInAutoScaling(
350  getParamFromNonlinearSystemVectorParam<bool>("off_diagonals_in_auto_scaling", i_nl_sys));
351  if (isParamValid("scaling_group_variables"))
352  nl.scalingGroupVariables(
353  getParamFromNonlinearSystemVectorParam<std::vector<std::vector<std::string>>>(
354  "scaling_group_variables", i_nl_sys));
355  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  if (isParamValid("scaling_group_variables"))
360  {
361  const auto & ignore_variables_for_autoscaling =
362  getParamFromNonlinearSystemVectorParam<std::vector<std::string>>(
363  "ignore_variables_for_autoscaling", i_nl_sys);
364  const auto & scaling_group_variables =
365  getParamFromNonlinearSystemVectorParam<std::vector<std::vector<std::string>>>(
366  "scaling_group_variables", i_nl_sys);
367  for (const auto & group : scaling_group_variables)
368  for (const auto & var_name : group)
369  if (std::find(ignore_variables_for_autoscaling.begin(),
370  ignore_variables_for_autoscaling.end(),
371  var_name) != ignore_variables_for_autoscaling.end())
372  paramError("ignore_variables_for_autoscaling",
373  "Variables cannot be in a scaling grouping and also be ignored");
374  }
375  nl.ignoreVariablesForAutoscaling(
376  getParamFromNonlinearSystemVectorParam<std::vector<std::string>>(
377  "ignore_variables_for_autoscaling", i_nl_sys));
378  }
379  }
380 
381  // Multi-grid options
383 }
384 
385 template <typename T>
386 T
388  unsigned int index) const
389 {
390  const auto & param_vec = getParam<std::vector<T>>(param_name);
391  if (index > _num_nl_systems)
392  paramError(param_name,
393  "Vector parameter is requested at index (" + std::to_string(index) +
394  ") which is larger than number of nonlinear systems (" +
395  std::to_string(_num_nl_systems) + ").");
396  if (param_vec.size() == 0)
397  paramError(
398  param_name,
399  "This parameter was passed to a routine which cannot handle empty vector parameters");
400  if (param_vec.size() != 1 && param_vec.size() != _num_nl_systems)
401  paramError(param_name,
402  "Vector parameter size (" + std::to_string(param_vec.size()) +
403  ") is different than the number of nonlinear systems (" +
404  std::to_string(_num_nl_systems) + ").");
405 
406  // User passed only one parameter, assume it applies to all nonlinear systems
407  if (param_vec.size() == 1)
408  return param_vec[0];
409  else
410  return param_vec[index];
411 }
412 
413 void
415 {
418 }
419 
420 void
422 {
423  // nonlinear
424  const auto conv_names = _problem.getNonlinearConvergenceNames();
425  for (const auto & conv_name : conv_names)
426  {
427  auto & conv = _problem.getConvergence(conv_name);
429  }
430 
431  // linear
432  if (isParamValid("linear_convergence"))
433  {
434  const auto conv_names = getParam<std::vector<ConvergenceName>>("linear_convergence");
435  for (const auto & conv_name : conv_names)
436  {
437  auto & conv = _problem.getConvergence(conv_name);
439  }
440  }
441 
442  // multisystem fixed point
443  if (isParamValid("multi_system_fixed_point_convergence"))
444  {
446  &_problem.getConvergence(getParam<ConvergenceName>("multi_system_fixed_point_convergence"));
449  }
450 }
451 
452 bool
454 {
455  // Outer loop for multi-grid convergence
456  bool converged = false;
457  unsigned int num_fp_multisys_iters = 0;
458  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  num_fp_multisys_iters = 0;
464  converged = false;
465  while (!converged)
466  {
467  // Loop over each system
468  for (const auto sys : _systems)
469  {
470  const bool is_nonlinear = (dynamic_cast<NonlinearSystemBase *>(sys) != nullptr);
471 
472  // Call solve on the problem for that system
473  if (is_nonlinear)
474  _problem.solve(sys->number());
475  else
476  {
477  const auto linear_sys_number =
478  cast_int<unsigned int>(sys->number() - _problem.numNonlinearSystems());
479  _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  _problem.getLinearSystem(linear_sys_number).computeGradients();
486  }
487 
488  // Check convergence
489  const auto solve_name =
490  _systems.size() == 1 ? " Solve" : "System " + sys->name() + ": Solve";
491  if (_problem.shouldSolve())
492  {
493  if (_problem.converged(sys->number()))
494  _console << COLOR_GREEN << solve_name << " Converged!" << COLOR_DEFAULT << std::endl;
495  else
496  {
497  _console << COLOR_RED << solve_name << " Did NOT Converge!" << COLOR_DEFAULT
498  << std::endl;
499  return false;
500  }
501  }
502  else
503  _console << COLOR_GREEN << solve_name << " Skipped!" << COLOR_DEFAULT << std::endl;
504  }
505 
506  // Assess convergence of the multi-system fixed point iteration
508  converged = true;
509  else
510  {
511  converged = _multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
513  if (_multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
515  break;
516  }
517  num_fp_multisys_iters++;
518  }
519 
520  if (grid_step != _num_grid_steps)
522  }
523 
525  return (_multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
527  else
528  return converged;
529 }
static InputParameters validParams()
bool shouldSolve() const
const std::vector< ConvergenceName > & getNonlinearConvergenceNames() const
Gets the nonlinear system convergence object name(s).
FEProblemBase & _problem
Reference to FEProblem.
Definition: SolveObject.h:47
Moose::PetscSupport::PetscOptions & getPetscOptions()
Retrieve a writable reference the PETSc options (used by PetscSupport)
void computeGradients()
Compute the Green-Gauss gradients.
Definition: LinearSystem.C:167
std::string _prefix
Definition: SolverParams.h:35
virtual std::size_t numNonlinearSystems() const override
std::set< std::string > getPetscValidLineSearches()
Returns the valid petsc line search options as a set of strings.
Definition: PetscSupport.C:918
virtual MooseConvergenceStatus checkConvergence(unsigned int iter)=0
Returns convergence status.
static std::set< std::string > const _moose_line_searches
Moose provided line searches.
virtual bool onlyAllowDefaultNonlinearConvergence() const
Returns true if an error will result if the user supplies &#39;nonlinear_convergence&#39;.
virtual bool solve() override
Picard solve the FEProblem.
FEProblemSolve(Executioner &ex)
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void initialSetup()
Method that should be executed once, before any solve calls.
Definition: SolveObject.h:32
const unsigned int _num_grid_steps
The number of steps to perform in a grid sequencing algorithm.
virtual void solve(const unsigned int nl_sys_num)
MooseApp & getMooseApp() const
Get the MooseApp this class is associated with.
Definition: MooseBase.h:45
InputParameters emptyInputParameters()
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Nonlinear system to be solved.
void skipExceptionCheck(bool skip_exception_check)
Set a flag that indicates if we want to skip exception and stop solve.
virtual void addLineSearch(const InputParameters &)
add a MOOSE line search
void uniformRefine()
uniformly refine the problem mesh(es).
const std::string _type
The type of this class.
Definition: MooseBase.h:87
void numGridSteps(unsigned int num_grid_steps)
Set the number of steps in a grid sequences.
Convergence * _multi_sys_fp_convergence
Convergence object to assess the convergence of the multi-system fixed point iteration.
virtual Convergence & getConvergence(const std::string &name, const THREAD_ID tid=0) const
Gets a Convergence object.
Jacobian-Free Newton Krylov.
Definition: MooseTypes.h:845
virtual libMesh::EquationSystems & es() override
virtual bool converged(const unsigned int sys_num)
Eventually we want to convert this virtual over to taking a solver system number argument.
Definition: SubProblem.h:113
static InputParameters feProblemDefaultConvergenceParams()
const T & getParam(const std::string &name) const
Retrieve a parameter for the object.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
void setLinearConvergenceNames(const std::vector< ConvergenceName > &convergence_names)
Sets the linear convergence object name(s) if there is one.
Executioners are objects that do the actual work of solving your problem.
Definition: Executioner.h:30
virtual void checkIterationType(IterationType) const
Perform checks related to the iteration type.
Definition: Convergence.h:48
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
void setNonlinearConvergenceNames(const std::vector< ConvergenceName > &convergence_names)
Sets the nonlinear convergence object name(s) if there is one.
static InputParameters validParams()
std::string stringify(const T &t)
conversion to string
Definition: Conversion.h:64
void setNeedToAddDefaultNonlinearConvergence()
Sets _need_to_add_default_nonlinear_convergence to true.
LinearSystem & getLinearSystem(unsigned int sys_num)
Get non-constant reference to a linear system.
void setConvergedReasonFlags(FEProblemBase &fe_problem, const std::string &prefix)
Set flags that will instruct the user on the reason their simulation diverged from PETSc&#39;s perspectiv...
Definition: PetscSupport.C:713
bool isParamSetByUser(const std::string &name) const
Method returns true if the parameter was set by the user.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
InputParameters getPetscValidParams()
Returns the PETSc options that are common between Executioners and Preconditioners.
Definition: PetscSupport.C:924
T getParamFromNonlinearSystemVectorParam(const std::string &param_name, unsigned int index) const
Helper routine to get the nonlinear system parameter at the right index.
T & set(const std::string &)
unsigned int _num_nl_systems
Number of nonlinear systems.
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
const InputParameters & _pars
Parameters of this object, references the InputParameters stored in the InputParametersWarehouse.
SolverParams & solverParams(unsigned int solver_sys_num=0)
Get the solver parameters.
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
const bool _using_multi_sys_fp_iterations
Whether we are using fixed point iterations for multi-system.
void convergenceSetup()
Performs setup related to Convergence objects.
virtual std::size_t numLinearSystems() const override
const ConsoleStream _console
An instance of helper class to write streams to the Console objects.
void paramWarning(const std::string &param, Args... args) const
Emits a warning prefixed with the file and line number of the given param (from the input file) along...
void automaticScaling(bool automatic_scaling) override
Automatic scaling setter.
void storePetscOptions(FEProblemBase &fe_problem, const std::string &prefix, const ParallelParamObject &param_object)
Stores the PETSc options supplied from the parameter object on the problem.
Definition: PetscSupport.C:587
void setSNESMFReuseBase(bool reuse, bool set_by_user)
If or not to reuse the base vector for matrix-free calculation.
virtual void solveLinearSystem(const unsigned int linear_sys_num, const Moose::PetscSupport::PetscOptions *po=nullptr)
Build and solve a linear system.
void ErrorVector unsigned int
auto index_range(const T &sizable)
std::vector< SolverSystem * > _systems
Vector of pointers to the systems.
virtual void initialSetup() override
Method that should be executed once, before any solve calls.
void setPreSMOResidual(bool use)
Set whether to evaluate the pre-SMO residual and use it in the subsequent relative convergence checks...
Tnew cast_int(Told oldvar)
static const std::set< std::string > & mooseLineSearches()
A solve object for use when wanting to solve multiple systems.
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)
This method takes a space delimited list of parameter names and adds them to the specified group name...