https://mooseframework.inl.gov
SIMPLESolveNonlinearAssembly.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 
11 #include "FEProblem.h"
12 #include "SegregatedSolverUtils.h"
13 #include "NonlinearSystemBase.h"
14 
15 #include "libmesh/nonlinear_implicit_system.h"
16 
17 using namespace libMesh;
18 
21 {
23 
24  params.addParam<TagName>("pressure_gradient_tag",
25  "pressure_momentum_kernels",
26  "The name of the tags associated with the kernels in the momentum "
27  "equations which are not related to the pressure gradient.");
28 
29  /*
30  * The names of the different systems in the segregated solver
31  */
32  params.addParam<std::vector<SolverSystemName>>(
33  "turbulence_systems", {}, "The solver system(s) for the turbulence equation(s).");
34 
35  /*
36  * Relaxation parameters for the different system
37  */
38  params.addParam<std::vector<Real>>(
39  "turbulence_equation_relaxation",
40  std::vector<Real>(),
41  "The relaxation which should be used for the turbulence equations "
42  "equations. (=1 for no relaxation, "
43  "diagonal dominance will still be enforced)");
44 
45  params.addParam<std::vector<Real>>(
46  "turbulence_field_min_limit",
47  std::vector<Real>(),
48  "The lower limit imposed on turbulent quantities. The recommended value for robustness "
49  "is 1e-8.");
50 
51  params.addParamNamesToGroup("turbulence_equation_relaxation", "Relaxation");
52 
53  /*
54  * Petsc options for every equations in the system
55  */
56  params.addParam<MultiMooseEnum>("turbulence_petsc_options",
58  "Singleton PETSc options for the turbulence equation(s)");
59  params.addParam<MultiMooseEnum>("turbulence_petsc_options_iname",
61  "Names of PETSc name/value pairs for the turbulence equation(s)");
62  params.addParam<std::vector<std::string>>(
63  "turbulence_petsc_options_value",
64  "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\" for the "
65  "turbulence equation");
66 
67  params.addParamNamesToGroup(
68  "turbulence_petsc_options turbulence_petsc_options_iname turbulence_petsc_options_value",
69  "PETSc Control");
70 
71  /*
72  * Iteration tolerances for the different equations
73  */
74 
75  params.addParam<std::vector<Real>>(
76  "turbulence_absolute_tolerance",
77  std::vector<Real>(),
78  "The absolute tolerance(s) on the normalized residual(s) of the turbulence equation(s).");
79 
80  params.addParamNamesToGroup("solid_energy_absolute_tolerance turbulence_absolute_tolerance",
81  "Iteration Control");
82  /*
83  * Linear iteration tolerances for the different equations
84  */
85  params.addRangeCheckedParam<Real>("turbulence_l_tol",
86  1e-5,
87  "0.0<=turbulence_l_tol & turbulence_l_tol<1.0",
88  "The relative tolerance on the normalized residual in the "
89  "linear solver of the turbulence equation(s).");
90  params.addRangeCheckedParam<Real>("turbulence_l_abs_tol",
91  1e-10,
92  "0.0<turbulence_l_abs_tol",
93  "The absolute tolerance on the normalized residual in the "
94  "linear solver of the turbulence equation(s).");
95  params.addParam<unsigned int>(
96  "turbulence_l_max_its",
97  10000,
98  "The maximum allowed iterations in the linear solver of the turbulence equation(s).");
99 
100  params.addParamNamesToGroup("turbulence_l_tol "
101  "turbulence_l_abs_tol turbulence_l_max_its",
102  "Linear Iteration Control");
103 
104  return params;
105 }
106 
108  : SIMPLESolveBase(ex),
109  _pressure_sys_number(_problem.nlSysNum(getParam<SolverSystemName>("pressure_system"))),
110  _pressure_system(_problem.getNonlinearSystemBase(_pressure_sys_number)),
111  _has_turbulence_systems(!getParam<std::vector<SolverSystemName>>("turbulence_systems").empty()),
112  _energy_sys_number(_has_energy_system
113  ? _problem.nlSysNum(getParam<SolverSystemName>("energy_system"))
114  : libMesh::invalid_uint),
115  _energy_system(_has_energy_system ? &_problem.getNonlinearSystemBase(_energy_sys_number)
116  : nullptr),
117  _solid_energy_sys_number(
118  _has_solid_energy_system
119  ? _problem.nlSysNum(getParam<SolverSystemName>("solid_energy_system"))
120  : libMesh::invalid_uint),
121  _solid_energy_system(_has_solid_energy_system
122  ? &_problem.getNonlinearSystemBase(_solid_energy_sys_number)
123  : nullptr),
124  _turbulence_system_names(getParam<std::vector<SolverSystemName>>("turbulence_systems")),
125  _turbulence_equation_relaxation(getParam<std::vector<Real>>("turbulence_equation_relaxation")),
126  _turbulence_field_min_limit(getParam<std::vector<Real>>("turbulence_field_min_limit")),
127  _turbulence_l_abs_tol(getParam<Real>("turbulence_l_abs_tol")),
128  _turbulence_absolute_tolerance(getParam<std::vector<Real>>("turbulence_absolute_tolerance")),
129  _pressure_tag_name(getParam<TagName>("pressure_gradient_tag")),
130  _pressure_tag_id(_problem.addVectorTag(_pressure_tag_name))
131 {
132  // We disable this considering that this object passes petsc options a little differently
134 
135  // We fetch the system numbers for the momentum components plus add vectors
136  // for removing the contribution from the pressure gradient terms.
137  for (auto system_i : index_range(_momentum_system_names))
138  {
140  _momentum_systems.push_back(
142  _momentum_systems[system_i]->addVector(_pressure_tag_id, false, ParallelType::PARALLEL);
143 
144  // We disable this considering that this object passes petsc options a little differently
145  _momentum_systems[system_i]->system().prefix_with_name(false);
146  }
147 
149  for (auto system_i : index_range(_passive_scalar_system_names))
150  {
153  _passive_scalar_systems.push_back(
155 
156  // We disable this considering that this object passes petsc options a little differently
157  _passive_scalar_systems[system_i]->system().prefix_with_name(false);
158  }
159 
161  {
162  for (auto system_i : index_range(_turbulence_system_names))
163  {
165  _turbulence_systems.push_back(
167 
168  // We disable this considering that this object passes petsc options a little differently
169  _turbulence_systems[system_i]->system().prefix_with_name(false);
170  }
171 
172  // We check for input errors with regards to the turbulence equations. At the same time, we
173  // set up the corresponding system numbers
175  paramError("turbulence_equation_relaxation",
176  "The number of equation relaxation parameters does not match the number of "
177  "turbulence scalar equations!");
179  paramError("turbulence_absolute_tolerance",
180  "The number of absolute tolerances does not match the number of "
181  "turbulence equations!");
182  if (_turbulence_field_min_limit.empty())
183  // If no minimum bounds are given, initialize to default value 1e-8
185  else if (_turbulence_system_names.size() != _turbulence_field_min_limit.size())
186  paramError("turbulence_field_min_limit",
187  "The number of lower bounds for turbulent quantities does not match the "
188  "number of turbulence equations!");
189  }
190 
191  if (isParamValid("solid_energy_system") && !_has_energy_system)
192  paramError(
193  "solid_energy_system",
194  "We cannot solve a solid energy system without solving for the fluid energy as well!");
195 
197  {
198  const auto & turbulence_petsc_options = getParam<MultiMooseEnum>("turbulence_petsc_options");
199  const auto & turbulence_petsc_pair_options = getParam<MooseEnumItem, std::string>(
200  "turbulence_petsc_options_iname", "turbulence_petsc_options_value");
202  turbulence_petsc_options, "-", *this, _turbulence_petsc_options);
203  Moose::PetscSupport::addPetscPairsToPetscOptions(turbulence_petsc_pair_options,
204  _problem.mesh().dimension(),
205  "-",
206  *this,
208 
209  _turbulence_linear_control.real_valued_data["rel_tol"] = getParam<Real>("turbulence_l_tol");
210  _turbulence_linear_control.real_valued_data["abs_tol"] = getParam<Real>("turbulence_l_abs_tol");
212  getParam<unsigned int>("turbulence_l_max_its");
213  }
214  else
215  checkDependentParameterError("turbulence_system",
216  {"turbulence_petsc_options",
217  "turbulence_petsc_options_iname",
218  "turbulence_petsc_options_value",
219  "turbulence_l_tol",
220  "turbulence_l_abs_tol",
221  "turbulence_l_max_its",
222  "turbulence_equation_relaxation",
223  "turbulence_absolute_tolerance"},
224  false);
225 }
226 
227 void
229 {
230  // Fetch the segregated rhie-chow object and transfer some information about the momentum
231  // system(s)
233  &getUserObject<INSFVRhieChowInterpolatorSegregated>("rhie_chow_user_object"));
235 
236  // Initialize the face velocities in the RC object
238 }
239 
240 std::vector<std::pair<unsigned int, Real>>
242 {
243  // Temporary storage for the (flux-normalized) residuals form
244  // different momentum components
245  std::vector<std::pair<unsigned int, Real>> its_normalized_residuals;
246 
247  // We can create this here with the assumption that every momentum component has the same number
248  // of dofs
249  auto zero_solution = _momentum_systems[0]->system().current_local_solution->zero_clone();
250 
251  // Solve the momentum equations.
252  // TO DO: These equations are VERY similar. If we can store the differences (things coming from
253  // BCs for example) separately, it is enough to construct one matrix.
254  for (const auto system_i : index_range(_momentum_systems))
255  {
257 
258  // We will need the right hand side and the solution of the next component
259  NonlinearImplicitSystem & momentum_system =
260  libMesh::cast_ref<NonlinearImplicitSystem &>(_momentum_systems[system_i]->system());
261 
262  PetscLinearSolver<Real> & momentum_solver =
263  libMesh::cast_ref<PetscLinearSolver<Real> &>(*momentum_system.get_linear_solver());
264 
265  NumericVector<Number> & solution = *(momentum_system.solution);
266  NumericVector<Number> & rhs = *(momentum_system.rhs);
267  SparseMatrix<Number> & mmat = *(momentum_system.matrix);
268 
269  auto diff_diagonal = solution.zero_clone();
270 
271  // We plug zero in this to get the system matrix and the right hand side of the linear problem
272  _problem.computeResidualAndJacobian(*zero_solution, rhs, mmat);
273  // Sadly, this returns -b so we multiply with -1
274  rhs.scale(-1.0);
275 
276  // Still need to relax the right hand side with the same vector
277  NS::FV::relaxMatrix(mmat, _momentum_equation_relaxation, *diff_diagonal);
278  NS::FV::relaxRightHandSide(rhs, solution, *diff_diagonal);
279 
280  // The normalization factor depends on the right hand side so we need to recompute it for this
281  // component
282  Real norm_factor = NS::FV::computeNormalizationFactor(solution, mmat, rhs);
283 
284  // Very important, for deciding the convergence, we need the unpreconditioned
285  // norms in the linear solve
286  LibmeshPetscCall(KSPSetNormType(momentum_solver.ksp(), KSP_NORM_UNPRECONDITIONED));
287  // Solve this component. We don't update the ghosted solution yet, that will come at the end
288  // of the corrector step. Also setting the linear tolerances and maximum iteration counts.
290  momentum_solver.set_solver_configuration(_momentum_linear_control);
291 
292  // We solve the equation
293  auto its_resid_pair = momentum_solver.solve(mmat, mmat, solution, rhs);
294  momentum_system.update();
295 
296  // Save the normalized residual
297  its_normalized_residuals.push_back(
298  std::make_pair(its_resid_pair.first, momentum_solver.get_initial_residual() / norm_factor));
299 
300  if (_print_fields)
301  {
302  _console << " matrix when we solve " << std::endl;
303  mmat.print();
304  _console << " rhs when we solve " << std::endl;
305  rhs.print();
306  _console << " velocity solution component " << system_i << std::endl;
307  solution.print();
308  _console << "Norm factor " << norm_factor << std::endl;
309  _console << Moose::stringify(momentum_solver.get_initial_residual()) << std::endl;
310  }
311  }
312 
313  for (const auto system_i : index_range(_momentum_systems))
314  {
315  NonlinearImplicitSystem & momentum_system =
316  libMesh::cast_ref<NonlinearImplicitSystem &>(_momentum_systems[system_i]->system());
317  _momentum_systems[system_i]->setSolution(*(momentum_system.current_local_solution));
318  _momentum_systems[system_i]->copyPreviousNonlinearSolutions();
319  }
320 
321  return its_normalized_residuals;
322 }
323 
324 std::pair<unsigned int, Real>
326 {
328 
329  // We will need some members from the implicit nonlinear system
330  NonlinearImplicitSystem & pressure_system =
331  libMesh::cast_ref<NonlinearImplicitSystem &>(_pressure_system.system());
332 
333  // We will need the solution, the right hand side and the matrix
334  NumericVector<Number> & current_local_solution = *(pressure_system.current_local_solution);
335  NumericVector<Number> & solution = *(pressure_system.solution);
336  SparseMatrix<Number> & mmat = *(pressure_system.matrix);
337  NumericVector<Number> & rhs = *(pressure_system.rhs);
338 
339  // Fetch the linear solver from the system
340  PetscLinearSolver<Real> & pressure_solver =
341  libMesh::cast_ref<PetscLinearSolver<Real> &>(*pressure_system.get_linear_solver());
342 
343  // We need a zero vector to be able to emulate the Ax=b system by evaluating the
344  // residual and jacobian. Unfortunately, this will leave us with the -b on the right hand side
345  // so we correct it by multiplying it with (-1)
346  auto zero_solution = current_local_solution.zero_clone();
347  _problem.computeResidualAndJacobian(*zero_solution, rhs, mmat);
348  rhs.scale(-1.0);
349 
350  if (_print_fields)
351  {
352  _console << "Pressure matrix" << std::endl;
353  mmat.print();
354  }
355 
356  // We compute the normalization factors based on the fluxes
357  Real norm_factor = NS::FV::computeNormalizationFactor(solution, mmat, rhs);
358 
359  // We need the non-preconditioned norm to be consistent with the norm factor
360  LibmeshPetscCall(KSPSetNormType(pressure_solver.ksp(), KSP_NORM_UNPRECONDITIONED));
361 
362  // Setting the linear tolerances and maximum iteration counts
365 
366  if (_pin_pressure)
368 
369  auto its_res_pair = pressure_solver.solve(mmat, mmat, solution, rhs);
370  pressure_system.update();
371 
372  if (_print_fields)
373  {
374  _console << " rhs when we solve pressure " << std::endl;
375  rhs.print();
376  _console << " Pressure " << std::endl;
377  solution.print();
378  _console << "Norm factor " << norm_factor << std::endl;
379  }
380 
381  _pressure_system.setSolution(current_local_solution);
382 
383  return std::make_pair(its_res_pair.first, pressure_solver.get_initial_residual() / norm_factor);
384 }
385 
386 std::pair<unsigned int, Real>
388  NonlinearSystemBase & system,
389  const Real relaxation_factor,
390  libMesh::SolverConfiguration & solver_config,
391  const Real absolute_tol)
392 {
394 
395  // We will need some members from the implicit nonlinear system
396  NonlinearImplicitSystem & ni_system =
397  libMesh::cast_ref<NonlinearImplicitSystem &>(system.system());
398 
399  // We will need the solution, the right hand side and the matrix
400  NumericVector<Number> & current_local_solution = *(ni_system.current_local_solution);
401  NumericVector<Number> & solution = *(ni_system.solution);
402  SparseMatrix<Number> & mmat = *(ni_system.matrix);
403  NumericVector<Number> & rhs = *(ni_system.rhs);
404 
405  // We need a vector that stores the (diagonal_relaxed-original_diagonal) vector
406  auto diff_diagonal = solution.zero_clone();
407 
408  // Fetch the linear solver from the system
409  PetscLinearSolver<Real> & linear_solver =
410  libMesh::cast_ref<PetscLinearSolver<Real> &>(*ni_system.get_linear_solver());
411 
412  // We need a zero vector to be able to emulate the Ax=b system by evaluating the
413  // residual and jacobian. Unfortunately, this will leave us with the -b on the right hand side
414  // so we correct it by multiplying it with (-1)
415  auto zero_solution = current_local_solution.zero_clone();
416  _problem.computeResidualAndJacobian(*zero_solution, rhs, mmat);
417  rhs.scale(-1.0);
418 
419  // Go and relax the system matrix and the right hand side
420  NS::FV::relaxMatrix(mmat, relaxation_factor, *diff_diagonal);
421  NS::FV::relaxRightHandSide(rhs, solution, *diff_diagonal);
422 
423  if (_print_fields)
424  {
425  _console << system.name() << " system matrix" << std::endl;
426  mmat.print();
427  }
428 
429  // We compute the normalization factors based on the fluxes
430  Real norm_factor = NS::FV::computeNormalizationFactor(solution, mmat, rhs);
431 
432  // We need the non-preconditioned norm to be consistent with the norm factor
433  LibmeshPetscCall(KSPSetNormType(linear_solver.ksp(), KSP_NORM_UNPRECONDITIONED));
434 
435  // Setting the linear tolerances and maximum iteration counts
436  solver_config.real_valued_data["abs_tol"] = absolute_tol * norm_factor;
437  linear_solver.set_solver_configuration(solver_config);
438 
439  // Solve the system and update current local solution
440  auto its_res_pair = linear_solver.solve(mmat, mmat, solution, rhs);
441  ni_system.update();
442 
443  if (_print_fields)
444  {
445  _console << " rhs when we solve " << system.name() << std::endl;
446  rhs.print();
447  _console << system.name() << " solution " << std::endl;
448  solution.print();
449  _console << " Norm factor " << norm_factor << std::endl;
450  }
451 
452  system.setSolution(current_local_solution);
453 
454  return std::make_pair(its_res_pair.first, linear_solver.get_initial_residual() / norm_factor);
455 }
456 
457 std::pair<unsigned int, Real>
459 {
461 
462  // We will need some members from the implicit nonlinear system
463  NonlinearImplicitSystem & se_system =
464  libMesh::cast_ref<NonlinearImplicitSystem &>(_solid_energy_system->system());
465 
466  // We will need the solution, the right hand side and the matrix
467  NumericVector<Number> & current_local_solution = *(se_system.current_local_solution);
468  NumericVector<Number> & solution = *(se_system.solution);
469  SparseMatrix<Number> & mat = *(se_system.matrix);
470  NumericVector<Number> & rhs = *(se_system.rhs);
471 
472  // Fetch the linear solver from the system
473  PetscLinearSolver<Real> & se_solver =
474  libMesh::cast_ref<PetscLinearSolver<Real> &>(*se_system.get_linear_solver());
475 
476  // We need a zero vector to be able to emulate the Ax=b system by evaluating the
477  // residual and jacobian. Unfortunately, this will leave us with the -b on the righ hand side
478  // so we correct it by multiplying it with (-1)
479  auto zero_solution = current_local_solution.zero_clone();
480  _problem.computeResidualAndJacobian(*zero_solution, rhs, mat);
481  rhs.scale(-1.0);
482 
483  if (_print_fields)
484  {
485  _console << "Solid energy matrix" << std::endl;
486  mat.print();
487  }
488 
489  // We compute the normalization factors based on the fluxes
490  Real norm_factor = NS::FV::computeNormalizationFactor(solution, mat, rhs);
491 
492  // We need the non-preconditioned norm to be consistent with the norm factor
493  LibmeshPetscCall(KSPSetNormType(se_solver.ksp(), KSP_NORM_UNPRECONDITIONED));
494 
495  // Setting the linear tolerances and maximum iteration counts
498 
499  auto its_res_pair = se_solver.solve(mat, mat, solution, rhs);
500  se_system.update();
501 
502  if (_print_fields)
503  {
504  _console << " Solid energy rhs " << std::endl;
505  rhs.print();
506  _console << " Solid temperature " << std::endl;
507  solution.print();
508  _console << "Norm factor " << norm_factor << std::endl;
509  }
510 
511  _solid_energy_system->setSolution(current_local_solution);
512 
513  return std::make_pair(its_res_pair.first, se_solver.get_initial_residual() / norm_factor);
514 }
515 
516 bool
518 {
519  // Dummy solver parameter file which is needed for switching petsc options
520  SolverParams solver_params;
521  solver_params._type = Moose::SolveType::ST_LINEAR;
522  solver_params._line_search = Moose::LineSearchType::LS_NONE;
523 
524  // Initialize the quantities which matter in terms of the iteration
525  unsigned int iteration_counter = 0;
526 
527  // Assign residuals to general residual vector
528  unsigned int no_systems =
531  no_systems += _turbulence_systems.size();
532  std::vector<std::pair<unsigned int, Real>> ns_its_residuals(no_systems, std::make_pair(0, 1.0));
533  std::vector<Real> ns_abs_tols(_momentum_systems.size(), _momentum_absolute_tolerance);
534  ns_abs_tols.push_back(_pressure_absolute_tolerance);
535  if (_has_energy_system)
536  {
537  ns_abs_tols.push_back(_energy_absolute_tolerance);
539  ns_abs_tols.push_back(_solid_energy_absolute_tolerance);
540  }
542  for (auto system_i : index_range(_turbulence_absolute_tolerance))
543  ns_abs_tols.push_back(_turbulence_absolute_tolerance[system_i]);
544 
545  bool converged = false;
546  // Loop until converged or hit the maximum allowed iteration number
547  while (iteration_counter < _num_iterations && !converged)
548  {
549  iteration_counter++;
550  // Resdiual index
551  size_t residual_index = 0;
552 
553  // Execute all objects tagged as nonlinear
554  // This will execute everything in the problem at nonlinear, including the aux kernels.
555  // This way we compute the aux kernels before the momentum equations are solved.
557 
558  // We clear the caches in the momentum and pressure variables
559  for (auto system_i : index_range(_momentum_systems))
560  _momentum_systems[system_i]->residualSetup();
562 
563  // If we solve for energy, we clear the caches there too
564  if (_has_energy_system)
565  {
569  }
570 
571  // If we solve for turbulence, we clear the caches there too
573  for (auto system_i : index_range(_turbulence_systems))
574  _turbulence_systems[system_i]->residualSetup();
575 
576  // We set the preconditioner/controllable parameters through petsc options. Linear
577  // tolerances will be overridden within the solver. In case of a segregated momentum
578  // solver, we assume that every velocity component uses the same preconditioner
580 
581  // Solve the momentum predictor step
582  auto momentum_residual = solveMomentumPredictor();
583  for (const auto system_i : index_range(momentum_residual))
584  ns_its_residuals[system_i] = momentum_residual[system_i];
585 
586  // Compute the coupling fields between the momentum and pressure equations
588 
589  // We set the preconditioner/controllable parameters for the pressure equations through
590  // petsc options. Linear tolerances will be overridden within the solver.
592 
593  // Solve the pressure corrector
594  ns_its_residuals[momentum_residual.size()] = solvePressureCorrector();
595  // We need this to make sure we evaluate cell gradients for the nonorthogonal correction in
596  // the face velocity update
598 
599  // Compute the face velocity which is used in the advection terms
601 
602  auto & pressure_current_solution = *(_pressure_system.system().current_local_solution.get());
603  auto & pressure_old_solution = *(_pressure_system.solutionPreviousNewton());
604  // Relax the pressure update for the next momentum predictor
606  pressure_current_solution, pressure_old_solution, _pressure_variable_relaxation);
607 
608  // Overwrite old solution
609  pressure_old_solution = pressure_current_solution;
610  _pressure_system.setSolution(pressure_current_solution);
611 
612  // We clear out the caches so that the gradients can be computed with the relaxed solution
614 
615  // Reconstruct the cell velocity as well to accelerate convergence
617 
618  // Update residual index
619  residual_index = momentum_residual.size();
620 
621  // If we have an energy equation, solve it here. We assume the material properties in the
622  // Navier-Stokes equations depend on temperature, therefore we can not solve for temperature
623  // outside of the velocity-pressure loop
624  if (_has_energy_system)
625  {
626  // We set the preconditioner/controllable parameters through petsc options. Linear
627  // tolerances will be overridden within the solver.
629  residual_index += 1;
630  ns_its_residuals[residual_index] = solveAdvectedSystem(_energy_sys_number,
635 
637  {
638  // We set the preconditioner/controllable parameters through petsc options. Linear
639  // tolerances will be overridden within the solver.
641  residual_index += 1;
642  ns_its_residuals[residual_index] = solveSolidEnergySystem();
643  }
644  }
645 
646  // If we have an turbulence equations, we solve it here. We solve it inside the
647  // momentum-pressure loop because it affects the turbulent viscosity
649  {
651 
652  for (auto system_i : index_range(_turbulence_systems))
653  {
654  residual_index += 1;
655  ns_its_residuals[residual_index] =
657  *_turbulence_systems[system_i],
661 
662  auto & current_solution =
663  *(_turbulence_systems[system_i]->system().current_local_solution.get());
664  NS::FV::limitSolutionUpdate(current_solution, _turbulence_field_min_limit[system_i]);
665 
666  // Relax the turbulence update for the next momentum predictor
667  auto & old_solution = *(_turbulence_systems[system_i]->solutionPreviousNewton());
668 
669  // Relax the pressure update for the next momentum predictor
671  current_solution, old_solution, _turbulence_equation_relaxation[system_i]);
672 
673  // Overwrite old solution
674  old_solution = current_solution;
675  _turbulence_systems[system_i]->setSolution(current_solution);
676 
677  // We clear out the caches so that the gradients can be computed with the relaxed solution
678  _turbulence_systems[system_i]->residualSetup();
679  }
680  }
681 
682  // Printing residuals
683  residual_index = 0;
684  _console << "Iteration " << iteration_counter << " Initial residual norms:" << std::endl;
685  for (auto system_i : index_range(_momentum_systems))
686  _console << " Momentum equation:"
687  << (_momentum_systems.size() > 1
688  ? std::string(" Component ") + std::to_string(system_i + 1) +
689  std::string(" ")
690  : std::string(" "))
691  << COLOR_GREEN << ns_its_residuals[system_i].second << COLOR_DEFAULT << std::endl;
692  _console << " Pressure equation: " << COLOR_GREEN
693  << ns_its_residuals[momentum_residual.size()].second << COLOR_DEFAULT << std::endl;
694  residual_index = momentum_residual.size();
695 
696  if (_has_energy_system)
697  {
698  residual_index += 1;
699  _console << " Energy equation: " << COLOR_GREEN << ns_its_residuals[residual_index].second
700  << COLOR_DEFAULT << std::endl;
702  {
703  residual_index += 1;
704  _console << " Solid energy equation: " << COLOR_GREEN
705  << ns_its_residuals[residual_index].second << COLOR_DEFAULT << std::endl;
706  }
707  }
708 
710  {
711  _console << "Turbulence Iteration " << std::endl;
712  for (auto system_i : index_range(_turbulence_systems))
713  {
714  residual_index += 1;
715  _console << _turbulence_systems[system_i]->name() << " " << COLOR_GREEN
716  << ns_its_residuals[residual_index].second << COLOR_DEFAULT << std::endl;
717  }
718  }
719 
720  converged = NS::FV::converged(ns_its_residuals, ns_abs_tols);
721  }
722 
724 
725  // Now we solve for the passive scalar equations, they should not influence the solution of the
726  // system above. The reason why we need more than one iteration is due to the matrix relaxation
727  // which can be used to stabilize the equations
729  {
730  _console << " Passive Scalar Iteration " << iteration_counter << std::endl;
731 
732  // We set the options used by Petsc (preconditioners etc). We assume that every passive
733  // scalar equation uses the same options for now.
735 
736  iteration_counter = 0;
737  std::vector<std::pair<unsigned int, Real>> passive_scalar_residuals(
738  _passive_scalar_systems.size(), std::make_pair(0, 1.0));
739 
740  bool passive_scalar_converged =
741  NS::FV::converged(passive_scalar_residuals, _passive_scalar_absolute_tolerance);
742  while (iteration_counter < _num_iterations && !passive_scalar_converged)
743  {
744  // We clear the caches in the passive scalar variables
745  for (auto system_i : index_range(_passive_scalar_systems))
746  _passive_scalar_systems[system_i]->residualSetup();
747 
748  iteration_counter++;
749 
750  // Solve the passive scalar equations
751  for (auto system_i : index_range(_passive_scalar_systems))
752  passive_scalar_residuals[system_i] =
754  *_passive_scalar_systems[system_i],
758 
759  _console << "Iteration " << iteration_counter << " Initial residual norms:" << std::endl;
760  for (auto system_i : index_range(_passive_scalar_systems))
761  _console << _passive_scalar_systems[system_i]->name() << " " << COLOR_GREEN
762  << passive_scalar_residuals[system_i].second << COLOR_DEFAULT << std::endl;
763 
764  passive_scalar_converged =
765  NS::FV::converged(passive_scalar_residuals, _passive_scalar_absolute_tolerance);
766  }
767 
768  converged = _continue_on_max_its ? true : passive_scalar_converged;
769  }
770 
771  return converged;
772 }
773 
774 void
776 {
777  // check to make sure that we don't have any time kernels in this simulation (Steady State)
778  for (const auto system : _momentum_systems)
779  checkTimeKernels(*system);
780 
782 
783  if (_has_energy_system)
784  {
788  }
789 
791  for (const auto system : _passive_scalar_systems)
792  checkTimeKernels(*system);
793 
795  for (const auto system : _turbulence_systems)
796  checkTimeKernels(*system);
797 }
798 
799 void
801 {
802  // check to make sure that we don't have any time kernels in this simulation (Steady State)
803  if (system.containsTimeKernel())
804  mooseError("You have specified time kernels in your steady state simulation in system",
805  system.name());
806 }
const unsigned int _energy_sys_number
The number of the system corresponding to the energy equation.
MultiMooseEnum getCommonPetscKeys()
void relaxMatrix(SparseMatrix< Number > &matrix_in, const Real relaxation_parameter, NumericVector< Number > &diff_diagonal)
Relax the matrix to ensure diagonal dominance, we hold onto the difference in diagonals for later use...
FEProblemBase & _problem
SIMPLESolverConfiguration _momentum_linear_control
Options for the linear solver of the momentum equation.
const Real _energy_absolute_tolerance
The user-defined absolute tolerance for determining the convergence in energy.
const TagID _pressure_tag_id
The ID of the tag which corresponds to the pressure gradient terms in the momentum equation...
const unsigned int _num_iterations
The maximum number of momentum-pressure iterations.
SIMPLESolverConfiguration _energy_linear_control
Options for the linear solver of the energy equation.
const unsigned int invalid_uint
const bool _print_fields
Debug parameter which allows printing the coupling and solution vectors/matrices. ...
INSFVRhieChowInterpolatorSegregated * _rc_uo
Pointer to the segregated RhieChow interpolation object.
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
const Real _solid_energy_l_abs_tol
Absolute linear tolerance for the energy equations.
const bool _has_energy_system
Boolean for easy check if a fluid energy system shall be solved or not.
Moose::LineSearchType _line_search
Real computeNormalizationFactor(const NumericVector< Number > &solution, const SparseMatrix< Number > &mat, const NumericVector< Number > &rhs)
Compute a normalization factor which is applied to the linear residual to determine convergence...
const Real _passive_scalar_l_abs_tol
Absolute linear tolerance for the passive scalar equation(s).
std::vector< Real > _turbulence_field_min_limit
The user-defined lower limit for turbulent quantities e.g. k, eps/omega, etc..
const Real _momentum_equation_relaxation
The user-defined relaxation parameter for the momentum equation.
virtual std::unique_ptr< NumericVector< T > > zero_clone() const=0
std::vector< unsigned int > _turbulence_system_numbers
virtual void checkTimeKernels(NonlinearSystemBase &system)
Check if the system contains time kernels.
void checkDependentParameterError(const std::string &main_parameter, const std::vector< std::string > &dependent_parameters, const bool should_be_defined)
static InputParameters validParams()
NumericVector< Number > * rhs
const unsigned int _solid_energy_sys_number
The number of the system corresponding to the solid energy equation.
const Real _momentum_l_abs_tol
Absolute linear tolerance for the momentum equation(s).
void setSolution(const NumericVector< Number > &soln)
virtual LinearSolver< Number > * get_linear_solver() const
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &matrix_in, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const std::optional< double > tol=std::nullopt, const std::optional< unsigned int > m_its=std::nullopt) override
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
virtual void linkRhieChowUserObject() override
Fetch the Rhie Chow user object that is reponsible for determining face velocities and mass flux...
const Real _turbulence_l_abs_tol
Absolute linear tolerance for the turbulence equation(s).
NonlinearSystemBase & _pressure_system
Reference to the nonlinear system corresponding to the pressure equation.
void relaxRightHandSide(NumericVector< Number > &rhs_in, const NumericVector< Number > &solution_in, const NumericVector< Number > &diff_diagonal)
Relax the right hand side of an equation, this needs to be called once and the system matrix has been...
std::map< std::string, Real > real_valued_data
bool converged(const std::vector< std::pair< unsigned int, Real >> &residuals, const std::vector< Real > &abs_tolerances)
Based on the residuals, determine if the iterative process converged or not.
const Real _pressure_absolute_tolerance
The user-defined absolute tolerance for determining the convergence in pressure.
bool isParamValid(const std::string &name) const
const Real _pressure_variable_relaxation
The user-defined relaxation parameter for the pressure variable.
const Real _momentum_absolute_tolerance
The user-defined absolute tolerance for determining the convergence in momentum.
virtual const std::string & name() const
const std::vector< SolverSystemName > & _passive_scalar_system_names
The names of the passive scalar systems.
virtual bool containsTimeKernel() override
SIMPLESolverConfiguration _turbulence_linear_control
Options for the linear solver of the turbulence equation(s)
virtual void scale(const T factor)=0
void computeFaceVelocity()
Update the values of the face velocities in the containers.
SIMPLESolverConfiguration _passive_scalar_linear_control
Options for the linear solver of the passive scalar equation(s)
virtual void execute(const ExecFlagType &exec_type)
Moose::PetscSupport::PetscOptions _energy_petsc_options
Options which hold the petsc settings for the fluid energy equation.
SIMPLESolverConfiguration _pressure_linear_control
Options for the linear solver of the pressure equation.
MultiMooseEnum getCommonPetscFlags()
void setCurrentNonlinearSystem(const unsigned int nl_sys_num)
Moose::PetscSupport::PetscOptions _passive_scalar_petsc_options
Options which hold the petsc settings for the passive scalar equation(s)
virtual unsigned int dimension() const
std::map< std::string, int > int_valued_data
virtual std::unique_ptr< SparseMatrix< T > > zero_clone() const=0
const bool _has_turbulence_systems
Boolean for easy check if turbulence systems shall be solved or not.
Moose::PetscSupport::PetscOptions _momentum_petsc_options
Options which hold the petsc settings for the momentum equation.
const Real _pressure_l_abs_tol
Absolute linear tolerance for the pressure equation.
const std::vector< Real > _passive_scalar_absolute_tolerance
The user-defined absolute tolerance for determining the convergence in passive scalars.
std::unique_ptr< NumericVector< Number > > solution
dof_id_type _pressure_pin_dof
The dof ID where the pressure needs to be pinned.
Moose::SolveType _type
const std::vector< SolverSystemName > & _turbulence_system_names
The names of the turbulence scalar systems.
virtual void print(std::ostream &os=libMesh::out) const
const bool _pin_pressure
If the pressure needs to be pinned.
NonlinearSystemBase & getNonlinearSystemBase(const unsigned int sys_num)
void paramError(const std::string &param, Args... args) const
const Real _energy_l_abs_tol
Absolute linear tolerance for the energy equations.
std::string stringify(const T &t)
std::pair< unsigned int, Real > solveSolidEnergySystem()
Solve the solid energy conservation equation.
void petscSetOptions(const PetscOptions &po, const SolverParams &solver_params, FEProblemBase *const problem=nullptr)
void set_solver_configuration(SolverConfiguration &solver_configuration)
const std::vector< SolverSystemName > & _momentum_system_names
The names of the momentum systems.
const ExecFlagType EXEC_NONLINEAR
NonlinearSystemBase * _solid_energy_system
Pointer to the nonlinear system corresponding to the solid energy equation.
const std::vector< Real > _passive_scalar_equation_relaxation
The user-defined relaxation parameter(s) for the passive scalar equation(s)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool solve() override
Performs the momentum pressure coupling.
NonlinearSystemBase * _energy_system
Pointer to the nonlinear system corresponding to the fluid energy equation.
void addPetscPairsToPetscOptions(const std::vector< std::pair< MooseEnumItem, std::string >> &petsc_pair_options, const unsigned int mesh_dimension, const std::string &prefix, const ParallelParamObject &param_object, PetscOptions &petsc_options)
const Real _energy_equation_relaxation
The user-defined relaxation parameter for the energy equation.
SparseMatrix< Number > * matrix
std::vector< NonlinearSystemBase * > _momentum_systems
Pointer(s) to the system(s) corresponding to the momentum equation(s)
const Real _pressure_pin_value
The value we want to enforce for pressure.
void computeResidualAndJacobian(const NumericVector< libMesh::Number > &soln, NumericVector< libMesh::Number > &residual, libMesh::SparseMatrix< libMesh::Number > &jacobian)
virtual MooseMesh & mesh() override
void computeCellVelocity()
Update the cell values of the velocity variables.
void mooseError(Args &&... args) const
std::pair< unsigned int, Real > solveAdvectedSystem(const unsigned int system_num, NonlinearSystemBase &system, const Real relaxation_factor, libMesh::SolverConfiguration &solver_config, const Real abs_tol)
Solve an equation which contains an advection term that depends on the solution of the segregated Nav...
std::unique_ptr< NumericVector< Number > > current_local_solution
virtual const NumericVector< Number > * solutionPreviousNewton() const
Solve class serving as a base class for the two SIMPLE solvers that operate with different assembly a...
const bool _has_passive_scalar_systems
Boolean for easy check if a passive scalar systems shall be solved or not.
void addPetscFlagsToPetscOptions(const MultiMooseEnum &petsc_flags, const std::string &prefix, const ParallelParamObject &param_object, PetscOptions &petsc_options)
void addRangeCheckedParam(const std::string &name, const T &value, const std::string &parsed_function, const std::string &doc_string)
const ConsoleStream _console
const Real _solid_energy_absolute_tolerance
The user-defined absolute tolerance for determining the convergence in solid energy.
SIMPLESolverConfiguration _solid_energy_linear_control
Options for the linear solver of the energy equation.
const bool _continue_on_max_its
If solve should continue if maximum number of iterations is hit.
Moose::PetscSupport::PetscOptions _turbulence_petsc_options
Options which hold the petsc settings for the turbulence equation(s)
const unsigned int _pressure_sys_number
The number of the system corresponding to the pressure equation.
void constrainSystem(SparseMatrix< Number > &mx, NumericVector< Number > &rhs, const Real desired_value, const dof_id_type dof_id)
Implicitly constrain the system by adding a factor*(u-u_desired) to it at a desired dof value...
void prefix_with_name(bool value)
const std::vector< Real > _turbulence_absolute_tolerance
The user-defined absolute tolerance for determining the convergence in turbulence equations...
void limitSolutionUpdate(NumericVector< Number > &solution, const Real min_limit=std::numeric_limits< Real >::epsilon(), const Real max_limit=1e10)
Limit a solution to its minimum and maximum bounds: $u = min(max(u, min_limit), max_limit)$.
void relaxSolutionUpdate(NumericVector< Number > &vec_new, const NumericVector< Number > &vec_old, const Real relaxation_factor)
Relax the update on a solution field using the following approach: $u = u_{old}+ (u - u_{old})$...
void computeHbyA(bool verbose)
Computes the inverse of the digaonal (1/A) of the system matrix plus the H/A components for the press...
const std::vector< Real > _turbulence_equation_relaxation
The user-defined relaxation parameter(s) for the turbulence equation(s)
std::vector< NonlinearSystemBase * > _passive_scalar_systems
Pointer(s) to the system(s) corresponding to the passive scalar equation(s)
auto index_range(const T &sizable)
std::vector< NonlinearSystemBase * > _turbulence_systems
Pointer(s) to the system(s) corresponding to the turbulence equation(s)
void print(std::ostream &os=libMesh::out, const bool sparse=false) const
std::vector< unsigned int > _momentum_system_numbers
The number(s) of the system(s) corresponding to the momentum equation(s)
virtual void checkIntegrity() override
Check if the user defined time kernels.
std::vector< unsigned int > _passive_scalar_system_numbers
virtual unsigned int nlSysNum(const NonlinearSystemName &nl_sys_name) const override
Moose::PetscSupport::PetscOptions _pressure_petsc_options
Options which hold the petsc settings for the pressure equation.
A user object which implements the Rhie Chow interpolation for segregated momentum-pressure systems...
void linkMomentumSystem(std::vector< NonlinearSystemBase *> momentum_systems, const std::vector< unsigned int > &momentum_system_numbers, const TagID pressure_gradient_tag)
Update the momentum system-related information.
const bool _has_solid_energy_system
Boolean for easy check if a solid energy system shall be solved or not.
virtual std::vector< std::pair< unsigned int, Real > > solveMomentumPredictor() override
Solve a momentum predictor step with a fixed pressure field.
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)
virtual std::pair< unsigned int, Real > solvePressureCorrector() override
Solve a pressure corrector step.
void initFaceVelocities()
Initialize the container for face velocities.
virtual void residualSetup() override
virtual libMesh::System & system() override