https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
LinearAssemblySegregatedSolve Class Reference

Common base class for segregated solvers for the Navier-Stokes equations with linear FV assembly routines. More...

#include <LinearAssemblySegregatedSolve.h>

Inheritance diagram for LinearAssemblySegregatedSolve:
[legend]

Public Types

typedef DataFileName DataFileParameterType
 

Public Member Functions

 LinearAssemblySegregatedSolve (Executioner &ex)
 
virtual void linkRhieChowUserObject () override
 Fetch the Rhie Chow user object that is reponsible for determining face velocities and mass flux. More...
 
virtual bool solve () override
 Performs the momentum pressure coupling. More...
 
const std::vector< LinearSystem * > systemsToSolve () const
 Return pointers to the systems which are solved for within this object. More...
 
virtual void setInnerSolve (SolveObject &) override
 
void setupPressurePin ()
 Setup pressure pin if there is need for one. More...
 
virtual void checkIntegrity ()
 Check if the user defined time kernels. More...
 
virtual void initialSetup ()
 
virtual bool enabled () const
 
std::shared_ptr< MooseObjectgetSharedPtr ()
 
std::shared_ptr< const MooseObjectgetSharedPtr () const
 
MooseAppgetMooseApp () const
 
const std::string & type () const
 
virtual const std::string & name () const
 
std::string typeAndName () const
 
std::string errorPrefix (const std::string &error_type) const
 
void callMooseError (std::string msg, const bool with_prefix) const
 
MooseObjectParameterName uniqueParameterName (const std::string &parameter_name) const
 
const InputParametersparameters () const
 
MooseObjectName uniqueName () const
 
const T & getParam (const std::string &name) const
 
std::vector< std::pair< T1, T2 > > getParam (const std::string &param1, const std::string &param2) const
 
const T * queryParam (const std::string &name) const
 
const T & getRenamedParam (const std::string &old_name, const std::string &new_name) const
 
getCheckedPointerParam (const std::string &name, const std::string &error_string="") const
 
bool isParamValid (const std::string &name) const
 
bool isParamSetByUser (const std::string &nm) const
 
void paramError (const std::string &param, Args... args) const
 
void paramWarning (const std::string &param, Args... args) const
 
void paramInfo (const std::string &param, Args... args) const
 
void connectControllableParams (const std::string &parameter, const std::string &object_type, const std::string &object_name, const std::string &object_parameter) const
 
void mooseError (Args &&... args) const
 
void mooseErrorNonPrefixed (Args &&... args) const
 
