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  std::set<std::string> line_searches = mooseLineSearches();
68 
69  std::set<std::string> alias_line_searches = {"default", "none", "basic"};
70  line_searches.insert(alias_line_searches.begin(), alias_line_searches.end());
71  std::set<std::string> petsc_line_searches = Moose::PetscSupport::getPetscValidLineSearches();
72  line_searches.insert(petsc_line_searches.begin(), petsc_line_searches.end());
73  std::string line_search_string = Moose::stringify(line_searches, " ");
74  MooseEnum line_search(line_search_string, "default");
75  std::string addtl_doc_str(" (Note: none = basic)");
76  params.addParam<MooseEnum>(
77  "line_search", line_search, "Specifies the line search type" + addtl_doc_str);
78  MooseEnum line_search_package("petsc moose", "petsc");
79  params.addParam<MooseEnum>("line_search_package",
80  line_search_package,
81  "The solver package to use to conduct the line-search");
82 
83  params.addParam<unsigned>("contact_line_search_allowed_lambda_cuts",
84  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  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 
94  params.addParam<Real>("l_tol", 1.0e-5, "Linear Relative Tolerance");
95  params.addParam<Real>("l_abs_tol", 1.0e-50, "Linear Absolute Tolerance");
96  params.addParam<unsigned int>("l_max_its", 10000, "Max Linear Iterations");
97  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  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  params.addParam<bool>(
107  "snesmf_reuse_base",
108  true,
109  "Specifies whether or not to reuse the base vector for matrix-free calculation");
110  params.addParam<bool>(
111  "skip_exception_check", false, "Specifies whether or not to skip exception check");
112  params.addParam<bool>("compute_initial_residual_before_preset_bcs",
113  false,
114  "Use the residual norm computed *before* solution modifying objects like "
115  "preset BCs are imposed in relative convergence check.");
116  params.deprecateParam(
117  "compute_initial_residual_before_preset_bcs", "use_pre_SMO_residual", "12/31/2024");
118  params.addParam<bool>(
119  "use_pre_SMO_residual",
120  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  params.addParam<bool>("automatic_scaling", "Whether to use automatic scaling for the variables.");
126  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  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  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  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  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  params.addRangeCheckedParam<unsigned int>(
156  "num_grids",
157  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  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  params.addParam<bool>("reuse_preconditioner",
167  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  params.addParam<unsigned int>("reuse_preconditioner_max_linear_its",
175  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  params.addParam<bool>(
185  "multi_system_fixed_point",
186  false,
187  "Whether to perform fixed point (Picard) iterations between the nonlinear systems.");
188  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  params.addParamNamesToGroup("l_tol l_abs_tol l_max_its reuse_preconditioner "
195  "reuse_preconditioner_max_linear_its",
196  "Linear Solver");
197  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  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  params.addParamNamesToGroup("line_search line_search_package contact_line_search_ltol "
206  "contact_line_search_allowed_lambda_cuts",
207  "Solver line search");
208  params.addParamNamesToGroup("multi_system_fixed_point multi_system_fixed_point_convergence",
209  "Multiple solver system");
210  params.addParamNamesToGroup("skip_exception_check", "Advanced");
211 
212  return params;
213 }
214 
217  _num_grid_steps(cast_int<unsigned int>(getParam<unsigned int>("num_grids") - 1)),
218  _using_multi_sys_fp_iterations(getParam<bool>("multi_system_fixed_point")),
219  _multi_sys_fp_convergence(nullptr) // has not been created yet
220 {
221  if (_moose_line_searches.find(getParam<MooseEnum>("line_search").operator std::string()) !=
222  _moose_line_searches.end())
224 
225  auto set_solver_params = [this, &ex](const SolverSystem & sys)
226  {
227  const auto prefix = sys.prefix();
230 
231  // Set solver parameter prefix and system number
232  auto & solver_params = _problem.solverParams(sys.number());
233  solver_params._prefix = prefix;
234  solver_params._solver_sys_num = sys.number();
235  };
236 
237  // Extract and store PETSc related settings on FEProblemBase
238  for (const auto * const sys : _systems)
239  set_solver_params(*sys);
240 
241  // Set linear solve parameters in the equation system
242  // Nonlinear solve parameters are added in the DefaultNonlinearConvergence
243  EquationSystems & es = _problem.es();
244  es.parameters.set<Real>("linear solver tolerance") = getParam<Real>("l_tol");
245  es.parameters.set<Real>("linear solver absolute tolerance") = getParam<Real>("l_abs_tol");
246  es.parameters.set<unsigned int>("linear solver maximum iterations") =
247  getParam<unsigned int>("l_max_its");
248  es.parameters.set<bool>("reuse preconditioner") = getParam<bool>("reuse_preconditioner");
249  es.parameters.set<unsigned int>("reuse preconditioner maximum linear iterations") =
250  getParam<unsigned int>("reuse_preconditioner_max_linear_its");
251 
252  // Transfer to the Problem misc nonlinear solve optimization parameters
253  _problem.setSNESMFReuseBase(getParam<bool>("snesmf_reuse_base"),
254  _pars.isParamSetByUser("snesmf_reuse_base"));
255  _problem.skipExceptionCheck(getParam<bool>("skip_exception_check"));
256 
257  if (isParamValid("nonlinear_convergence"))
258  {
260  mooseError("The selected problem does not allow 'nonlinear_convergence' to be set.");
262  getParam<std::vector<ConvergenceName>>("nonlinear_convergence"));
263  }
264  else
266  if (isParamValid("linear_convergence"))
267  {
268  if (_problem.numLinearSystems() == 0)
269  paramError(
270  "linear_convergence",
271  "Setting 'linear_convergence' is currently only possible for solving linear systems");
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  if ((_pars.isParamSetByUser("automatic_scaling") && getParam<bool>("automatic_scaling")) &&
279  std::all_of(_systems.begin(),
280  _systems.end(),
281  [this](const auto & solver_sys)
282  { return _problem.solverParams(solver_sys->number())._type == Moose::ST_JFNK; }))
283  {
284  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  _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
293  isParamValid("automatic_scaling")
294  ? getParam<bool>("automatic_scaling")
295  : (getMooseApp().defaultAutomaticScaling() &&
296  std::any_of(_systems.begin(),
297  _systems.end(),
298  [this](const auto & solver_sys) {
299  return _problem.solverParams(solver_sys->number())._type !=
301  })));
302 
303  if (!_using_multi_sys_fp_iterations && isParamValid("multi_system_fixed_point_convergence"))
304  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  if (_using_multi_sys_fp_iterations && !isParamValid("multi_system_fixed_point_convergence"))
308  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  int i_nl_sys = -1;
314  for (const auto i_sys : index_range(_systems))
315  {
316  auto nl_ptr = dynamic_cast<NonlinearSystemBase *>(_systems[i_sys]);
317  // Linear systems have very different parameters at the moment
318  if (!nl_ptr)
319  continue;
320  auto & nl = *nl_ptr;
321  i_nl_sys++;
322 
323  nl.setPreSMOResidual(getParam<bool>("use_pre_SMO_residual"));
324 
325  const auto res_and_jac =
326  getParamFromNonlinearSystemVectorParam<bool>("residual_and_jacobian_together", i_nl_sys);
327  if (res_and_jac)
328  nl.residualAndJacobianTogether();
329 
330  // Automatic scaling parameters
331  nl.computeScalingOnce(
332  getParamFromNonlinearSystemVectorParam<bool>("compute_scaling_once", i_nl_sys));
333  nl.autoScalingParam(
334  getParamFromNonlinearSystemVectorParam<Real>("resid_vs_jac_scaling_param", i_nl_sys));
335  nl.offDiagonalsInAutoScaling(
336  getParamFromNonlinearSystemVectorParam<bool>("off_diagonals_in_auto_scaling", i_nl_sys));
337  if (isParamValid("scaling_group_variables"))
338  nl.scalingGroupVariables(
339  getParamFromNonlinearSystemVectorParam<std::vector<std::vector<std::string>>>(
340  "scaling_group_variables", i_nl_sys));
341  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  if (isParamValid("scaling_group_variables"))
346  {
347  const auto & ignore_variables_for_autoscaling =
348  getParamFromNonlinearSystemVectorParam<std::vector<std::string>>(
349  "ignore_variables_for_autoscaling", i_nl_sys);
350  const auto & scaling_group_variables =
351  getParamFromNonlinearSystemVectorParam<std::vector<std::vector<std::string>>>(
352  "scaling_group_variables", i_nl_sys);
353  for (const auto & group : scaling_group_variables)
354  for (const auto & var_name : group)
355  if (std::find(ignore_variables_for_autoscaling.begin(),
356  ignore_variables_for_autoscaling.end(),
357  var_name) != ignore_variables_for_autoscaling.end())
358  paramError("ignore_variables_for_autoscaling",
359  "Variables cannot be in a scaling grouping and also be ignored");
360  }
361  nl.ignoreVariablesForAutoscaling(
362  getParamFromNonlinearSystemVectorParam<std::vector<std::string>>(
363  "ignore_variables_for_autoscaling", i_nl_sys));
364  }
365  }
366 
367  // Multi-grid options
369 }
370 
371 template <typename T>
372 T
374  unsigned int index) const
375 {
376  const auto & param_vec = getParam<std::vector<T>>(param_name);
377  if (index > _num_nl_systems)
378  paramError(param_name,
379  "Vector parameter is requested at index (" + std::to_string(index) +
380  ") which is larger than number of nonlinear systems (" +
381  std::to_string(_num_nl_systems) + ").");
382  if (param_vec.size() == 0)
383  paramError(
384  param_name,
385  "This parameter was passed to a routine which cannot handle empty vector parameters");
386  if (param_vec.size() != 1 && param_vec.size() != _num_nl_systems)
387  paramError(param_name,
388  "Vector parameter size (" + std::to_string(param_vec.size()) +
389  ") is different than the number of nonlinear systems (" +
390  std::to_string(_num_nl_systems) + ").");
391 
392  // User passed only one parameter, assume it applies to all nonlinear systems
393  if (param_vec.size() == 1)
394  return param_vec[0];
395  else
396  return param_vec[index];
397 }
398 
399 void
401 {
404 }
405 
406 void
408 {
409  // nonlinear
410  const auto conv_names = _problem.getNonlinearConvergenceNames();
411  for (const auto & conv_name : conv_names)
412  {
413  auto & conv = _problem.getConvergence(conv_name);
415  }
416 
417  // linear
418  if (isParamValid("linear_convergence"))
419  {
420  const auto conv_names = getParam<std::vector<ConvergenceName>>("linear_convergence");
421  for (const auto & conv_name : conv_names)
422  {
423  auto & conv = _problem.getConvergence(conv_name);
425  }
426  }
427 
428  // multisystem fixed point
429  if (isParamValid("multi_system_fixed_point_convergence"))
430  {
432  &_problem.getConvergence(getParam<ConvergenceName>("multi_system_fixed_point_convergence"));
435  }
436 }
437 
438 bool
440 {
441  // Outer loop for multi-grid convergence
442  bool converged = false;
443  unsigned int num_fp_multisys_iters = 0;
444  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  num_fp_multisys_iters = 0;
450  converged = false;
451  while (!converged)
452  {
453  // Loop over each system
454  for (const auto sys : _systems)
455  {
456  const bool is_nonlinear = (dynamic_cast<NonlinearSystemBase *>(sys) != nullptr);
457 
458  // Call solve on the problem for that system
459  if (is_nonlinear)
460  _problem.solve(sys->number());
461  else
462  {
463  const auto linear_sys_number =
464  cast_int<unsigned int>(sys->number() - _problem.numNonlinearSystems());
465  _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  _problem.getLinearSystem(linear_sys_number).computeGradients();
472  }
473 
474  // Check convergence
475  const auto solve_name =
476  _systems.size() == 1 ? " Solve" : "System " + sys->name() + ": Solve";
477  if (_problem.shouldSolve())
478  {
479  if (_problem.converged(sys->number()))
480  _console << COLOR_GREEN << solve_name << " Converged!" << COLOR_DEFAULT << std::endl;
481  else
482  {
483  _console << COLOR_RED << solve_name << " Did NOT Converge!" << COLOR_DEFAULT
484  << std::endl;
485  return false;
486  }
487  }
488  else
489  _console << COLOR_GREEN << solve_name << " Skipped!" << COLOR_DEFAULT << std::endl;
490  }
491 
492  // Assess convergence of the multi-system fixed point iteration
494  converged = true;
495  else
496  {
497  converged = _multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
499  if (_multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
501  break;
502  }
503  num_fp_multisys_iters++;
504  }
505 
506  if (grid_step != _num_grid_steps)
508  }
509 
511  return (_multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
513  else
514  return converged;
515 }
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)
const InputParameters & _pars
The object&#39;s parameters.
Definition: MooseBase.h:362
void computeGradients()
Compute the Green-Gauss gradients.
Definition: LinearSystem.C:167
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 ...
Definition: MooseBase.h:435
std::string _prefix
Definition: SolverParams.h:35
const T & getParam(const std::string &name) const
Retrieve a parameter for the object.
Definition: MooseBase.h:384
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:879
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:83
InputParameters emptyInputParameters()
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 setConvergedReasonFlags(FEProblemBase &fe_problem, std::string prefix)
Set flags that will instruct the user on the reason their simulation diverged from PETSc&#39;s perspectiv...
Definition: PetscSupport.C:672
void uniformRefine()
uniformly refine the problem mesh(es).
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:846
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()
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 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.
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:885
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 and optionally a file path to the top-level block p...
Definition: MooseBase.h:267
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
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Definition: MooseBase.h:195
const ConsoleStream _console
An instance of helper class to write streams to the Console objects.
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:545
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...
Definition: MooseBase.h:442
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)
const std::string & _type
The type of this class.
Definition: MooseBase.h:356
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...