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>(
113  "use_pre_SMO_residual",
114  false,
115  "Compute the pre-SMO residual norm and use it in the relative convergence check. The "
116  "pre-SMO residual is computed at the begining of the time step before solution-modifying "
117  "objects are executed. Solution-modifying objects include preset BCs, constraints, "
118  "predictors, etc.");
119  params.addParam<bool>("automatic_scaling", "Whether to use automatic scaling for the variables.");
120  params.addParam<std::vector<bool>>(
121  "compute_scaling_once",
122  {true},
123  "Whether the scaling factors should only be computed once at the beginning of the simulation "
124  "through an extra Jacobian evaluation. If this is set to false, then the scaling factors "
125  "will be computed during an extra Jacobian evaluation at the beginning of every time step. "
126  "Vector entries correspond to each nonlinear system.");
127  params.addParam<std::vector<bool>>(
128  "off_diagonals_in_auto_scaling",
129  {false},
130  "Whether to consider off-diagonals when determining automatic scaling factors. Vector "
131  "entries correspond to each nonlinear system.");
132  params.addRangeCheckedParam<std::vector<Real>>(
133  "resid_vs_jac_scaling_param",
134  {0},
135  "0<=resid_vs_jac_scaling_param<=1",
136  "A parameter that indicates the weighting of the residual vs the Jacobian in determining "
137  "variable scaling parameters. A value of 1 indicates pure residual-based scaling. A value of "
138  "0 indicates pure Jacobian-based scaling. Vector entries correspond to each nonlinear "
139  "system.");
140  params.addParam<std::vector<std::vector<std::vector<std::string>>>>(
141  "scaling_group_variables",
142  "Name of variables that are grouped together for determining scale factors. (Multiple "
143  "groups can be provided, separated by semicolon). Vector entries correspond to each "
144  "nonlinear system.");
145  params.addParam<std::vector<std::vector<std::string>>>(
146  "ignore_variables_for_autoscaling",
147  "List of variables that do not participate in autoscaling. Vector entries correspond to each "
148  "nonlinear system.");
149  params.addRangeCheckedParam<unsigned int>(
150  "num_grids",
151  1,
152  "num_grids>0",
153  "The number of grids to use for a grid sequencing algorithm. This includes the final grid, "
154  "so num_grids = 1 indicates just one solve in a time-step");
155  params.addParam<std::vector<bool>>("residual_and_jacobian_together",
156  {false},
157  "Whether to compute the residual and Jacobian together. "
158  "Vector entries correspond to each nonlinear system.");
159 
160  params.addParam<bool>("reuse_preconditioner",
161  false,
162  "If true reuse the previously calculated "
163  "preconditioner for the linearized "
164  "system across multiple solves "
165  "spanning nonlinear iterations and time steps. "
166  "The preconditioner resets as controlled by "
167  "reuse_preconditioner_max_linear_its");
168  params.addParam<unsigned int>("reuse_preconditioner_max_linear_its",
169  25,
170  "Reuse the previously calculated "
171  "preconditioner for the linear system "
172  "until the number of linear iterations "
173  "exceeds this number");
174 
175  // Multi-system fixed point
176  // Defaults to false because of the difficulty of defining a good multi-system convergence
177  // criterion, unless we add a default one to the simulation?
178  params.addParam<bool>(
179  "multi_system_fixed_point",
180  false,
181  "Whether to perform fixed point (Picard) iterations between the nonlinear systems.");
182  params.addParam<ConvergenceName>(
183  "multi_system_fixed_point_convergence",
184  "Convergence object to determine the convergence of the multi-system fixed point iteration. "
185  "If unspecified, defaults to checking that every system is converged (based on their own "
186  "convergence criterion)");
187 
188  params.addParamNamesToGroup("l_tol l_abs_tol l_max_its reuse_preconditioner "
189  "reuse_preconditioner_max_linear_its",
190  "Linear Solver");
191  params.addParamNamesToGroup(
192  "solve_type snesmf_reuse_base use_pre_SMO_residual "
193  "num_grids residual_and_jacobian_together nonlinear_convergence linear_convergence",
194  "Nonlinear Solver");
195  params.addParamNamesToGroup(
196  "automatic_scaling compute_scaling_once off_diagonals_in_auto_scaling "
197  "scaling_group_variables resid_vs_jac_scaling_param ignore_variables_for_autoscaling",
198  "Solver variable scaling");
199  params.addParamNamesToGroup("line_search line_search_package contact_line_search_ltol "
200  "contact_line_search_allowed_lambda_cuts",
201  "Solver line search");
202  params.addParamNamesToGroup("multi_system_fixed_point multi_system_fixed_point_convergence",
203  "Multiple solver system");
204  params.addParamNamesToGroup("skip_exception_check", "Advanced");
205 
206  return params;
207 }
208 
211  _num_grid_steps(cast_int<unsigned int>(getParam<unsigned int>("num_grids") - 1)),
212  _using_multi_sys_fp_iterations(getParam<bool>("multi_system_fixed_point")),
213  _multi_sys_fp_convergence(nullptr) // has not been created yet
214 {
215  if (_moose_line_searches.find(getParam<MooseEnum>("line_search").operator std::string()) !=
216  _moose_line_searches.end())
218 
219  auto set_solver_params = [this, &ex](const SolverSystem & sys)
220  {
221  const auto prefix = sys.prefix();
224 
225  // Set solver parameter prefix and system number
226  auto & solver_params = _problem.solverParams(sys.number());
227  solver_params._prefix = prefix;
228  solver_params._solver_sys_num = sys.number();
229  };
230 
231  // Extract and store PETSc related settings on FEProblemBase
232  for (const auto * const sys : _systems)
233  set_solver_params(*sys);
234 
235  // Set linear solve parameters in the equation system
236  // Nonlinear solve parameters are added in the DefaultNonlinearConvergence
237  EquationSystems & es = _problem.es();
238  es.parameters.set<Real>("linear solver tolerance") = getParam<Real>("l_tol");
239  es.parameters.set<Real>("linear solver absolute tolerance") = getParam<Real>("l_abs_tol");
240  es.parameters.set<unsigned int>("linear solver maximum iterations") =
241  getParam<unsigned int>("l_max_its");
242  es.parameters.set<bool>("reuse preconditioner") = getParam<bool>("reuse_preconditioner");
243  es.parameters.set<unsigned int>("reuse preconditioner maximum linear iterations") =
244  getParam<unsigned int>("reuse_preconditioner_max_linear_its");
245 
246  // Transfer to the Problem misc nonlinear solve optimization parameters
247  _problem.setSNESMFReuseBase(getParam<bool>("snesmf_reuse_base"),
248  _pars.isParamSetByUser("snesmf_reuse_base"));
249  _problem.skipExceptionCheck(getParam<bool>("skip_exception_check"));
250 
251  if (isParamValid("nonlinear_convergence"))
252  {
254  mooseError("The selected problem does not allow 'nonlinear_convergence' to be set.");
256  getParam<std::vector<ConvergenceName>>("nonlinear_convergence"));
257  }
258  else
260  if (isParamValid("linear_convergence"))
261  {
262  if (_problem.numLinearSystems() == 0)
263  paramError(
264  "linear_convergence",
265  "Setting 'linear_convergence' is currently only possible for solving linear systems");
267  getParam<std::vector<ConvergenceName>>("linear_convergence"));
268  }
269 
270  // Check whether the user has explicitly requested automatic scaling and is using a solve type
271  // without a matrix. If so, then we warn them
272  if ((_pars.isParamSetByUser("automatic_scaling") && getParam<bool>("automatic_scaling")) &&
273  std::all_of(_systems.begin(),
274  _systems.end(),
275  [this](const auto & solver_sys)
276  { return _problem.solverParams(solver_sys->number())._type == Moose::ST_JFNK; }))
277  {
278  paramWarning("automatic_scaling",
279  "Automatic scaling isn't implemented for the case where you do not have a "
280  "preconditioning matrix. No scaling will be applied");
281  _problem.automaticScaling(false);
282  }
283  else
284  // Check to see whether automatic_scaling has been specified anywhere, including at the
285  // application level. No matter what: if we don't have a matrix, we don't do scaling
287  isParamValid("automatic_scaling")
288  ? getParam<bool>("automatic_scaling")
289  : (getMooseApp().defaultAutomaticScaling() &&
290  std::any_of(_systems.begin(),
291  _systems.end(),
292  [this](const auto & solver_sys) {
293  return _problem.solverParams(solver_sys->number())._type !=
295  })));
296 
297  if (!_using_multi_sys_fp_iterations && isParamValid("multi_system_fixed_point_convergence"))
298  paramError("multi_system_fixed_point_convergence",
299  "Cannot set a convergence object for multi-system fixed point iterations if "
300  "'multi_system_fixed_point' is set to false");
301  if (_using_multi_sys_fp_iterations && !isParamValid("multi_system_fixed_point_convergence"))
302  paramError("multi_system_fixed_point_convergence",
303  "Must set a convergence object for multi-system fixed point iterations if using "
304  "multi-system fixed point iterations");
305 
306  // Set the same parameters to every nonlinear system by default
307  int i_nl_sys = -1;
308  for (const auto i_sys : index_range(_systems))
309  {
310  auto nl_ptr = dynamic_cast<NonlinearSystemBase *>(_systems[i_sys]);
311  // Linear systems have very different parameters at the moment
312  if (!nl_ptr)
313  continue;
314  auto & nl = *nl_ptr;
315  i_nl_sys++;
316 
317  nl.setPreSMOResidual(getParam<bool>("use_pre_SMO_residual"));
318 
319  const auto res_and_jac =
320  getParamFromNonlinearSystemVectorParam<bool>("residual_and_jacobian_together", i_nl_sys);
321  if (res_and_jac)
322  nl.residualAndJacobianTogether();
323 
324  // Automatic scaling parameters
325  nl.computeScalingOnce(
326  getParamFromNonlinearSystemVectorParam<bool>("compute_scaling_once", i_nl_sys));
327  nl.autoScalingParam(
328  getParamFromNonlinearSystemVectorParam<Real>("resid_vs_jac_scaling_param", i_nl_sys));
329  nl.offDiagonalsInAutoScaling(
330  getParamFromNonlinearSystemVectorParam<bool>("off_diagonals_in_auto_scaling", i_nl_sys));
331  if (isParamValid("scaling_group_variables"))
332  nl.scalingGroupVariables(
333  getParamFromNonlinearSystemVectorParam<std::vector<std::vector<std::string>>>(
334  "scaling_group_variables", i_nl_sys));
335  if (isParamValid("ignore_variables_for_autoscaling"))
336  {
337  // Before setting ignore_variables_for_autoscaling, check that they are not present in
338  // scaling_group_variables
339  if (isParamValid("scaling_group_variables"))
340  {
341  const auto & ignore_variables_for_autoscaling =
342  getParamFromNonlinearSystemVectorParam<std::vector<std::string>>(
343  "ignore_variables_for_autoscaling", i_nl_sys);
344  const auto & scaling_group_variables =
345  getParamFromNonlinearSystemVectorParam<std::vector<std::vector<std::string>>>(
346  "scaling_group_variables", i_nl_sys);
347  for (const auto & group : scaling_group_variables)
348  for (const auto & var_name : group)
349  if (std::find(ignore_variables_for_autoscaling.begin(),
350  ignore_variables_for_autoscaling.end(),
351  var_name) != ignore_variables_for_autoscaling.end())
352  paramError("ignore_variables_for_autoscaling",
353  "Variables cannot be in a scaling grouping and also be ignored");
354  }
355  nl.ignoreVariablesForAutoscaling(
356  getParamFromNonlinearSystemVectorParam<std::vector<std::string>>(
357  "ignore_variables_for_autoscaling", i_nl_sys));
358  }
359  }
360 
361  // Multi-grid options
363 }
364 
365 template <typename T>
366 T
368  unsigned int index) const
369 {
370  const auto & param_vec = getParam<std::vector<T>>(param_name);
371  if (index > _num_nl_systems)
372  paramError(param_name,
373  "Vector parameter is requested at index (" + std::to_string(index) +
374  ") which is larger than number of nonlinear systems (" +
375  std::to_string(_num_nl_systems) + ").");
376  if (param_vec.size() == 0)
377  paramError(
378  param_name,
379  "This parameter was passed to a routine which cannot handle empty vector parameters");
380  if (param_vec.size() != 1 && param_vec.size() != _num_nl_systems)
381  paramError(param_name,
382  "Vector parameter size (" + std::to_string(param_vec.size()) +
383  ") is different than the number of nonlinear systems (" +
384  std::to_string(_num_nl_systems) + ").");
385 
386  // User passed only one parameter, assume it applies to all nonlinear systems
387  if (param_vec.size() == 1)
388  return param_vec[0];
389  else
390  return param_vec[index];
391 }
392 
393 void
395 {
398  // Keep track of the solution warnings from the setup
399  if (!_app.isRecovering())
400  {
404  }
405 }
406 
407 void
409 {
410  // nonlinear
411  const auto conv_names = _problem.getNonlinearConvergenceNames();
412  for (const auto & conv_name : conv_names)
413  {
414  auto & conv = _problem.getConvergence(conv_name);
416  }
417 
418  // linear
419  if (isParamValid("linear_convergence"))
420  {
421  const auto conv_names = getParam<std::vector<ConvergenceName>>("linear_convergence");
422  for (const auto & conv_name : conv_names)
423  {
424  auto & conv = _problem.getConvergence(conv_name);
426  }
427  }
428 
429  // multisystem fixed point
430  if (isParamValid("multi_system_fixed_point_convergence"))
431  {
433  &_problem.getConvergence(getParam<ConvergenceName>("multi_system_fixed_point_convergence"));
436  }
437 }
438 
439 bool
441 {
442  // Outer loop for multi-grid convergence
443  bool converged = false;
444  unsigned int num_fp_multisys_iters = 0;
445  for (MooseIndex(_num_grid_steps) grid_step = 0; grid_step <= _num_grid_steps; ++grid_step)
446  {
447  // Multi-system fixed point loop
448  // Use a convergence object if provided, if not, use a reasonable default of every nested system
449  // being converged
450  num_fp_multisys_iters = 0;
451  converged = false;
452  while (!converged)
453  {
454  // Loop over each system
455  for (const auto sys : _systems)
456  {
457  const bool is_nonlinear = (dynamic_cast<NonlinearSystemBase *>(sys) != nullptr);
458 
459  // Call solve on the problem for that system
460  if (is_nonlinear)
461  _problem.solve(sys->number());
462  else
463  {
464  const auto linear_sys_number =
465  cast_int<unsigned int>(sys->number() - _problem.numNonlinearSystems());
466  _problem.solveLinearSystem(linear_sys_number, &_problem.getPetscOptions());
467 
468  // This is for postprocessing purposes in case none of the objects
469  // request the gradients.
470  // TODO: Somehow collect information if the postprocessors
471  // need gradients and if nothing needs this, just skip it
472  _problem.getLinearSystem(linear_sys_number).computeGradients();
473  }
474 
475  // Check convergence
476  const auto solve_name =
477  _systems.size() == 1 ? " Solve" : "System " + sys->name() + ": Solve";
478  if (_problem.shouldSolve())
479  {
480  if (_problem.converged(sys->number()))
481  _console << COLOR_GREEN << solve_name << " Converged!" << COLOR_DEFAULT << std::endl;
482  else
483  {
484  _console << COLOR_RED << solve_name << " Did NOT Converge!" << COLOR_DEFAULT
485  << std::endl;
486  return false;
487  }
488  }
489  else
490  _console << COLOR_GREEN << solve_name << " Skipped!" << COLOR_DEFAULT << std::endl;
491  }
492 
493  // Assess convergence of the multi-system fixed point iteration
495  converged = true;
496  else
497  {
498  converged = _multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
500  if (_multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
502  break;
503  }
504  num_fp_multisys_iters++;
505  }
506 
507  if (grid_step != _num_grid_steps)
509  }
510 
512  return (_multi_sys_fp_convergence->checkConvergence(num_fp_multisys_iters) ==
514  else
515  return converged;
516 }
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)
KOKKOS_INLINE_FUNCTION const T * find(const T &target, const T *const begin, const T *const end)
Find a value in an array.
Definition: KokkosUtils.h:30
const InputParameters & _pars
The object&#39;s parameters.
Definition: MooseBase.h:366
void solutionInvalidAccumulationTimeStep(const unsigned int timestep_index)
Pass the number of solution invalid occurrences from current iteration to cumulative time iteration c...
void computeGradients()
Compute the Green-Gauss gradients.
Definition: LinearSystem.C:180
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:439
std::string _prefix
Definition: SolverParams.h:35
const T & getParam(const std::string &name) const
Retrieve a parameter for the object.
Definition: MooseBase.h:388
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:883
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:87
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 solutionInvalidAccumulation()
Pass the number of solution invalid occurrences from current iteration to cumulative counters...
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:676
void syncIteration()
Sync iteration counts to main processor.
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.
SolutionInvalidity & solutionInvalidity()
Get the SolutionInvalidity for this app.
Definition: MooseApp.h:179
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
MooseApp & _app
The MOOSE application this is associated with.
Definition: MooseBase.h:357
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:889
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:271
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:199
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 paramWarning(const std::string &param, Args... args) const
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:549
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.
bool isRecovering() const
Whether or not this is a "recover" calculation.
Definition: MooseApp.C:1874
void ErrorVector unsigned int
auto index_range(const T &sizable)
const std::string & _type
The type of this class.
Definition: MooseBase.h:360
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...