void mooseDocumentedError (const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
 
void mooseWarning (Args &&... args) const
 
void mooseWarningNonPrefixed (Args &&... args) const
 
void mooseDeprecated (Args &&... args) const
 
void mooseInfo (Args &&... args) const
 
std::string getDataFileName (const std::string &param) const
 
std::string getDataFileNameByName (const std::string &relative_path) const
 
std::string getDataFilePath (const std::string &relative_path) const
 
PerfGraphperfGraph ()
 
bool isDefaultPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
bool hasPostprocessor (const std::string &param_name, const unsigned int index=0) const
 
bool hasPostprocessorByName (const PostprocessorName &name) const
 
std::size_t coupledPostprocessors (const std::string &param_name) const
 
const PostprocessorName & getPostprocessorName (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOld (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOld (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOlder (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOlder (const std::string &param_name, const unsigned int index=0) const
 
virtual const PostprocessorValuegetPostprocessorValueByName (const PostprocessorName &name) const
 
virtual const PostprocessorValuegetPostprocessorValueByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOldByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOldByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOlderByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOlderByName (const PostprocessorName &name) const
 
const Parallel::Communicator & comm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 
UserObjectName getUserObjectName (const std::string &param_name) const
 
const T & getUserObject (const std::string &param_name, bool is_dependency=true) const
 
const T & getUserObjectByName (const UserObjectName &object_name, bool is_dependency=true) const
 
const UserObjectgetUserObjectBase (const std::string &param_name, bool is_dependency=true) const
 
const UserObjectgetUserObjectBaseByName (const UserObjectName &object_name, bool is_dependency=true) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 

Static Public Member Functions

static InputParameters validParams ()
 

Public Attributes

const ConsoleStream _console
 

Protected Member Functions

virtual std::vector< std::pair< unsigned int, Real > > solveMomentumPredictor () override
 Solve a momentum predictor step with a fixed pressure field. More...
 
virtual std::pair< unsigned int, RealsolvePressureCorrector () override
 Solve a pressure corrector step. More...
 
virtual std::pair< unsigned int, RealcorrectVelocity (const bool subtract_updated_pressure, const bool recompute_face_mass_flux, const SolverParams &solver_params)
 Computes new velocity field based on computed pressure gradients. More...
 
std::pair< unsigned int, RealsolveAdvectedSystem (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 Navier-Stokes equations. More...
 
std::pair< unsigned int, RealsolveSolidEnergy ()
 Solve an equation which contains the solid energy conservation. More...
 
void checkDependentParameterError (const std::string &main_parameter, const std::vector< std::string > &dependent_parameters, const bool should_be_defined)
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level) const
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level, const std::string &live_message, const bool print_dots=true) const
 
std::string timedSectionName (const std::string &section_name) const
 
virtual void addPostprocessorDependencyHelper (const PostprocessorName &) const
 
virtual void addUserObjectDependencyHelper (const UserObject &) const
 

Protected Attributes

std::vector< unsigned int_momentum_system_numbers
 The number(s) of the system(s) corresponding to the momentum equation(s) More...
 
std::vector< LinearSystem * > _momentum_systems
 Pointer(s) to the system(s) corresponding to the momentum equation(s) More...
 
const unsigned int _pressure_sys_number
 The number of the system corresponding to the pressure equation. More...
 
LinearSystem_pressure_system
 Reference to the nonlinear system corresponding to the pressure equation. More...
 
const unsigned int _energy_sys_number
 The number of the system corresponding to the energy equation. More...
 
LinearSystem_energy_system
 Pointer to the nonlinear system corresponding to the fluid energy equation. More...
 
const unsigned int _solid_energy_sys_number
 The number of the system corresponding to the solid energy equation. More...
 
LinearSystem_solid_energy_system
 Pointer to the nonlinear system corresponding to the solid energy equation. More...
 
std::vector< LinearSystem * > _passive_scalar_systems
 Pointer(s) to the system(s) corresponding to the passive scalar equation(s) More...
 
std::vector< LinearSystem * > _active_scalar_systems
 Pointer(s) to the system(s) corresponding to the active scalar equation(s) More...
 
RhieChowMassFlux_rc_uo
 Pointer to the segregated RhieChow interpolation object. More...
 
std::vector< LinearSystem * > _systems_to_solve
 Shortcut to every linear system that we solve for here. More...
 
const std::vector< SolverSystemName > & _active_scalar_system_names
 The names of the active scalar systems. More...
 
const bool _has_active_scalar_systems
 Boolean for easy check if a active scalar systems shall be solved or not. More...
 
std::vector< unsigned int_active_scalar_system_numbers
 
const std::vector< Real_active_scalar_equation_relaxation
 The user-defined relaxation parameter(s) for the active scalar equation(s) More...
 
Moose::PetscSupport::PetscOptions _active_scalar_petsc_options
 Options which hold the petsc settings for the active scalar equation(s) More...
 
SIMPLESolverConfiguration _active_scalar_linear_control
 Options for the linear solver of the active scalar equation(s) More...
 
const Real _active_scalar_l_abs_tol
 Absolute linear tolerance for the active scalar equation(s). More...
 
const std::vector< Real_active_scalar_absolute_tolerance
 The user-defined absolute tolerance for determining the convergence in active scalars. More...
 
const std::vector< SolverSystemName > & _momentum_system_names
 The names of the momentum systems. More...
 
SIMPLESolverConfiguration _momentum_linear_control
 Options for the linear solver of the momentum equation. More...
 
const Real _momentum_l_abs_tol
 Absolute linear tolerance for the momentum equation(s). More...
 
Moose::PetscSupport::PetscOptions _momentum_petsc_options
 Options which hold the petsc settings for the momentum equation. More...
 
const Real _momentum_equation_relaxation
 The user-defined relaxation parameter for the momentum equation. More...
 
const SolverSystemName & _pressure_system_name
 The name of the pressure system. More...
 
SIMPLESolverConfiguration _pressure_linear_control
 Options for the linear solver of the pressure equation. More...
 
const Real _pressure_l_abs_tol
 Absolute linear tolerance for the pressure equation. More...
 
Moose::PetscSupport::PetscOptions _pressure_petsc_options
 Options which hold the petsc settings for the pressure equation. More...
 
const Real _pressure_variable_relaxation
 The user-defined relaxation parameter for the pressure variable. More...
 
const bool _pin_pressure
 If the pressure needs to be pinned. More...
 
const Real _pressure_pin_value
 The value we want to enforce for pressure. More...
 
dof_id_type _pressure_pin_dof
 The dof ID where the pressure needs to be pinned. More...
 
const bool _has_energy_system
 Boolean for easy check if a fluid energy system shall be solved or not. More...
 
const Real _energy_equation_relaxation
 The user-defined relaxation parameter for the energy equation. More...
 
Moose::PetscSupport::PetscOptions _energy_petsc_options
 Options which hold the petsc settings for the fluid energy equation. More...
 
SIMPLESolverConfiguration _energy_linear_control
 Options for the linear solver of the energy equation. More...
 
const Real _energy_l_abs_tol
 Absolute linear tolerance for the energy equations. More...
 
const bool _has_solid_energy_system
 Boolean for easy check if a solid energy system shall be solved or not. More...
 
Moose::PetscSupport::PetscOptions _solid_energy_petsc_options
 Options which hold the petsc settings for the fluid energy equation. More...
 
SIMPLESolverConfiguration _solid_energy_linear_control
 Options for the linear solver of the energy equation. More...
 
const Real _solid_energy_l_abs_tol
 Absolute linear tolerance for the energy equations. More...
 
const std::vector< SolverSystemName > & _passive_scalar_system_names
 The names of the passive scalar systems. More...
 
const bool _has_passive_scalar_systems
 Boolean for easy check if a passive scalar systems shall be solved or not. More...
 
std::vector< unsigned int_passive_scalar_system_numbers
 
const std::vector< Real_passive_scalar_equation_relaxation
 The user-defined relaxation parameter(s) for the passive scalar equation(s) More...
 
Moose::PetscSupport::PetscOptions _passive_scalar_petsc_options
 Options which hold the petsc settings for the passive scalar equation(s) More...
 
SIMPLESolverConfiguration _passive_scalar_linear_control
 Options for the linear solver of the passive scalar equation(s) More...
 
const Real _passive_scalar_l_abs_tol
 Absolute linear tolerance for the passive scalar equation(s). More...
 
const Real _momentum_absolute_tolerance
 The user-defined absolute tolerance for determining the convergence in momentum. More...
 
const Real _pressure_absolute_tolerance
 The user-defined absolute tolerance for determining the convergence in pressure. More...
 
const Real _energy_absolute_tolerance
 The user-defined absolute tolerance for determining the convergence in energy. More...
 
const Real _solid_energy_absolute_tolerance
 The user-defined absolute tolerance for determining the convergence in solid energy. More...
 
const std::vector< Real_passive_scalar_absolute_tolerance
 The user-defined absolute tolerance for determining the convergence in passive scalars. More...
 
const unsigned int _num_iterations
 The maximum number of momentum-pressure iterations. More...
 
const bool _continue_on_max_its
 If solve should continue if maximum number of iterations is hit. More...
 
const bool _print_fields
 Debug parameter which allows printing the coupling and solution vectors/matrices. More...
 
Executioner_executioner
 
FEProblemBase_problem
 
DisplacedProblem_displaced_problem
 
MooseMesh_mesh
 
MooseMesh_displaced_mesh
 
SystemBase_solver_sys
 
AuxiliarySystem_aux
 
SolveObject_inner_solve
 
const bool & _enabled
 
MooseApp_app
 
const std::string _type
 
const std::string _name
 
const InputParameters_pars
 
Factory_factory
 
ActionFactory_action_factory
 
MooseApp_pg_moose_app
 
const std::string _prefix
 
const Parallel::Communicator & _communicator
 

Detailed Description

Common base class for segregated solvers for the Navier-Stokes equations with linear FV assembly routines.

Once the nonlinear assembly-based routines are retired, this will be the primary base class instead of SIMPLESolveBase.

Definition at line 22 of file LinearAssemblySegregatedSolve.h.

Constructor & Destructor Documentation

◆ LinearAssemblySegregatedSolve()

LinearAssemblySegregatedSolve::LinearAssemblySegregatedSolve ( Executioner ex)

Definition at line 75 of file LinearAssemblySegregatedSolve.C.

76  : SIMPLESolveBase(ex),
77  _pressure_sys_number(_problem.linearSysNum(getParam<SolverSystemName>("pressure_system"))),
80  ? _problem.linearSysNum(getParam<SolverSystemName>("energy_system"))
85  ? _problem.linearSysNum(getParam<SolverSystemName>("solid_energy_system"))
89  _active_scalar_system_names(getParam<std::vector<SolverSystemName>>("active_scalar_systems")),
92  getParam<std::vector<Real>>("active_scalar_equation_relaxation")),
93  _active_scalar_l_abs_tol(getParam<Real>("active_scalar_l_abs_tol")),
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 }
const std::vector< Real > _active_scalar_equation_relaxation
The user-defined relaxation parameter(s) for the active scalar equation(s)
FEProblemBase & _problem
const unsigned int invalid_uint
const bool _has_energy_system
Boolean for easy check if a fluid energy system shall be solved or not.
SIMPLESolverConfiguration _active_scalar_linear_control
Options for the linear solver of the active scalar equation(s)
std::vector< LinearSystem * > _passive_scalar_systems
Pointer(s) to the system(s) corresponding to the passive scalar equation(s)
LinearSystem * _solid_energy_system
Pointer to the nonlinear system corresponding to the solid energy equation.
const bool _has_active_scalar_systems
Boolean for easy check if a active scalar systems shall be solved or not.
const Real _active_scalar_l_abs_tol
Absolute linear tolerance for the active scalar equation(s).
std::map< std::string, Real > real_valued_data
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)
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)
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.
virtual unsigned int dimension() const
std::map< std::string, int > int_valued_data
const T & getParam(const std::string &name) const
void paramError(const std::string &param, Args... args) const
LinearSystem & getLinearSystem(unsigned int sys_num)
const std::vector< SolverSystemName > & _momentum_system_names
The names of the momentum systems.
SIMPLESolveBase(Executioner &ex)
const std::vector< SolverSystemName > & _active_scalar_system_names
The names of the active scalar systems.
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)
std::vector< LinearSystem * > _systems_to_solve
Shortcut to every linear system that we solve for here.
virtual MooseMesh & mesh() override
unsigned int linearSysNum(const LinearSystemName &linear_sys_name) const override
Moose::PetscSupport::PetscOptions _active_scalar_petsc_options
Options which hold the petsc settings for the active scalar equation(s)
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)
std::vector< unsigned int > _active_scalar_system_numbers
std::vector< LinearSystem * > _momentum_systems
Pointer(s) to the system(s) corresponding to the momentum equation(s)
auto index_range(const T &sizable)
std::vector< unsigned int > _passive_scalar_system_numbers
const bool _has_solid_energy_system
Boolean for easy check if a solid energy system shall be solved or not.
const unsigned int _pressure_sys_number
The number of the system corresponding to the pressure equation.

