https://mooseframework.inl.gov
LinearAssemblySegregatedSolve.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 "LinearSystem.h"
14 
15 using namespace libMesh;
16 
19 {
21 
22  params.addParam<std::vector<SolverSystemName>>(
23  "active_scalar_systems", {}, "The solver system for each active scalar advection equation.");
24 
25  /*
26  * Parameters to control the solution of each scalar advection system
27  */
28  params.addParam<std::vector<Real>>("active_scalar_equation_relaxation",
29  std::vector<Real>(),
30  "The relaxation which should be used for the active scalar "
31  "equations. (=1 for no relaxation, "
32  "diagonal dominance will still be enforced)");
33 
34  params.addParam<MultiMooseEnum>("active_scalar_petsc_options",
36  "Singleton PETSc options for the active scalar equation(s)");
37  params.addParam<MultiMooseEnum>(
38  "active_scalar_petsc_options_iname",
40  "Names of PETSc name/value pairs for the active scalar equation(s)");
41  params.addParam<std::vector<std::string>>(
42  "active_scalar_petsc_options_value",
43  "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\" for the "
44  "active scalar equation(s)");
45  params.addParam<std::vector<Real>>(
46  "active_scalar_absolute_tolerance",
47  std::vector<Real>(),
48  "The absolute tolerance(s) on the normalized residual(s) of the active scalar equation(s).");
49  params.addRangeCheckedParam<Real>("active_scalar_l_tol",
50  1e-5,
51  "0.0<=active_scalar_l_tol & active_scalar_l_tol<1.0",
52  "The relative tolerance on the normalized residual in the "
53  "linear solver of the active scalar equation(s).");
54  params.addRangeCheckedParam<Real>("active_scalar_l_abs_tol",
55  1e-10,
56  "0.0<active_scalar_l_abs_tol",
57  "The absolute tolerance on the normalized residual in the "
58  "linear solver of the active scalar equation(s).");
59  params.addParam<unsigned int>(
60  "active_scalar_l_max_its",
61  10000,
62  "The maximum allowed iterations in the linear solver of the turbulence equation.");
63 
64  params.addParamNamesToGroup(
65  "active_scalar_systems active_scalar_equation_relaxation active_scalar_petsc_options "
66  "active_scalar_petsc_options_iname "
67  "active_scalar_petsc_options_value active_scalar_petsc_options_value "
68  "active_scalar_absolute_tolerance "
69  "active_scalar_l_tol active_scalar_l_abs_tol active_scalar_l_max_its",
70  "Active Scalars Equations");
71 
72  return params;
73 }
74 
76  : SIMPLESolveBase(ex),
77  _pressure_sys_number(_problem.linearSysNum(getParam<SolverSystemName>("pressure_system"))),
78  _pressure_system(_problem.getLinearSystem(_pressure_sys_number)),
79  _energy_sys_number(_has_energy_system
80  ? _problem.linearSysNum(getParam<SolverSystemName>("energy_system"))
82  _energy_system(_has_energy_system ? &_problem.getLinearSystem(_energy_sys_number) : nullptr),
83  _solid_energy_sys_number(
84  _has_solid_energy_system
85  ? _problem.linearSysNum(getParam<SolverSystemName>("solid_energy_system"))
87  _solid_energy_system(
88  _has_solid_energy_system ? &_problem.getLinearSystem(_solid_energy_sys_number) : nullptr),
89  _active_scalar_system_names(getParam<std::vector<SolverSystemName>>("active_scalar_systems")),
90  _has_active_scalar_systems(!_active_scalar_system_names.empty()),
91  _active_scalar_equation_relaxation(
92  getParam<std::vector<Real>>("active_scalar_equation_relaxation")),
93  _active_scalar_l_abs_tol(getParam<Real>("active_scalar_l_abs_tol")),
94  _active_scalar_absolute_tolerance(
95  getParam<std::vector<Real>>("active_scalar_absolute_tolerance"))
96 {
97  // We fetch the systems and their numbers for the momentum equations.
98  for (auto system_i : index_range(_momentum_system_names))
99  {
102  _systems_to_solve.push_back(_momentum_systems.back());
103  }
104 
106 
107  if (_has_energy_system)
109 
112 
113  // and for the passive scalar equations
115  for (auto system_i : index_range(_passive_scalar_system_names))
116  {
119  _passive_scalar_systems.push_back(
121  _systems_to_solve.push_back(_passive_scalar_systems.back());
122  }
123 
124  // We disable the prefix here for the time being, the segregated solvers use a different approach
125  // for setting the petsc parameters
126  for (auto & system : _systems_to_solve)
127  system->system().prefix_with_name(false);
128 
129  // and for the active scalar equations
131  for (auto system_i : index_range(_active_scalar_system_names))
132  {
135  _active_scalar_systems.push_back(
137 
138  const auto & active_scalar_petsc_options =
139  getParam<MultiMooseEnum>("active_scalar_petsc_options");
140  const auto & active_scalar_petsc_pair_options = getParam<MooseEnumItem, std::string>(
141  "active_scalar_petsc_options_iname", "active_scalar_petsc_options_value");
143  active_scalar_petsc_options, "-", *this, _active_scalar_petsc_options);
144  Moose::PetscSupport::addPetscPairsToPetscOptions(active_scalar_petsc_pair_options,
145  _problem.mesh().dimension(),
146  "-",
147  *this,
149 
151  getParam<Real>("active_scalar_l_tol");
153  getParam<Real>("active_scalar_l_abs_tol");
155  getParam<unsigned int>("active_scalar_l_max_its");
156  }
157 
159  paramError("active_scalar_equation_relaxation",
160  "Should be the same size as the number of systems");
161 }
162 
163 void
165 {
166  _rc_uo =
167  const_cast<RhieChowMassFlux *>(&getUserObject<RhieChowMassFlux>("rhie_chow_user_object"));
170 
171  // Initialize the face velocities in the RC object
172  if (!_app.isRecovering())
175 }
176 
177 std::vector<std::pair<unsigned int, Real>>
179 {
180  // Temporary storage for the (flux-normalized) residuals from
181  // different momentum components
182  std::vector<std::pair<unsigned int, Real>> its_normalized_residuals;
183 
184  LinearImplicitSystem & momentum_system_0 =
185  libMesh::cast_ref<LinearImplicitSystem &>(_momentum_systems[0]->system());
186 
187  PetscLinearSolver<Real> & momentum_solver =
188  libMesh::cast_ref<PetscLinearSolver<Real> &>(*momentum_system_0.get_linear_solver());
189 
190  // Solve the momentum equations.
191  // TO DO: These equations are VERY similar. If we can store the differences (things coming from
192  // BCs for example) separately, it is enough to construct one matrix.
193  for (const auto system_i : index_range(_momentum_systems))
194  {
196 
197  // We will need the right hand side and the solution of the next component
198  LinearImplicitSystem & momentum_system =
199  libMesh::cast_ref<LinearImplicitSystem &>(_momentum_systems[system_i]->system());
200 
201  NumericVector<Number> & solution = *(momentum_system.solution);
202  NumericVector<Number> & rhs = *(momentum_system.rhs);
203  SparseMatrix<Number> & mmat = *(momentum_system.matrix);
204 
205  auto diff_diagonal = solution.zero_clone();
206 
207  // We assemble the matrix and the right hand side
208  _problem.computeLinearSystemSys(momentum_system, mmat, rhs, /*compute_grads*/ true);
209 
210  // Still need to relax the right hand side with the same vector
211  NS::FV::relaxMatrix(mmat, _momentum_equation_relaxation, *diff_diagonal);
212  NS::FV::relaxRightHandSide(rhs, solution, *diff_diagonal);
213 
214  // The normalization factor depends on the right hand side so we need to recompute it for this
215  // component
216  Real norm_factor = NS::FV::computeNormalizationFactor(solution, mmat, rhs);
217 
218  // Very important, for deciding the convergence, we need the unpreconditioned
219  // norms in the linear solve
220  LibmeshPetscCall(KSPSetNormType(momentum_solver.ksp(), KSP_NORM_UNPRECONDITIONED));
221  // Solve this component. We don't update the ghosted solution yet, that will come at the end
222  // of the corrector step. Also setting the linear tolerances and maximum iteration counts.
224  momentum_solver.set_solver_configuration(_momentum_linear_control);
225 
226  // We solve the equation
227  auto its_resid_pair = momentum_solver.solve(mmat, mmat, solution, rhs);
228  momentum_system.update();
229 
230  // We will reuse the preconditioner for every momentum system
231  if (system_i == 0)
232  momentum_solver.reuse_preconditioner(true);
233 
234  // Save the normalized residual
235  its_normalized_residuals.push_back(
236  std::make_pair(its_resid_pair.first, momentum_solver.get_initial_residual() / norm_factor));
237 
238  if (_print_fields)
239  {
240  _console << " solution after solve " << std::endl;
241  solution.print();
242  _console << " matrix when we solve " << std::endl;
243  mmat.print();
244  _console << " rhs when we solve " << std::endl;
245  rhs.print();
246  _console << " velocity solution component " << system_i << std::endl;
247  solution.print();
248  _console << "Norm factor " << norm_factor << std::endl;
249  _console << Moose::stringify(momentum_solver.get_initial_residual()) << std::endl;
250  }
251 
252  // Printing residuals
253  _console << " Momentum equation:"
254  << (_momentum_systems.size() > 1
255  ? std::string(" Component ") + std::to_string(system_i + 1) + std::string(" ")
256  : std::string(" "))
257  << COLOR_GREEN << its_normalized_residuals[system_i].second << COLOR_DEFAULT
258  << " Linear its: " << its_normalized_residuals[system_i].first << std::endl;
259  }
260 
261  for (const auto system_i : index_range(_momentum_systems))
262  {
263  LinearImplicitSystem & momentum_system =
264  libMesh::cast_ref<LinearImplicitSystem &>(_momentum_systems[system_i]->system());
265  _momentum_systems[system_i]->setSolution(*(momentum_system.current_local_solution));
266  _momentum_systems[system_i]->copyPreviousNonlinearSolutions();
267  }
268 
269  // We reset this to ensure the preconditioner is recomputed new time we go to the momentum
270  // predictor
271  momentum_solver.reuse_preconditioner(false);
272 
273  return its_normalized_residuals;
274 }
275 
276 std::pair<unsigned int, Real>
278 {
280 
281  // We will need some members from the linear system
282  LinearImplicitSystem & pressure_system =
283  libMesh::cast_ref<LinearImplicitSystem &>(_pressure_system.system());
284 
285  // We will need the solution, the right hand side and the matrix
286  NumericVector<Number> & current_local_solution = *(pressure_system.current_local_solution);
287  NumericVector<Number> & solution = *(pressure_system.solution);
288  SparseMatrix<Number> & mmat = *(pressure_system.matrix);
289  NumericVector<Number> & rhs = *(pressure_system.rhs);
290 
291  // Fetch the linear solver from the system
292  PetscLinearSolver<Real> & pressure_solver =
293  libMesh::cast_ref<PetscLinearSolver<Real> &>(*pressure_system.get_linear_solver());
294 
295  _problem.computeLinearSystemSys(pressure_system, mmat, rhs, false);
296 
297  if (_print_fields)
298  {
299  _console << "Pressure matrix" << std::endl;
300  mmat.print();
301  }
302 
303  // We compute the normalization factors based on the fluxes
304  Real norm_factor = NS::FV::computeNormalizationFactor(solution, mmat, rhs);
305 
306  // We need the non-preconditioned norm to be consistent with the norm factor
307  LibmeshPetscCall(KSPSetNormType(pressure_solver.ksp(), KSP_NORM_UNPRECONDITIONED));
308 
309  // Setting the linear tolerances and maximum iteration counts
312 
313  if (_pin_pressure)
315  pressure_system.update();
316 
317  auto its_res_pair = pressure_solver.solve(mmat, mmat, solution, rhs);
318  pressure_system.update();
319 
320  if (_print_fields)
321  {
322  _console << " rhs when we solve pressure " << std::endl;
323  rhs.print();
324  _console << " Pressure " << std::endl;
325  solution.print();
326  _console << "Norm factor " << norm_factor << std::endl;
327  }
328 
329  _pressure_system.setSolution(current_local_solution);
330 
331  const auto residuals =
332  std::make_pair(its_res_pair.first, pressure_solver.get_initial_residual() / norm_factor);
333 
334  _console << " Pressure equation: " << COLOR_GREEN << residuals.second << COLOR_DEFAULT
335  << " Linear its: " << residuals.first << std::endl;
336 
337  return residuals;
338 }
339 
340 std::pair<unsigned int, Real>
342 {
344 
345  // We will need some members from the linear system
346  LinearImplicitSystem & system =
347  libMesh::cast_ref<LinearImplicitSystem &>(_solid_energy_system->system());
348 
349  // We will need the solution, the right hand side and the matrix
350  NumericVector<Number> & current_local_solution = *(system.current_local_solution);
351  NumericVector<Number> & solution = *(system.solution);
352  SparseMatrix<Number> & mmat = *(system.matrix);
353  NumericVector<Number> & rhs = *(system.rhs);
354 
355  // Fetch the linear solver from the system
356  PetscLinearSolver<Real> & solver =
357  libMesh::cast_ref<PetscLinearSolver<Real> &>(*system.get_linear_solver());
358 
359  _problem.computeLinearSystemSys(system, mmat, rhs, false);
360 
361  if (_print_fields)
362  {
363  _console << "Solid energy matrix" << std::endl;
364  mmat.print();
365  }
366 
367  // We compute the normalization factors based on the fluxes
368  Real norm_factor = NS::FV::computeNormalizationFactor(solution, mmat, rhs);
369 
370  // We need the non-preconditioned norm to be consistent with the norm factor
371  LibmeshPetscCall(KSPSetNormType(solver.ksp(), KSP_NORM_UNPRECONDITIONED));
372 
373  // Setting the linear tolerances and maximum iteration counts
376 
377  auto its_res_pair = solver.solve(mmat, mmat, solution, rhs);
378  system.update();
379 
380  if (_print_fields)
381  {
382  _console << " rhs when we solve solid energy " << std::endl;
383  rhs.print();
384  _console << " Solid energy " << std::endl;
385  solution.print();
386  _console << "Norm factor " << norm_factor << std::endl;
387  }
388 
389  _solid_energy_system->setSolution(current_local_solution);
390 
391  const auto residuals =
392  std::make_pair(its_res_pair.first, solver.get_initial_residual() / norm_factor);
393 
394  _console << " Solid energy equation: " << COLOR_GREEN << residuals.second << COLOR_DEFAULT
395  << " Linear its: " << residuals.first << std::endl;
396 
397  return residuals;
398 }
399 
400 std::pair<unsigned int, Real>
401 LinearAssemblySegregatedSolve::correctVelocity(const bool subtract_updated_pressure,
402  const bool recompute_face_mass_flux,
403  const SolverParams & solver_params)
404 {
405  // Compute the coupling fields between the momentum and pressure equations.
406  // The first argument makes sure the pressure gradient is staged at the first
407  // iteration
408  _rc_uo->computeHbyA(subtract_updated_pressure, _print_fields);
409 
410  // We set the preconditioner/controllable parameters for the pressure equations through
411  // petsc options. Linear tolerances will be overridden within the solver.
413 
414  // Solve the pressure corrector
415  const auto residuals = solvePressureCorrector();
416 
417  // Compute the face velocity which is used in the advection terms. In certain
418  // segregated solver algorithms (like PISO) this is only done on the last iteration.
419  if (recompute_face_mass_flux)
421 
422  auto & pressure_current_solution = *(_pressure_system.system().current_local_solution.get());
423  auto & pressure_old_solution = *(_pressure_system.solutionPreviousNewton());
424 
425  // Relax the pressure update for the next momentum predictor
427  pressure_current_solution, pressure_old_solution, _pressure_variable_relaxation);
428 
429  // Overwrite old solution
430  pressure_old_solution = pressure_current_solution;
431  _pressure_system.setSolution(pressure_current_solution);
432 
433  // We recompute the updated pressure gradient
435 
436  // Reconstruct the cell velocity as well to accelerate convergence
438 
439  return residuals;
440 }
441 
442 std::pair<unsigned int, Real>
444  LinearSystem & system,
445  const Real relaxation_factor,
446  SolverConfiguration & solver_config,
447  const Real absolute_tol)
448 {
449  _problem.setCurrentLinearSystem(system_num);
450 
451  // We will need some members from the implicit linear system
452  LinearImplicitSystem & li_system = libMesh::cast_ref<LinearImplicitSystem &>(system.system());
453 
454  // We will need the solution, the right hand side and the matrix
455  NumericVector<Number> & current_local_solution = *(li_system.current_local_solution);
456  NumericVector<Number> & solution = *(li_system.solution);
457  SparseMatrix<Number> & mmat = *(li_system.matrix);
458  NumericVector<Number> & rhs = *(li_system.rhs);
459 
460  // We need a vector that stores the (diagonal_relaxed-original_diagonal) vector
461  auto diff_diagonal = solution.zero_clone();
462 
463  // Fetch the linear solver from the system
464  PetscLinearSolver<Real> & linear_solver =
465  libMesh::cast_ref<PetscLinearSolver<Real> &>(*li_system.get_linear_solver());
466 
467  _problem.computeLinearSystemSys(li_system, mmat, rhs, true);
468 
469  // Go and relax the system matrix and the right hand side
470  NS::FV::relaxMatrix(mmat, relaxation_factor, *diff_diagonal);
471  NS::FV::relaxRightHandSide(rhs, solution, *diff_diagonal);
472 
473  if (_print_fields)
474  {
475  _console << system.name() << " system matrix" << std::endl;
476  mmat.print();
477  }
478 
479  // We compute the normalization factors based on the fluxes
480  Real norm_factor = NS::FV::computeNormalizationFactor(solution, mmat, rhs);
481 
482  // We need the non-preconditioned norm to be consistent with the norm factor
483  LibmeshPetscCall(KSPSetNormType(linear_solver.ksp(), KSP_NORM_UNPRECONDITIONED));
484 
485  // Setting the linear tolerances and maximum iteration counts
486  solver_config.real_valued_data["abs_tol"] = absolute_tol * norm_factor;
487  linear_solver.set_solver_configuration(solver_config);
488 
489  // Solve the system and update current local solution
490  auto its_res_pair = linear_solver.solve(mmat, mmat, solution, rhs);
491  li_system.update();
492 
493  if (_print_fields)
494  {
495  _console << " rhs when we solve " << system.name() << std::endl;
496  rhs.print();
497  _console << system.name() << " solution " << std::endl;
498  solution.print();
499  _console << " Norm factor " << norm_factor << std::endl;
500  }
501 
502  system.setSolution(current_local_solution);
503 
504  const auto residuals =
505  std::make_pair(its_res_pair.first, linear_solver.get_initial_residual() / norm_factor);
506 
507  _console << " Advected system: " << system.name() << " " << COLOR_GREEN << residuals.second
508  << COLOR_DEFAULT << " Linear its: " << residuals.first << std::endl;
509 
510  return residuals;
511 }
512 
513 bool
515 {
516  // Do not solve if problem is set not to
517  if (!_problem.shouldSolve())
518  return true;
519 
520  // Dummy solver parameter file which is needed for switching petsc options
521  SolverParams solver_params;
522  solver_params._type = Moose::SolveType::ST_LINEAR;
523  solver_params._line_search = Moose::LineSearchType::LS_NONE;
524 
525  // Initialize the SIMPLE iteration counter
526  unsigned int simple_iteration_counter = 0;
527 
528  // Assign residuals to general residual vector
529  const unsigned int no_systems = _momentum_systems.size() + 1 + _has_energy_system +
531  std::vector<std::pair<unsigned int, Real>> ns_residuals(no_systems, std::make_pair(0, 1.0));
532  std::vector<Real> ns_abs_tols(_momentum_systems.size(), _momentum_absolute_tolerance);
533  ns_abs_tols.push_back(_pressure_absolute_tolerance);
534  if (_has_energy_system)
535  ns_abs_tols.push_back(_energy_absolute_tolerance);
537  ns_abs_tols.push_back(_solid_energy_absolute_tolerance);
539  for (const auto scalar_tol : _active_scalar_absolute_tolerance)
540  ns_abs_tols.push_back(scalar_tol);
541 
542  bool converged = false;
543  // Loop until converged or hit the maximum allowed iteration number
544  while (simple_iteration_counter < _num_iterations && !converged)
545  {
546  simple_iteration_counter++;
547 
548  // We set the preconditioner/controllable parameters through petsc options. Linear
549  // tolerances will be overridden within the solver. In case of a segregated momentum
550  // solver, we assume that every velocity component uses the same preconditioner
552 
553  // Initialize pressure gradients, after this we just reuse the last ones from each
554  // iteration
555  if (simple_iteration_counter == 1)
557 
558  _console << "Iteration " << simple_iteration_counter << " Initial residual norms:" << std::endl;
559 
560  // Solve the momentum predictor step
561  auto momentum_residual = solveMomentumPredictor();
562  for (const auto system_i : index_range(momentum_residual))
563  ns_residuals[system_i] = momentum_residual[system_i];
564 
565  // Now we correct the velocity, this function depends on the method, it differs for
566  // SIMPLE/PIMPLE, this returns the pressure errors
567  ns_residuals[momentum_residual.size()] = correctVelocity(true, true, solver_params);
568 
569  // If we have an energy equation, solve it here.We assume the material properties in the
570  // Navier-Stokes equations depend on temperature, therefore we can not solve for temperature
571  // outside of the velocity-pressure loop
572  if (_has_energy_system)
573  {
574  // We set the preconditioner/controllable parameters through petsc options. Linear
575  // tolerances will be overridden within the solver.
577  ns_residuals[momentum_residual.size() + _has_energy_system] =
583  }
585  {
586  // We set the preconditioner/controllable parameters through petsc options. Linear
587  // tolerances will be overridden within the solver.
589  ns_residuals[momentum_residual.size() + _has_solid_energy_system + _has_energy_system] =
591  }
592  // If we have active scalar equations, solve them here in case they depend on temperature
593  // or they affect the fluid properties such that they must be solved concurrently with pressure
594  // and velocity
596  {
598 
599  // We set the preconditioner/controllable parameters through petsc options. Linear
600  // tolerances will be overridden within the solver.
602  for (const auto i : index_range(_active_scalar_system_names))
603  ns_residuals[momentum_residual.size() + 1 + _has_energy_system + i] =
609  }
610  else
612 
613  converged = NS::FV::converged(ns_residuals, ns_abs_tols);
614  }
615 
616  // If we have passive scalar equations, solve them here. We assume the material properties in the
617  // Navier-Stokes equations do not depend on passive scalars, as they are passive, therefore we
618  // solve outside of the velocity-pressure loop
620  {
621  // The reason why we need more than one iteration is due to the matrix relaxation
622  // which can be used to stabilize the equations
623  bool passive_scalar_converged = false;
624  unsigned int ps_iteration_counter = 0;
625 
626  _console << "Passive scalar iteration " << ps_iteration_counter
627  << " Initial residual norms:" << std::endl;
628 
629  while (ps_iteration_counter < _num_iterations && !passive_scalar_converged)
630  {
631  ps_iteration_counter++;
632  std::vector<std::pair<unsigned int, Real>> scalar_residuals(
633  _passive_scalar_system_names.size(), std::make_pair(0, 1.0));
634  std::vector<Real> scalar_abs_tols;
635  for (const auto scalar_tol : _passive_scalar_absolute_tolerance)
636  scalar_abs_tols.push_back(scalar_tol);
637 
638  // We set the preconditioner/controllable parameters through petsc options. Linear
639  // tolerances will be overridden within the solver.
641  for (const auto i : index_range(_passive_scalar_system_names))
642  scalar_residuals[i] = solveAdvectedSystem(_passive_scalar_system_numbers[i],
647 
648  passive_scalar_converged = NS::FV::converged(scalar_residuals, scalar_abs_tols);
649  }
650 
651  // Both flow and scalars must converge
652  converged = passive_scalar_converged && converged;
653  }
654 
656 
657  return converged;
658 }
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...
bool shouldSolve() const
const std::vector< Real > _active_scalar_equation_relaxation
The user-defined relaxation parameter(s) for the active scalar equation(s)
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 unsigned int _num_iterations
The maximum number of momentum-pressure iterations.
User object responsible for determining the face fluxes using the Rhie-Chow interpolation in a segreg...
SIMPLESolverConfiguration _energy_linear_control
Options for the linear solver of the energy equation.
const unsigned int invalid_uint
void computeGradients()
const bool _print_fields
Debug parameter which allows printing the coupling and solution vectors/matrices. ...
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...
SIMPLESolverConfiguration _active_scalar_linear_control
Options for the linear solver of the active scalar equation(s)
const Real _passive_scalar_l_abs_tol
Absolute linear tolerance for the passive scalar equation(s).
std::pair< unsigned int, Real > solveSolidEnergy()
Solve an equation which contains the solid energy conservation.
std::vector< LinearSystem * > _passive_scalar_systems
Pointer(s) to the system(s) corresponding to the passive scalar equation(s)
Moose::PetscSupport::PetscOptions _solid_energy_petsc_options
Options which hold the petsc settings for the fluid energy equation.
const Real _momentum_equation_relaxation
The user-defined relaxation parameter for the momentum equation.
virtual std::unique_ptr< NumericVector< T > > zero_clone() const=0
void computeFaceMassFlux()
Update the values of the face velocities in the containers.
static InputParameters validParams()
LinearSystem * _solid_energy_system
Pointer to the nonlinear system corresponding to the solid energy equation.
NumericVector< Number > * rhs
const Real _momentum_l_abs_tol
Absolute linear tolerance for the momentum equation(s).
void setSolution(const NumericVector< Number > &soln)
const bool _has_active_scalar_systems
Boolean for easy check if a active scalar systems shall be solved or not.
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
const Real _active_scalar_l_abs_tol
Absolute linear tolerance for the active scalar equation(s).
virtual LinearSolver< Number > * get_linear_solver() const override
RhieChowMassFlux * _rc_uo
Pointer to the segregated RhieChow interpolation object.
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
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
virtual std::pair< unsigned int, Real > correctVelocity(const bool subtract_updated_pressure, const bool recompute_face_mass_flux, const SolverParams &solver_params)
Computes new velocity field based on computed pressure gradients.
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.
const Real _pressure_variable_relaxation
The user-defined relaxation parameter for the pressure variable.
virtual void computeLinearSystemSys(libMesh::LinearImplicitSystem &sys, libMesh::SparseMatrix< libMesh::Number > &system_matrix, NumericVector< libMesh::Number > &rhs, const bool compute_gradients=true)
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.
const unsigned int _solid_energy_sys_number
The number of the system corresponding to the solid energy equation.
std::vector< LinearSystem * > _active_scalar_systems
Pointer(s) to the system(s) corresponding to the active scalar equation(s)
SIMPLESolverConfiguration _passive_scalar_linear_control
Options for the linear solver of the passive scalar equation(s)
void computeHbyA(const bool with_updated_pressure, const bool verbose)
Computes the inverse of the diagonal (1/A) of the system matrix plus the H/A components for the press...
virtual void execute(const ExecFlagType &exec_type)
Moose::PetscSupport::PetscOptions _energy_petsc_options
Options which hold the petsc settings for the fluid energy equation.
LinearSystem & _pressure_system
Reference to the nonlinear system corresponding to the pressure equation.
std::vector< unsigned int > _momentum_system_numbers
The number(s) of the system(s) corresponding to the momentum equation(s)
SIMPLESolverConfiguration _pressure_linear_control
Options for the linear solver of the pressure equation.
MultiMooseEnum getCommonPetscFlags()
LinearSystem * _energy_system
Pointer to the nonlinear system corresponding to the fluid energy equation.
const unsigned int _energy_sys_number
The number of the system corresponding to the energy equation.
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
void linkMomentumPressureSystems(const std::vector< LinearSystem *> &momentum_systems, const LinearSystem &pressure_system, const std::vector< unsigned int > &momentum_system_numbers)
Update the momentum system-related information.
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.
void computeCellVelocity()
Update the cell values of the velocity variables.
Moose::SolveType _type
virtual bool solve() override
Performs the momentum pressure coupling.
virtual void print(std::ostream &os=libMesh::out) const
const bool _pin_pressure
If the pressure needs to be pinned.
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)
void petscSetOptions(const PetscOptions &po, const SolverParams &solver_params, FEProblemBase *const problem=nullptr)
virtual void linkRhieChowUserObject() override
Fetch the Rhie Chow user object that is reponsible for determining face velocities and mass flux...
LinearSystem & getLinearSystem(unsigned int sys_num)
void set_solver_configuration(SolverConfiguration &solver_configuration)
const std::vector< SolverSystemName > & _momentum_system_names
The names of the momentum systems.
const ExecFlagType EXEC_NONLINEAR
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
const std::vector< SolverSystemName > & _active_scalar_system_names
The names of the active scalar systems.
MooseApp & _app
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
void setCurrentLinearSystem(unsigned int sys_num)
std::vector< LinearSystem * > _systems_to_solve
Shortcut to every linear system that we solve for here.
const Real _pressure_pin_value
The value we want to enforce for pressure.
virtual std::vector< std::pair< unsigned int, Real > > solveMomentumPredictor() override
Solve a momentum predictor step with a fixed pressure field.
virtual MooseMesh & mesh() override
unsigned int linearSysNum(const LinearSystemName &linear_sys_name) const override
virtual std::pair< unsigned int, Real > solvePressureCorrector() override
Solve a pressure corrector step.
Moose::PetscSupport::PetscOptions _active_scalar_petsc_options
Options which hold the petsc settings for the active scalar equation(s)
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.
const std::vector< Real > _active_scalar_absolute_tolerance
The user-defined absolute tolerance for determining the convergence in active scalars.
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.
std::vector< unsigned int > _active_scalar_system_numbers
SIMPLESolverConfiguration _solid_energy_linear_control
Options for the linear solver of the energy equation.
std::vector< LinearSystem * > _momentum_systems
Pointer(s) to the system(s) corresponding to the momentum equation(s)
const bool _continue_on_max_its
If solve should continue if maximum number of iterations is hit.
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...
bool isRecovering() const
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})$...
auto index_range(const T &sizable)
void print(std::ostream &os=libMesh::out, const bool sparse=false) const
std::vector< unsigned int > _passive_scalar_system_numbers
void initCouplingField()
Initialize the coupling fields (HbyA and Ainv)
Moose::PetscSupport::PetscOptions _pressure_petsc_options
Options which hold the petsc settings for the pressure equation.
virtual System & system() override
const bool _has_solid_energy_system
Boolean for easy check if a solid energy system shall be solved or not.
void initFaceMassFlux()
Initialize the container for face velocities.
std::pair< unsigned int, Real > solveAdvectedSystem(const unsigned int system_num, LinearSystem &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...
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)
const unsigned int _pressure_sys_number
The number of the system corresponding to the pressure equation.