Member Function Documentation

◆ checkDependentParameterError()

void SIMPLESolveBase::checkDependentParameterError ( const std::string &  main_parameter,
const std::vector< std::string > &  dependent_parameters,
const bool  should_be_defined 
)
protectedinherited

Definition at line 487 of file SIMPLESolveBase.C.

Referenced by SIMPLESolveBase::SIMPLESolveBase(), and SIMPLESolveNonlinearAssembly::SIMPLESolveNonlinearAssembly().

490 {
491  for (const auto & param : dependent_parameters)
492  if (parameters().isParamSetByUser(param) == !should_be_defined)
493  paramError(param,
494  "This parameter should " + std::string(should_be_defined ? "" : "not") +
495  " be given by the user with the corresponding " + main_parameter +
496  " setting!");
497 }
void paramError(const std::string &param, Args... args) const
bool isParamSetByUser(const std::string &nm) const
const InputParameters & parameters() const

◆ checkIntegrity()

virtual void SIMPLESolveBase::checkIntegrity ( )
inlinevirtualinherited

Check if the user defined time kernels.

Reimplemented in SIMPLESolve, and SIMPLESolveNonlinearAssembly.

Definition at line 61 of file SIMPLESolveBase.h.

61 {}

◆ correctVelocity()

std::pair< unsigned int, Real > LinearAssemblySegregatedSolve::correctVelocity ( const bool  subtract_updated_pressure,
const bool  recompute_face_mass_flux,
const SolverParams solver_params 
)
protectedvirtual

Computes new velocity field based on computed pressure gradients.

Parameters
subtract_updated_pressureIf we need to subtract the updated pressure gradient from the right hand side of the system
recompute_face_mass_fluxIf we want to recompute the face flux too
solver_paramsDummy solver parameter object for the linear solve

Reimplemented in PIMPLESolve.

Definition at line 401 of file LinearAssemblySegregatedSolve.C.

Referenced by PIMPLESolve::correctVelocity(), and solve().

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 }
void computeGradients()
const bool _print_fields
Debug parameter which allows printing the coupling and solution vectors/matrices. ...
void computeFaceMassFlux()
Update the values of the face velocities in the containers.
void setSolution(const NumericVector< Number > &soln)
RhieChowMassFlux * _rc_uo
Pointer to the segregated RhieChow interpolation object.
const Real _pressure_variable_relaxation
The user-defined relaxation parameter for the pressure variable.
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...
LinearSystem & _pressure_system
Reference to the nonlinear system corresponding to the pressure equation.
void computeCellVelocity()
Update the cell values of the velocity variables.
void petscSetOptions(const PetscOptions &po, const SolverParams &solver_params, FEProblemBase *const problem=nullptr)
virtual std::pair< unsigned int, Real > solvePressureCorrector() override
Solve a pressure corrector step.
std::unique_ptr< NumericVector< Number > > current_local_solution
virtual const NumericVector< Number > * solutionPreviousNewton() 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})$...
Moose::PetscSupport::PetscOptions _pressure_petsc_options
Options which hold the petsc settings for the pressure equation.
virtual System & system() override

◆ linkRhieChowUserObject()

void LinearAssemblySegregatedSolve::linkRhieChowUserObject ( )
overridevirtual

Fetch the Rhie Chow user object that is reponsible for determining face velocities and mass flux.

Implements SIMPLESolveBase.

Definition at line 164 of file LinearAssemblySegregatedSolve.C.

Referenced by SIMPLE::init(), and PIMPLE::init().

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 }
User object responsible for determining the face fluxes using the Rhie-Chow interpolation in a segreg...
RhieChowMassFlux * _rc_uo
Pointer to the segregated RhieChow interpolation object.
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)
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.
MooseApp & _app
std::vector< LinearSystem * > _momentum_systems
Pointer(s) to the system(s) corresponding to the momentum equation(s)
bool isRecovering() const
void initCouplingField()
Initialize the coupling fields (HbyA and Ainv)
void initFaceMassFlux()
Initialize the container for face velocities.

◆ setInnerSolve()

virtual void SIMPLESolveBase::setInnerSolve ( SolveObject )
inlineoverridevirtualinherited

Reimplemented from SolveObject.

Definition at line 48 of file SIMPLESolveBase.h.

49  {
50  mooseError("Cannot set inner solve object for solves that inherit from SIMPLESolveBase");
51  }
void mooseError(Args &&... args) const

◆ setupPressurePin()

void SIMPLESolveBase::setupPressurePin ( )
inherited

Setup pressure pin if there is need for one.

Definition at line 478 of file SIMPLESolveBase.C.

Referenced by SIMPLE::init(), SIMPLENonlinearAssembly::init(), and PIMPLE::init().

479 {
480  if (_pin_pressure)
482  _problem.mesh(),
483  getParam<Point>("pressure_pin_point"));
484 }
FEProblemBase & _problem
virtual const MooseVariableFieldBase & getVariable(const THREAD_ID tid, const std::string &var_name, Moose::VarKindType expected_var_type=Moose::VarKindType::VAR_ANY, Moose::VarFieldType expected_var_field_type=Moose::VarFieldType::VAR_FIELD_ANY) const override
dof_id_type _pressure_pin_dof
The dof ID where the pressure needs to be pinned.
const bool _pin_pressure
If the pressure needs to be pinned.
dof_id_type findPointDoFID(const MooseVariableFieldBase &variable, const MooseMesh &mesh, const Point &point)
Find the ID of the degree of freedom which corresponds to the variable and a given point on the mesh...
virtual MooseMesh & mesh() override

◆ solve()

bool LinearAssemblySegregatedSolve::solve ( )
overridevirtual

Performs the momentum pressure coupling.

Returns
True if solver is converged.

Implements SolveObject.

Definition at line 514 of file LinearAssemblySegregatedSolve.C.

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  {
597  _problem.execute(EXEC_NONLINEAR);
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
611  _problem.execute(EXEC_NONLINEAR);
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 }
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
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.
SIMPLESolverConfiguration _energy_linear_control
Options for the linear solver of the energy equation.
void computeGradients()
const bool _has_energy_system
Boolean for easy check if a fluid energy system shall be solved or not.
Moose::LineSearchType _line_search
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 bool _has_active_scalar_systems
Boolean for easy check if a active scalar systems shall be solved or not.
const Real _active_scalar_l_abs_tol
Absolute linear tolerance for the active scalar equation(s).
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 _momentum_absolute_tolerance
The user-defined absolute tolerance for determining the convergence in momentum.
const std::vector< SolverSystemName > & _passive_scalar_system_names
The names of the passive scalar systems.
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)
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.
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)
Moose::PetscSupport::PetscOptions _momentum_petsc_options
Options which hold the petsc settings for the momentum equation.
const std::vector< Real > _passive_scalar_absolute_tolerance
The user-defined absolute tolerance for determining the convergence in passive scalars.
Moose::SolveType _type
const Real _energy_l_abs_tol
Absolute linear tolerance for the energy equations.
void petscSetOptions(const PetscOptions &po, const SolverParams &solver_params, FEProblemBase *const problem=nullptr)
const std::vector< Real > _passive_scalar_equation_relaxation
The user-defined relaxation parameter(s) for the passive scalar equation(s)
const std::vector< SolverSystemName > & _active_scalar_system_names
The names of the active scalar systems.
const Real _energy_equation_relaxation
The user-defined relaxation parameter for the energy equation.
virtual std::vector< std::pair< unsigned int, Real > > solveMomentumPredictor() override
Solve a momentum predictor step with a fixed pressure field.
Moose::PetscSupport::PetscOptions _active_scalar_petsc_options
Options which hold the petsc settings for the active scalar equation(s)
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.
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
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.
auto index_range(const T &sizable)
std::vector< unsigned int > _passive_scalar_system_numbers
const bool _has_solid_energy_system
Boolean for easy check if a solid energy system shall be solved or not.
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...

◆ solveAdvectedSystem()

std::pair< unsigned int, Real > LinearAssemblySegregatedSolve::solveAdvectedSystem ( const unsigned int  system_num,
LinearSystem system,
const Real  relaxation_factor,
libMesh::SolverConfiguration solver_config,
const Real  abs_tol 
)
protected

Solve an equation which contains an advection term that depends on the solution of the segregated Navier-Stokes equations.

Parameters
system_numThe number of the system which is solved
systemReference to the system which is solved
relaxation_factorThe relaxation factor for matrix relaxation
solver_configThe solver configuration object for the linear solve
abs_tolThe scaled absolute tolerance for the linear solve
Returns
The normalized residual norm of the equation.

Definition at line 443 of file LinearAssemblySegregatedSolve.C.

Referenced by solve().

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 }
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
const bool _print_fields
Debug parameter which allows printing the coupling and solution vectors/matrices. ...
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...
virtual std::unique_ptr< NumericVector< T > > zero_clone() const=0
NumericVector< Number > * rhs
void setSolution(const NumericVector< Number > &soln)
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
virtual LinearSolver< Number > * get_linear_solver() const override
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 void computeLinearSystemSys(libMesh::LinearImplicitSystem &sys, libMesh::SparseMatrix< libMesh::Number > &system_matrix, NumericVector< libMesh::Number > &rhs, const bool compute_gradients=true)
virtual const std::string & name() const
std::unique_ptr< NumericVector< Number > > solution
virtual void print(std::ostream &os=libMesh::out) const
void set_solver_configuration(SolverConfiguration &solver_configuration)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
SparseMatrix< Number > * matrix
void setCurrentLinearSystem(unsigned int sys_num)
std::unique_ptr< NumericVector< Number > > current_local_solution
const ConsoleStream _console
void print(std::ostream &os=libMesh::out, const bool sparse=false) const
virtual System & system() override

◆ solveMomentumPredictor()

std::vector< std::pair< unsigned int, Real > > LinearAssemblySegregatedSolve::solveMomentumPredictor ( )
overrideprotectedvirtual

Solve a momentum predictor step with a fixed pressure field.

Returns
A vector of (number of linear iterations, normalized residual norm) pairs for the momentum equations. The length of the vector equals the dimensionality of the domain.

Implements SIMPLESolveBase.

Definition at line 178 of file LinearAssemblySegregatedSolve.C.

Referenced by solve().

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 }
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 bool _print_fields
Debug parameter which allows printing the coupling and solution vectors/matrices. ...
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 _momentum_equation_relaxation
The user-defined relaxation parameter for the momentum equation.
NumericVector< Number > * rhs
const Real _momentum_l_abs_tol
Absolute linear tolerance for the momentum equation(s).
virtual LinearSolver< Number > * get_linear_solver() const override
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 void computeLinearSystemSys(libMesh::LinearImplicitSystem &sys, libMesh::SparseMatrix< libMesh::Number > &system_matrix, NumericVector< libMesh::Number > &rhs, const bool compute_gradients=true)
std::vector< unsigned int > _momentum_system_numbers
The number(s) of the system(s) corresponding to the momentum equation(s)
virtual std::unique_ptr< SparseMatrix< T > > zero_clone() const=0
std::unique_ptr< NumericVector< Number > > solution
virtual void print(std::ostream &os=libMesh::out) const
std::string stringify(const T &t)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
SparseMatrix< Number > * matrix
void setCurrentLinearSystem(unsigned int sys_num)
std::unique_ptr< NumericVector< Number > > current_local_solution
const ConsoleStream _console
std::vector< LinearSystem * > _momentum_systems
Pointer(s) to the system(s) corresponding to the momentum equation(s)
auto index_range(const T &sizable)
void print(std::ostream &os=libMesh::out, const bool sparse=false) const

◆ solvePressureCorrector()

std::pair< unsigned int, Real > LinearAssemblySegregatedSolve::solvePressureCorrector ( )
overrideprotectedvirtual

Solve a pressure corrector step.

Returns
The number of linear iterations and the normalized residual norm of the pressure equation.

Implements SIMPLESolveBase.

Definition at line 277 of file LinearAssemblySegregatedSolve.C.

Referenced by correctVelocity().

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 }
FEProblemBase & _problem
const bool _print_fields
Debug parameter which allows printing the coupling and solution vectors/matrices. ...
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...
NumericVector< Number > * rhs
void setSolution(const NumericVector< Number > &soln)
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
virtual LinearSolver< Number > * get_linear_solver() const override
std::map< std::string, Real > real_valued_data
virtual void computeLinearSystemSys(libMesh::LinearImplicitSystem &sys, libMesh::SparseMatrix< libMesh::Number > &system_matrix, NumericVector< libMesh::Number > &rhs, const bool compute_gradients=true)
LinearSystem & _pressure_system
Reference to the nonlinear system corresponding to the pressure equation.
SIMPLESolverConfiguration _pressure_linear_control
Options for the linear solver of the pressure equation.
const Real _pressure_l_abs_tol
Absolute linear tolerance for the pressure equation.
std::unique_ptr< NumericVector< Number > > solution
dof_id_type _pressure_pin_dof
The dof ID where the pressure needs to be pinned.
virtual void print(std::ostream &os=libMesh::out) const
const bool _pin_pressure
If the pressure needs to be pinned.
void set_solver_configuration(SolverConfiguration &solver_configuration)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
SparseMatrix< Number > * matrix
void setCurrentLinearSystem(unsigned int sys_num)
const Real _pressure_pin_value
The value we want to enforce for pressure.
std::unique_ptr< NumericVector< Number > > current_local_solution
const ConsoleStream _console
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 print(std::ostream &os=libMesh::out, const bool sparse=false) const
virtual System & system() override
const unsigned int _pressure_sys_number
The number of the system corresponding to the pressure equation.

◆ solveSolidEnergy()

std::pair< unsigned int, Real > LinearAssemblySegregatedSolve::solveSolidEnergy ( )
protected

Solve an equation which contains the solid energy conservation.

Definition at line 341 of file LinearAssemblySegregatedSolve.C.

Referenced by solve().

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 }
FEProblemBase & _problem
const bool _print_fields
Debug parameter which allows printing the coupling and solution vectors/matrices. ...
const Real _solid_energy_l_abs_tol
Absolute linear tolerance for the energy equations.
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...
LinearSystem * _solid_energy_system
Pointer to the nonlinear system corresponding to the solid energy equation.
NumericVector< Number > * rhs
void setSolution(const NumericVector< Number > &soln)
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
virtual LinearSolver< Number > * get_linear_solver() const override
std::map< std::string, Real > real_valued_data
virtual void computeLinearSystemSys(libMesh::LinearImplicitSystem &sys, libMesh::SparseMatrix< libMesh::Number > &system_matrix, NumericVector< libMesh::Number > &rhs, const bool compute_gradients=true)
const unsigned int _solid_energy_sys_number
The number of the system corresponding to the solid energy equation.
std::unique_ptr< NumericVector< Number > > solution
virtual void print(std::ostream &os=libMesh::out) const
void set_solver_configuration(SolverConfiguration &solver_configuration)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
SparseMatrix< Number > * matrix
void setCurrentLinearSystem(unsigned int sys_num)
std::unique_ptr< NumericVector< Number > > current_local_solution
const ConsoleStream _console
SIMPLESolverConfiguration _solid_energy_linear_control
Options for the linear solver of the energy equation.
void print(std::ostream &os=libMesh::out, const bool sparse=false) const
virtual System & system() override

◆ systemsToSolve()

const std::vector<LinearSystem *> LinearAssemblySegregatedSolve::systemsToSolve ( ) const
inline

Return pointers to the systems which are solved for within this object.

Definition at line 38 of file LinearAssemblySegregatedSolve.h.

Referenced by PIMPLE::getTimeIntegrators(), and PIMPLE::relativeSolutionDifferenceNorm().

38 { return _systems_to_solve; }
std::vector< LinearSystem * > _systems_to_solve
Shortcut to every linear system that we solve for here.

◆ validParams()

InputParameters LinearAssemblySegregatedSolve::validParams ( )
static

Definition at line 18 of file LinearAssemblySegregatedSolve.C.

Referenced by SIMPLESolve::validParams(), and PIMPLESolve::validParams().

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 }
MultiMooseEnum getCommonPetscKeys()
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
static InputParameters validParams()
MultiMooseEnum getCommonPetscFlags()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

Member Data Documentation

◆ _active_scalar_absolute_tolerance

const std::vector<Real> LinearAssemblySegregatedSolve::_active_scalar_absolute_tolerance
protected

The user-defined absolute tolerance for determining the convergence in active scalars.

Definition at line 131 of file LinearAssemblySegregatedSolve.h.

Referenced by solve().

◆ _active_scalar_equation_relaxation

const std::vector<Real> LinearAssemblySegregatedSolve::_active_scalar_equation_relaxation
protected

The user-defined relaxation parameter(s) for the active scalar equation(s)

Definition at line 118 of file LinearAssemblySegregatedSolve.h.

Referenced by LinearAssemblySegregatedSolve(), and solve().

◆ _active_scalar_l_abs_tol

const Real LinearAssemblySegregatedSolve::_active_scalar_l_abs_tol
protected

Absolute linear tolerance for the active scalar equation(s).

We need to store this, because it needs to be scaled with a representative flux.

Definition at line 128 of file LinearAssemblySegregatedSolve.h.

Referenced by solve().

◆ _active_scalar_linear_control

SIMPLESolverConfiguration LinearAssemblySegregatedSolve::_active_scalar_linear_control
protected

Options for the linear solver of the active scalar equation(s)

Definition at line 124 of file LinearAssemblySegregatedSolve.h.

Referenced by LinearAssemblySegregatedSolve(), and solve().

◆ _active_scalar_petsc_options

Moose::PetscSupport::PetscOptions LinearAssemblySegregatedSolve::_active_scalar_petsc_options
protected

Options which hold the petsc settings for the active scalar equation(s)

Definition at line 121 of file LinearAssemblySegregatedSolve.h.

Referenced by LinearAssemblySegregatedSolve(), and solve().

◆ _active_scalar_system_names

const std::vector<SolverSystemName>& LinearAssemblySegregatedSolve::_active_scalar_system_names
protected

The names of the active scalar systems.

Definition at line 109 of file LinearAssemblySegregatedSolve.h.

Referenced by LinearAssemblySegregatedSolve(), and solve().

◆ _active_scalar_system_numbers

std::vector<unsigned int> LinearAssemblySegregatedSolve::_active_scalar_system_numbers
protected

Definition at line 115 of file LinearAssemblySegregatedSolve.h.

Referenced by LinearAssemblySegregatedSolve(), and solve().

◆ _active_scalar_systems

std::vector<LinearSystem *> LinearAssemblySegregatedSolve::_active_scalar_systems
protected

Pointer(s) to the system(s) corresponding to the active scalar equation(s)

Definition at line 98 of file LinearAssemblySegregatedSolve.h.

Referenced by SIMPLESolve::checkIntegrity(), LinearAssemblySegregatedSolve(), and solve().

◆ _continue_on_max_its

const bool SIMPLESolveBase::_continue_on_max_its
protectedinherited

If solve should continue if maximum number of iterations is hit.

Definition at line 202 of file SIMPLESolveBase.h.

Referenced by solve(), and SIMPLESolveNonlinearAssembly::solve().

◆ _energy_absolute_tolerance

const Real SIMPLESolveBase::_energy_absolute_tolerance
protectedinherited

The user-defined absolute tolerance for determining the convergence in energy.

Definition at line 190 of file SIMPLESolveBase.h.

Referenced by solve(), and SIMPLESolveNonlinearAssembly::solve().

◆ _energy_equation_relaxation

const Real SIMPLESolveBase::_energy_equation_relaxation
protectedinherited

The user-defined relaxation parameter for the energy equation.

Definition at line 130 of file SIMPLESolveBase.h.

Referenced by solve(), and SIMPLESolveNonlinearAssembly::solve().

◆ _energy_l_abs_tol

const Real SIMPLESolveBase::_energy_l_abs_tol
protectedinherited

Absolute linear tolerance for the energy equations.

We need to store this, because it needs to be scaled with a representative flux.

Definition at line 140 of file SIMPLESolveBase.h.

Referenced by solve(), and SIMPLESolveNonlinearAssembly::solve().

◆ _energy_linear_control

SIMPLESolverConfiguration SIMPLESolveBase::_energy_linear_control
protectedinherited

Options for the linear solver of the energy equation.

Definition at line 136 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::SIMPLESolveBase(), SIMPLESolveNonlinearAssembly::solve(), and solve().

◆ _energy_petsc_options

Moose::PetscSupport::PetscOptions SIMPLESolveBase::_energy_petsc_options
protectedinherited

Options which hold the petsc settings for the fluid energy equation.

Definition at line 133 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::SIMPLESolveBase(), SIMPLESolveNonlinearAssembly::solve(), and solve().

◆ _energy_sys_number

const unsigned int LinearAssemblySegregatedSolve::_energy_sys_number
protected

The number of the system corresponding to the energy equation.

Definition at line 83 of file LinearAssemblySegregatedSolve.h.

Referenced by solve().

◆ _energy_system

LinearSystem* LinearAssemblySegregatedSolve::_energy_system
protected

Pointer to the nonlinear system corresponding to the fluid energy equation.

Definition at line 86 of file LinearAssemblySegregatedSolve.h.

Referenced by SIMPLESolve::checkIntegrity(), LinearAssemblySegregatedSolve(), and solve().

◆ _has_active_scalar_systems

const bool LinearAssemblySegregatedSolve::_has_active_scalar_systems
protected

Boolean for easy check if a active scalar systems shall be solved or not.

Definition at line 112 of file LinearAssemblySegregatedSolve.h.

Referenced by SIMPLESolve::checkIntegrity(), LinearAssemblySegregatedSolve(), and solve().

◆ _has_energy_system

const bool SIMPLESolveBase::_has_energy_system
protectedinherited

◆ _has_passive_scalar_systems

const bool SIMPLESolveBase::_has_passive_scalar_systems
protectedinherited

◆ _has_solid_energy_system

const bool SIMPLESolveBase::_has_solid_energy_system
protectedinherited

Boolean for easy check if a solid energy system shall be solved or not.

Definition at line 145 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveNonlinearAssembly::checkIntegrity(), LinearAssemblySegregatedSolve(), SIMPLESolveBase::SIMPLESolveBase(), SIMPLESolveNonlinearAssembly::solve(), and solve().

◆ _momentum_absolute_tolerance

const Real SIMPLESolveBase::_momentum_absolute_tolerance
protectedinherited

The user-defined absolute tolerance for determining the convergence in momentum.

Definition at line 184 of file SIMPLESolveBase.h.

Referenced by solve(), and SIMPLESolveNonlinearAssembly::solve().

◆ _momentum_equation_relaxation

const Real SIMPLESolveBase::_momentum_equation_relaxation
protectedinherited

The user-defined relaxation parameter for the momentum equation.

Definition at line 95 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveNonlinearAssembly::solveMomentumPredictor(), and solveMomentumPredictor().

◆ _momentum_l_abs_tol

const Real SIMPLESolveBase::_momentum_l_abs_tol
protectedinherited

Absolute linear tolerance for the momentum equation(s).

We need to store this, because it needs to be scaled with a representative flux.

Definition at line 89 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveNonlinearAssembly::solveMomentumPredictor(), and solveMomentumPredictor().

◆ _momentum_linear_control

SIMPLESolverConfiguration SIMPLESolveBase::_momentum_linear_control
protectedinherited

Options for the linear solver of the momentum equation.

Definition at line 85 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::SIMPLESolveBase(), SIMPLESolveNonlinearAssembly::solveMomentumPredictor(), and solveMomentumPredictor().

◆ _momentum_petsc_options

Moose::PetscSupport::PetscOptions SIMPLESolveBase::_momentum_petsc_options
protectedinherited

Options which hold the petsc settings for the momentum equation.

Definition at line 92 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::SIMPLESolveBase(), SIMPLESolveNonlinearAssembly::solve(), and solve().

◆ _momentum_system_names

const std::vector<SolverSystemName>& SIMPLESolveBase::_momentum_system_names
protectedinherited

◆ _momentum_system_numbers

std::vector<unsigned int> LinearAssemblySegregatedSolve::_momentum_system_numbers
protected

The number(s) of the system(s) corresponding to the momentum equation(s)

Definition at line 71 of file LinearAssemblySegregatedSolve.h.

Referenced by LinearAssemblySegregatedSolve(), linkRhieChowUserObject(), and solveMomentumPredictor().

◆ _momentum_systems

std::vector<LinearSystem *> LinearAssemblySegregatedSolve::_momentum_systems
protected

Pointer(s) to the system(s) corresponding to the momentum equation(s)

Definition at line 74 of file LinearAssemblySegregatedSolve.h.

Referenced by SIMPLESolve::checkIntegrity(), LinearAssemblySegregatedSolve(), linkRhieChowUserObject(), solve(), and solveMomentumPredictor().

◆ _num_iterations

const unsigned int SIMPLESolveBase::_num_iterations
protectedinherited

The maximum number of momentum-pressure iterations.

Definition at line 199 of file SIMPLESolveBase.h.

Referenced by solve(), and SIMPLESolveNonlinearAssembly::solve().

◆ _passive_scalar_absolute_tolerance

const std::vector<Real> SIMPLESolveBase::_passive_scalar_absolute_tolerance
protectedinherited

The user-defined absolute tolerance for determining the convergence in passive scalars.

Definition at line 196 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::SIMPLESolveBase(), SIMPLESolveNonlinearAssembly::solve(), and solve().

◆ _passive_scalar_equation_relaxation

const std::vector<Real> SIMPLESolveBase::_passive_scalar_equation_relaxation
protectedinherited

The user-defined relaxation parameter(s) for the passive scalar equation(s)

Definition at line 169 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::SIMPLESolveBase(), SIMPLESolveNonlinearAssembly::solve(), and solve().

◆ _passive_scalar_l_abs_tol

const Real SIMPLESolveBase::_passive_scalar_l_abs_tol
protectedinherited

Absolute linear tolerance for the passive scalar equation(s).

We need to store this, because it needs to be scaled with a representative flux.

Definition at line 179 of file SIMPLESolveBase.h.

Referenced by solve(), and SIMPLESolveNonlinearAssembly::solve().

◆ _passive_scalar_linear_control

SIMPLESolverConfiguration SIMPLESolveBase::_passive_scalar_linear_control
protectedinherited

Options for the linear solver of the passive scalar equation(s)

Definition at line 175 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::SIMPLESolveBase(), SIMPLESolveNonlinearAssembly::solve(), and solve().

◆ _passive_scalar_petsc_options

Moose::PetscSupport::PetscOptions SIMPLESolveBase::_passive_scalar_petsc_options
protectedinherited

Options which hold the petsc settings for the passive scalar equation(s)

Definition at line 172 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::SIMPLESolveBase(), SIMPLESolveNonlinearAssembly::solve(), and solve().

◆ _passive_scalar_system_names

const std::vector<SolverSystemName>& SIMPLESolveBase::_passive_scalar_system_names
protectedinherited

◆ _passive_scalar_system_numbers

std::vector<unsigned int> SIMPLESolveBase::_passive_scalar_system_numbers
protectedinherited

◆ _passive_scalar_systems

std::vector<LinearSystem *> LinearAssemblySegregatedSolve::_passive_scalar_systems
protected

Pointer(s) to the system(s) corresponding to the passive scalar equation(s)

Definition at line 95 of file LinearAssemblySegregatedSolve.h.

Referenced by SIMPLESolve::checkIntegrity(), LinearAssemblySegregatedSolve(), and solve().

◆ _pin_pressure

const bool SIMPLESolveBase::_pin_pressure
protectedinherited

If the pressure needs to be pinned.

Definition at line 116 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::setupPressurePin(), SIMPLESolveNonlinearAssembly::solvePressureCorrector(), and solvePressureCorrector().

◆ _pressure_absolute_tolerance

const Real SIMPLESolveBase::_pressure_absolute_tolerance
protectedinherited

The user-defined absolute tolerance for determining the convergence in pressure.

Definition at line 187 of file SIMPLESolveBase.h.

Referenced by solve(), and SIMPLESolveNonlinearAssembly::solve().

◆ _pressure_l_abs_tol

const Real SIMPLESolveBase::_pressure_l_abs_tol
protectedinherited

Absolute linear tolerance for the pressure equation.

We need to store this, because it needs to be scaled with a representative flux.

Definition at line 107 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveNonlinearAssembly::solvePressureCorrector(), and solvePressureCorrector().

◆ _pressure_linear_control

SIMPLESolverConfiguration SIMPLESolveBase::_pressure_linear_control
protectedinherited

Options for the linear solver of the pressure equation.

Definition at line 103 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::SIMPLESolveBase(), SIMPLESolveNonlinearAssembly::solvePressureCorrector(), and solvePressureCorrector().

◆ _pressure_petsc_options

Moose::PetscSupport::PetscOptions SIMPLESolveBase::_pressure_petsc_options
protectedinherited

Options which hold the petsc settings for the pressure equation.

Definition at line 110 of file SIMPLESolveBase.h.

Referenced by correctVelocity(), SIMPLESolveBase::SIMPLESolveBase(), and SIMPLESolveNonlinearAssembly::solve().

◆ _pressure_pin_dof

dof_id_type SIMPLESolveBase::_pressure_pin_dof
protectedinherited

The dof ID where the pressure needs to be pinned.

Definition at line 122 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::setupPressurePin(), SIMPLESolveNonlinearAssembly::solvePressureCorrector(), and solvePressureCorrector().

◆ _pressure_pin_value

const Real SIMPLESolveBase::_pressure_pin_value
protectedinherited

The value we want to enforce for pressure.

Definition at line 119 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveNonlinearAssembly::solvePressureCorrector(), and solvePressureCorrector().

◆ _pressure_sys_number

const unsigned int LinearAssemblySegregatedSolve::_pressure_sys_number
protected

The number of the system corresponding to the pressure equation.

Definition at line 77 of file LinearAssemblySegregatedSolve.h.

Referenced by solvePressureCorrector().

◆ _pressure_system

LinearSystem& LinearAssemblySegregatedSolve::_pressure_system
protected

Reference to the nonlinear system corresponding to the pressure equation.

Definition at line 80 of file LinearAssemblySegregatedSolve.h.

Referenced by SIMPLESolve::checkIntegrity(), correctVelocity(), LinearAssemblySegregatedSolve(), linkRhieChowUserObject(), solve(), and solvePressureCorrector().

◆ _pressure_system_name

const SolverSystemName& SIMPLESolveBase::_pressure_system_name
protectedinherited

The name of the pressure system.

Definition at line 100 of file SIMPLESolveBase.h.

◆ _pressure_variable_relaxation

const Real SIMPLESolveBase::_pressure_variable_relaxation
protectedinherited

The user-defined relaxation parameter for the pressure variable.

Definition at line 113 of file SIMPLESolveBase.h.

Referenced by correctVelocity(), and SIMPLESolveNonlinearAssembly::solve().

◆ _print_fields

const bool SIMPLESolveBase::_print_fields
protectedinherited

◆ _rc_uo

RhieChowMassFlux* LinearAssemblySegregatedSolve::_rc_uo
protected

Pointer to the segregated RhieChow interpolation object.

Definition at line 101 of file LinearAssemblySegregatedSolve.h.

Referenced by correctVelocity(), and linkRhieChowUserObject().

◆ _solid_energy_absolute_tolerance

const Real SIMPLESolveBase::_solid_energy_absolute_tolerance
protectedinherited

The user-defined absolute tolerance for determining the convergence in solid energy.

Definition at line 193 of file SIMPLESolveBase.h.

Referenced by solve(), and SIMPLESolveNonlinearAssembly::solve().

◆ _solid_energy_l_abs_tol

const Real SIMPLESolveBase::_solid_energy_l_abs_tol
protectedinherited

Absolute linear tolerance for the energy equations.

We need to store this, because it needs to be scaled with a representative flux.

Definition at line 155 of file SIMPLESolveBase.h.

Referenced by solveSolidEnergy(), and SIMPLESolveNonlinearAssembly::solveSolidEnergySystem().

◆ _solid_energy_linear_control

SIMPLESolverConfiguration SIMPLESolveBase::_solid_energy_linear_control
protectedinherited

Options for the linear solver of the energy equation.

Definition at line 151 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::SIMPLESolveBase(), solveSolidEnergy(), and SIMPLESolveNonlinearAssembly::solveSolidEnergySystem().

◆ _solid_energy_petsc_options

Moose::PetscSupport::PetscOptions SIMPLESolveBase::_solid_energy_petsc_options
protectedinherited

Options which hold the petsc settings for the fluid energy equation.

Definition at line 148 of file SIMPLESolveBase.h.

Referenced by SIMPLESolveBase::SIMPLESolveBase(), and solve().

◆ _solid_energy_sys_number

const unsigned int LinearAssemblySegregatedSolve::_solid_energy_sys_number
protected

The number of the system corresponding to the solid energy equation.

Definition at line 89 of file LinearAssemblySegregatedSolve.h.

Referenced by solveSolidEnergy().

◆ _solid_energy_system

LinearSystem* LinearAssemblySegregatedSolve::_solid_energy_system
protected

Pointer to the nonlinear system corresponding to the solid energy equation.

Definition at line 92 of file LinearAssemblySegregatedSolve.h.

Referenced by LinearAssemblySegregatedSolve(), and solveSolidEnergy().

◆ _systems_to_solve

std::vector<LinearSystem *> LinearAssemblySegregatedSolve::_systems_to_solve
protected

Shortcut to every linear system that we solve for here.

Definition at line 104 of file LinearAssemblySegregatedSolve.h.

Referenced by LinearAssemblySegregatedSolve(), and systemsToSolve().


The documentation for this class was generated from the following files: