PIMPLE

Solves the transient Navier-Stokes equations using the PIMPLE algorithm and linear finite volume variables.

Overview

This transient executioner is based on the algorithm discussed in Greenshields and Weller (2022). It is a combination of SIMPLE and the PISO (Pressure-Implicit with Splitting of Operators) algorithm introduced by Issa (1986).

In PISO, the user iterates between the source term in the pressure equation and the corrected pressure and velocity fields without assembling the momentum system again. With the notation already used in SIMPLE, the PISO iteration is the following:

  1. Take a velocity guess and compute which is the offdiagonals of the momentum system matrix multiplied by the current guess minus the right hand side of the momentum matrix. Note that the face flux is computed using a different quess so this operation is just a matrix-vector multiplication and a vector-vector addition.

  2. Use and the diagonal of the momentum matrix to solve the pressure equation:

    The pressure solution might have to be relaxed in this iteration:

  3. Once the pressure solution is obtained, update the velocity to the next guess:

    (1)

    and return to (1) until the maximum number of iterations is reached which can be set using the "num_piso_iterations" parameter.

Example Input Syntax

The problem setup is exactly the same as discussed for SIMPLE, only the executioner block is different:

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  type = PIMPLE
  momentum_l_abs_tol = 1e-12
  pressure_l_abs_tol = 1e-12
  energy_l_abs_tol = 1e-12
  momentum_l_tol = 1e-12
  pressure_l_tol = 1e-12
  energy_l_tol = 1e-12
  rhie_chow_user_object = 'rc'
  momentum_systems = 'u_system v_system'
  pressure_system = 'pressure_system'
  energy_system = 'energy_system'
  momentum_equation_relaxation = 0.8
  pressure_variable_relaxation = 0.3
  energy_equation_relaxation = 0.9
  num_iterations = 100
  pressure_absolute_tolerance = 1e-11
  momentum_absolute_tolerance = 1e-11
  energy_absolute_tolerance = 1e-11
  momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
  momentum_petsc_options_value = 'hypre boomeramg'
  pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
  pressure_petsc_options_value = 'hypre boomeramg'
  energy_petsc_options_iname = '-pc_type -pc_hypre_type'
  energy_petsc_options_value = 'hypre boomeramg'
  print_fields = false
  continue_on_max_its = true
  dt = 0.01
  num_steps = 6
  num_piso_iterations = 0
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-boussinesq-transient.i)

Input Parameters

  • rhie_chow_user_objectThe rhie-chow user-object

    C++ Type:UserObjectName

    Controllable:No

    Description:The rhie-chow user-object

Required Parameters

  • continue_on_max_itsFalseIf solve should continue if maximum number of iterations is hit.

    Default:False

    C++ Type:bool

    Controllable:No

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

  • dt1The timestep size between solves

    Default:1

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The timestep size between solves

  • end_time1e+30The end time of the simulation

    Default:1e+30

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The end time of the simulation

  • energy_systemThe solver system for the energy equation.

    C++ Type:SolverSystemName

    Controllable:No

    Description:The solver system for the energy equation.

  • error_on_dtminTrueThrow error when timestep is less than dtmin instead of just aborting solve.

    Default:True

    C++ Type:bool

    Controllable:No

    Description:Throw error when timestep is less than dtmin instead of just aborting solve.

  • normalize_solution_diff_norm_by_dtTrueWhether to divide the solution difference norm by dt. If taking 'small' time steps you probably want this to be true. If taking very 'large' timesteps in an attempt to *reach* a steady-state, you probably want this parameter to be false.

    Default:True

    C++ Type:bool

    Controllable:No

    Description:Whether to divide the solution difference norm by dt. If taking 'small' time steps you probably want this to be true. If taking very 'large' timesteps in an attempt to *reach* a steady-state, you probably want this parameter to be false.

  • num_iterations1000The number of momentum-pressure-(other fields) iterations needed.

    Default:1000

    C++ Type:unsigned int

    Controllable:No

    Description:The number of momentum-pressure-(other fields) iterations needed.

  • num_piso_iterations0The number of PISO iterations without recomputing the momentum matrix.

    Default:0

    C++ Type:unsigned int

    Controllable:No

    Description:The number of PISO iterations without recomputing the momentum matrix.

  • num_steps4294967295The number of timesteps in a transient run

    Default:4294967295

    C++ Type:unsigned int

    Controllable:No

    Description:The number of timesteps in a transient run

  • print_fieldsFalseUse this to print the coupling and solution fields and matrices throughout the iteration.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Use this to print the coupling and solution fields and matrices throughout the iteration.

  • reset_dtFalseUse when restarting a calculation to force a change in dt.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Use when restarting a calculation to force a change in dt.

  • schemeimplicit-eulerTime integration scheme used.

    Default:implicit-euler

    C++ Type:MooseEnum

    Options:implicit-euler, explicit-euler, crank-nicolson, bdf2, explicit-midpoint, dirk, explicit-tvd-rk-2, newmark-beta

    Controllable:No

    Description:Time integration scheme used.

  • solid_energy_systemThe solver system for the solid energy equation.

    C++ Type:SolverSystemName

    Controllable:No

    Description:The solver system for the solid energy equation.

  • verboseFalseSet to true to print additional information

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Set to true to print additional information

Optional Parameters

  • abort_on_solve_failFalseabort if solve not converged rather than cut timestep

    Default:False

    C++ Type:bool

    Controllable:No

    Description:abort if solve not converged rather than cut timestep

  • control_tagsAdds user-defined labels for accessing object parameters via control logic.

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:Adds user-defined labels for accessing object parameters via control logic.

  • dtmax1e+30The maximum timestep size in an adaptive run

    Default:1e+30

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The maximum timestep size in an adaptive run

  • dtmin1e-12The minimum timestep size in an adaptive run

    Default:1e-12

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The minimum timestep size in an adaptive run

  • enableTrueSet the enabled status of the MooseObject.

    Default:True

    C++ Type:bool

    Controllable:No

    Description:Set the enabled status of the MooseObject.

  • n_startup_steps0The number of timesteps during startup

    Default:0

    C++ Type:int

    Controllable:No

    Description:The number of timesteps during startup

  • outputsVector of output names where you would like to restrict the output of variables(s) associated with this object

    C++ Type:std::vector<OutputName>

    Controllable:No

    Description:Vector of output names where you would like to restrict the output of variables(s) associated with this object

  • start_time0The start time of the simulation

    Default:0

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The start time of the simulation

  • timestep_tolerance1e-12the tolerance setting for final timestep size and sync times

    Default:1e-12

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:the tolerance setting for final timestep size and sync times

  • use_multiapp_dtFalseIf true then the dt for the simulation will be chosen by the MultiApps. If false (the default) then the minimum over the master dt and the MultiApps is used

    Default:False

    C++ Type:bool

    Controllable:No

    Description:If true then the dt for the simulation will be chosen by the MultiApps. If false (the default) then the minimum over the master dt and the MultiApps is used

Advanced Parameters

  • accept_on_max_fixed_point_iterationFalseTrue to treat reaching the maximum number of fixed point iterations as converged.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:True to treat reaching the maximum number of fixed point iterations as converged.

  • auto_advanceFalseWhether to automatically advance sub-applications regardless of whether their solve converges, for transient executioners only.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Whether to automatically advance sub-applications regardless of whether their solve converges, for transient executioners only.

  • custom_abs_tol1e-50The absolute nonlinear residual to shoot for during fixed point iterations. This check is performed based on postprocessor defined by the custom_pp residual.

    Default:1e-50

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The absolute nonlinear residual to shoot for during fixed point iterations. This check is performed based on postprocessor defined by the custom_pp residual.

  • custom_ppPostprocessor for custom fixed point convergence check.

    C++ Type:PostprocessorName

    Unit:(no unit assumed)

    Controllable:No

    Description:Postprocessor for custom fixed point convergence check.

  • custom_rel_tol1e-08The relative nonlinear residual drop to shoot for during fixed point iterations. This check is performed based on the postprocessor defined by custom_pp residual.

    Default:1e-08

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The relative nonlinear residual drop to shoot for during fixed point iterations. This check is performed based on the postprocessor defined by custom_pp residual.

  • direct_pp_valueFalseTrue to use direct postprocessor value (scaled by value on first iteration). False (default) to use difference in postprocessor value between fixed point iterations.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:True to use direct postprocessor value (scaled by value on first iteration). False (default) to use difference in postprocessor value between fixed point iterations.

  • disable_fixed_point_residual_norm_checkFalseDisable the residual norm evaluation thus the three parameters fixed_point_rel_tol, fixed_point_abs_tol and fixed_point_force_norms.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Disable the residual norm evaluation thus the three parameters fixed_point_rel_tol, fixed_point_abs_tol and fixed_point_force_norms.

  • fixed_point_abs_tol1e-50The absolute nonlinear residual to shoot for during fixed point iterations. This check is performed based on the main app's nonlinear residual.

    Default:1e-50

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The absolute nonlinear residual to shoot for during fixed point iterations. This check is performed based on the main app's nonlinear residual.

  • fixed_point_algorithmpicardThe fixed point algorithm to converge the sequence of problems.

    Default:picard

    C++ Type:MooseEnum

    Options:picard, secant, steffensen

    Controllable:No

    Description:The fixed point algorithm to converge the sequence of problems.

  • fixed_point_force_normsFalseForce the evaluation of both the TIMESTEP_BEGIN and TIMESTEP_END norms regardless of the existence of active MultiApps with those execute_on flags, default: false.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Force the evaluation of both the TIMESTEP_BEGIN and TIMESTEP_END norms regardless of the existence of active MultiApps with those execute_on flags, default: false.

  • fixed_point_max_its1Specifies the maximum number of fixed point iterations.

    Default:1

    C++ Type:unsigned int

    Controllable:No

    Description:Specifies the maximum number of fixed point iterations.

  • fixed_point_min_its1Specifies the minimum number of fixed point iterations.

    Default:1

    C++ Type:unsigned int

    Controllable:No

    Description:Specifies the minimum number of fixed point iterations.

  • fixed_point_rel_tol1e-08The relative nonlinear residual drop to shoot for during fixed point iterations. This check is performed based on the main app's nonlinear residual.

    Default:1e-08

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The relative nonlinear residual drop to shoot for during fixed point iterations. This check is performed based on the main app's nonlinear residual.

  • multiapp_fixed_point_convergenceName of the Convergence object to use to assess convergence of the MultiApp fixed point solve. If not provided, a default Convergence will be constructed internally from the executioner parameters.

    C++ Type:ConvergenceName

    Controllable:No

    Description:Name of the Convergence object to use to assess convergence of the MultiApp fixed point solve. If not provided, a default Convergence will be constructed internally from the executioner parameters.

  • relaxation_factor1Fraction of newly computed value to keep.Set between 0 and 2.

    Default:1

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:Fraction of newly computed value to keep.Set between 0 and 2.

  • transformed_postprocessorsList of main app postprocessors to transform during fixed point iterations

    C++ Type:std::vector<PostprocessorName>

    Unit:(no unit assumed)

    Controllable:No

    Description:List of main app postprocessors to transform during fixed point iterations

  • transformed_variablesList of main app variables to transform during fixed point iterations

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:List of main app variables to transform during fixed point iterations

Fixed Point Iterations Parameters

  • active_scalar_absolute_toleranceThe absolute tolerance(s) on the normalized residual(s) of the active scalar equation(s).

    C++ Type:std::vector<double>

    Unit:(no unit assumed)

    Controllable:No

    Description:The absolute tolerance(s) on the normalized residual(s) of the active scalar equation(s).

  • active_scalar_equation_relaxationThe relaxation which should be used for the active scalar equations. (=1 for no relaxation, diagonal dominance will still be enforced)

    C++ Type:std::vector<double>

    Unit:(no unit assumed)

    Controllable:No

    Description:The relaxation which should be used for the active scalar equations. (=1 for no relaxation, diagonal dominance will still be enforced)

  • active_scalar_l_abs_tol1e-10The absolute tolerance on the normalized residual in the linear solver of the active scalar equation(s).

    Default:1e-10

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The absolute tolerance on the normalized residual in the linear solver of the active scalar equation(s).

  • active_scalar_l_max_its10000The maximum allowed iterations in the linear solver of the turbulence equation.

    Default:10000

    C++ Type:unsigned int

    Controllable:No

    Description:The maximum allowed iterations in the linear solver of the turbulence equation.

  • active_scalar_l_tol1e-05The relative tolerance on the normalized residual in the linear solver of the active scalar equation(s).

    Default:1e-05

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The relative tolerance on the normalized residual in the linear solver of the active scalar equation(s).

  • active_scalar_petsc_optionsSingleton PETSc options for the active scalar equation(s)

    C++ Type:MultiMooseEnum

    Options:-dm_moose_print_embedding, -dm_view, -KSP_CONVERGED_REASON, -KSP_GMRES_MODIFIEDGRAMSCHMIDT, -KSP_MONITOR, -KSP_MONITOR_SNES_LG, -SNES_KSP_EW, -KSP_SNES_EW, -SNES_CONVERGED_REASON, -SNES_KSP, -SNES_LINESEARCH_MONITOR, -SNES_MF, -SNES_MF_OPERATOR, -SNES_MONITOR, -SNES_TEST_DISPLAY, -SNES_VIEW, -SNES_MONITOR_CANCEL

    Controllable:No

    Description:Singleton PETSc options for the active scalar equation(s)

  • active_scalar_petsc_options_inameNames of PETSc name/value pairs for the active scalar equation(s)

    C++ Type:MultiMooseEnum

    Options:-mat_fd_coloring_err, -mat_fd_type, -mat_mffd_type, -pc_asm_overlap, -pc_factor_levels, -pc_factor_mat_ordering_type, -pc_hypre_boomeramg_grid_sweeps_all, -pc_hypre_boomeramg_max_iter, -pc_hypre_boomeramg_strong_threshold, -pc_hypre_type, -pc_type, -sub_pc_type, -KSP_ATOL, -KSP_GMRES_RESTART, -KSP_MAX_IT, -KSP_PC_SIDE, -KSP_RTOL, -KSP_TYPE, -SUB_KSP_TYPE, -SNES_ATOL, -SNES_LINESEARCH_TYPE, -SNES_LS, -SNES_MAX_IT, -SNES_RTOL, -SNES_DIVERGENCE_TOLERANCE, -SNES_TYPE

    Controllable:No

    Description:Names of PETSc name/value pairs for the active scalar equation(s)

  • active_scalar_petsc_options_valueValues of PETSc name/value pairs (must correspond with "petsc_options_iname" for the active scalar equation(s)

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:Values of PETSc name/value pairs (must correspond with "petsc_options_iname" for the active scalar equation(s)

  • active_scalar_systemsThe solver system for each active scalar advection equation.

    C++ Type:std::vector<SolverSystemName>

    Controllable:No

    Description:The solver system for each active scalar advection equation.

Active Scalars Equations Parameters

  • check_auxFalseWhether to check the auxiliary system for convergence to steady-state. If false, then the nonlinear system is used.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Whether to check the auxiliary system for convergence to steady-state. If false, then the nonlinear system is used.

  • steady_state_detectionFalseWhether or not to check for steady state conditions

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Whether or not to check for steady state conditions

  • steady_state_start_time0Minimum amount of time to run before checking for steady state conditions.

    Default:0

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:Minimum amount of time to run before checking for steady state conditions.

  • steady_state_tolerance1e-08Whenever the relative residual changes by less than this the solution will be considered to be at steady state.

    Default:1e-08

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:Whenever the relative residual changes by less than this the solution will be considered to be at steady state.

Steady State Detection Parameters

  • energy_absolute_tolerance1e-05The absolute tolerance on the normalized residual of the energy equation.

    Default:1e-05

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The absolute tolerance on the normalized residual of the energy equation.

  • energy_equation_relaxation1The relaxation which should be used for the energy equation. (=1 for no relaxation, diagonal dominance will still be enforced)

    Default:1

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The relaxation which should be used for the energy equation. (=1 for no relaxation, diagonal dominance will still be enforced)

  • energy_l_abs_tol1e-10The absolute tolerance on the normalized residual in the linear solver of the energy equation.

    Default:1e-10

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The absolute tolerance on the normalized residual in the linear solver of the energy equation.

  • energy_l_max_its10000The maximum allowed iterations in the linear solver of the energy equation.

    Default:10000

    C++ Type:unsigned int

    Controllable:No

    Description:The maximum allowed iterations in the linear solver of the energy equation.

  • energy_l_tol1e-05The relative tolerance on the normalized residual in the linear solver of the energy equation.

    Default:1e-05

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The relative tolerance on the normalized residual in the linear solver of the energy equation.

  • energy_petsc_optionsSingleton PETSc options for the energy equation

    C++ Type:MultiMooseEnum

    Options:-dm_moose_print_embedding, -dm_view, -KSP_CONVERGED_REASON, -KSP_GMRES_MODIFIEDGRAMSCHMIDT, -KSP_MONITOR, -KSP_MONITOR_SNES_LG, -SNES_KSP_EW, -KSP_SNES_EW, -SNES_CONVERGED_REASON, -SNES_KSP, -SNES_LINESEARCH_MONITOR, -SNES_MF, -SNES_MF_OPERATOR, -SNES_MONITOR, -SNES_TEST_DISPLAY, -SNES_VIEW, -SNES_MONITOR_CANCEL

    Controllable:No

    Description:Singleton PETSc options for the energy equation

  • energy_petsc_options_inameNames of PETSc name/value pairs for the energy equation

    C++ Type:MultiMooseEnum

    Options:-mat_fd_coloring_err, -mat_fd_type, -mat_mffd_type, -pc_asm_overlap, -pc_factor_levels, -pc_factor_mat_ordering_type, -pc_hypre_boomeramg_grid_sweeps_all, -pc_hypre_boomeramg_max_iter, -pc_hypre_boomeramg_strong_threshold, -pc_hypre_type, -pc_type, -sub_pc_type, -KSP_ATOL, -KSP_GMRES_RESTART, -KSP_MAX_IT, -KSP_PC_SIDE, -KSP_RTOL, -KSP_TYPE, -SUB_KSP_TYPE, -SNES_ATOL, -SNES_LINESEARCH_TYPE, -SNES_LS, -SNES_MAX_IT, -SNES_RTOL, -SNES_DIVERGENCE_TOLERANCE, -SNES_TYPE

    Controllable:No

    Description:Names of PETSc name/value pairs for the energy equation

  • energy_petsc_options_valueValues of PETSc name/value pairs (must correspond with "petsc_options_iname" for the energy equation

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:Values of PETSc name/value pairs (must correspond with "petsc_options_iname" for the energy equation

Energy Equation Parameters

  • max_xfem_update4294967295Maximum number of times to update XFEM crack topology in a step due to evolving cracks

    Default:4294967295

    C++ Type:unsigned int

    Controllable:No

    Description:Maximum number of times to update XFEM crack topology in a step due to evolving cracks

  • update_xfem_at_timestep_beginFalseShould XFEM update the mesh at the beginning of the timestep

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Should XFEM update the mesh at the beginning of the timestep

Xfem Fixed Point Iterations Parameters

  • momentum_absolute_tolerance1e-05The absolute tolerance on the normalized residual of the momentum equation.

    Default:1e-05

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The absolute tolerance on the normalized residual of the momentum equation.

  • momentum_equation_relaxation1The relaxation which should be used for the momentum equation. (=1 for no relaxation, diagonal dominance will still be enforced)

    Default:1

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The relaxation which should be used for the momentum equation. (=1 for no relaxation, diagonal dominance will still be enforced)

  • momentum_l_abs_tol1e-50The absolute tolerance on the normalized residual in the linear solver of the momentum equation.

    Default:1e-50

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The absolute tolerance on the normalized residual in the linear solver of the momentum equation.

  • momentum_l_max_its10000The maximum allowed iterations in the linear solver of the momentum equation.

    Default:10000

    C++ Type:unsigned int

    Controllable:No

    Description:The maximum allowed iterations in the linear solver of the momentum equation.

  • momentum_l_tol1e-05The relative tolerance on the normalized residual in the linear solver of the momentum equation.

    Default:1e-05

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The relative tolerance on the normalized residual in the linear solver of the momentum equation.

  • momentum_petsc_optionsSingleton PETSc options for the momentum equation

    C++ Type:MultiMooseEnum

    Options:-dm_moose_print_embedding, -dm_view, -KSP_CONVERGED_REASON, -KSP_GMRES_MODIFIEDGRAMSCHMIDT, -KSP_MONITOR, -KSP_MONITOR_SNES_LG, -SNES_KSP_EW, -KSP_SNES_EW, -SNES_CONVERGED_REASON, -SNES_KSP, -SNES_LINESEARCH_MONITOR, -SNES_MF, -SNES_MF_OPERATOR, -SNES_MONITOR, -SNES_TEST_DISPLAY, -SNES_VIEW, -SNES_MONITOR_CANCEL

    Controllable:No

    Description:Singleton PETSc options for the momentum equation

  • momentum_petsc_options_inameNames of PETSc name/value pairs for the momentum equation

    C++ Type:MultiMooseEnum

    Options:-mat_fd_coloring_err, -mat_fd_type, -mat_mffd_type, -pc_asm_overlap, -pc_factor_levels, -pc_factor_mat_ordering_type, -pc_hypre_boomeramg_grid_sweeps_all, -pc_hypre_boomeramg_max_iter, -pc_hypre_boomeramg_strong_threshold, -pc_hypre_type, -pc_type, -sub_pc_type, -KSP_ATOL, -KSP_GMRES_RESTART, -KSP_MAX_IT, -KSP_PC_SIDE, -KSP_RTOL, -KSP_TYPE, -SUB_KSP_TYPE, -SNES_ATOL, -SNES_LINESEARCH_TYPE, -SNES_LS, -SNES_MAX_IT, -SNES_RTOL, -SNES_DIVERGENCE_TOLERANCE, -SNES_TYPE

    Controllable:No

    Description:Names of PETSc name/value pairs for the momentum equation

  • momentum_petsc_options_valueValues of PETSc name/value pairs (must correspond with "petsc_options_iname" for the momentum equation

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:Values of PETSc name/value pairs (must correspond with "petsc_options_iname" for the momentum equation

  • momentum_systemsThe solver system(s) for the momentum equation(s).

    C++ Type:std::vector<SolverSystemName>

    Controllable:No

    Description:The solver system(s) for the momentum equation(s).

Momentum Equation Parameters

  • passive_scalar_absolute_toleranceThe absolute tolerance(s) on the normalized residual(s) of the passive scalar equation(s).

    C++ Type:std::vector<double>

    Unit:(no unit assumed)

    Controllable:No

    Description:The absolute tolerance(s) on the normalized residual(s) of the passive scalar equation(s).

  • passive_scalar_equation_relaxationThe relaxation which should be used for the passive scalar equations. (=1 for no relaxation, diagonal dominance will still be enforced)

    C++ Type:std::vector<double>

    Unit:(no unit assumed)

    Controllable:No

    Description:The relaxation which should be used for the passive scalar equations. (=1 for no relaxation, diagonal dominance will still be enforced)

  • passive_scalar_l_abs_tol1e-10The absolute tolerance on the normalized residual in the linear solver of the passive scalar equation(s).

    Default:1e-10

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The absolute tolerance on the normalized residual in the linear solver of the passive scalar equation(s).

  • passive_scalar_l_max_its10000The maximum allowed iterations in the linear solver of the turbulence equation.

    Default:10000

    C++ Type:unsigned int

    Controllable:No

    Description:The maximum allowed iterations in the linear solver of the turbulence equation.

  • passive_scalar_l_tol1e-05The relative tolerance on the normalized residual in the linear solver of the passive scalar equation(s).

    Default:1e-05

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The relative tolerance on the normalized residual in the linear solver of the passive scalar equation(s).

  • passive_scalar_petsc_optionsSingleton PETSc options for the passive scalar equation(s)

    C++ Type:MultiMooseEnum

    Options:-dm_moose_print_embedding, -dm_view, -KSP_CONVERGED_REASON, -KSP_GMRES_MODIFIEDGRAMSCHMIDT, -KSP_MONITOR, -KSP_MONITOR_SNES_LG, -SNES_KSP_EW, -KSP_SNES_EW, -SNES_CONVERGED_REASON, -SNES_KSP, -SNES_LINESEARCH_MONITOR, -SNES_MF, -SNES_MF_OPERATOR, -SNES_MONITOR, -SNES_TEST_DISPLAY, -SNES_VIEW, -SNES_MONITOR_CANCEL

    Controllable:No

    Description:Singleton PETSc options for the passive scalar equation(s)

  • passive_scalar_petsc_options_inameNames of PETSc name/value pairs for the passive scalar equation(s)

    C++ Type:MultiMooseEnum

    Options:-mat_fd_coloring_err, -mat_fd_type, -mat_mffd_type, -pc_asm_overlap, -pc_factor_levels, -pc_factor_mat_ordering_type, -pc_hypre_boomeramg_grid_sweeps_all, -pc_hypre_boomeramg_max_iter, -pc_hypre_boomeramg_strong_threshold, -pc_hypre_type, -pc_type, -sub_pc_type, -KSP_ATOL, -KSP_GMRES_RESTART, -KSP_MAX_IT, -KSP_PC_SIDE, -KSP_RTOL, -KSP_TYPE, -SUB_KSP_TYPE, -SNES_ATOL, -SNES_LINESEARCH_TYPE, -SNES_LS, -SNES_MAX_IT, -SNES_RTOL, -SNES_DIVERGENCE_TOLERANCE, -SNES_TYPE

    Controllable:No

    Description:Names of PETSc name/value pairs for the passive scalar equation(s)

  • passive_scalar_petsc_options_valueValues of PETSc name/value pairs (must correspond with "petsc_options_iname" for the passive scalar equation(s)

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:Values of PETSc name/value pairs (must correspond with "petsc_options_iname" for the passive scalar equation(s)

  • passive_scalar_systemsThe solver system for each scalar advection equation.

    C++ Type:std::vector<SolverSystemName>

    Controllable:No

    Description:The solver system for each scalar advection equation.

Passive_Scalar Equation Parameters

  • pin_pressureFalseIf the pressure field needs to be pinned at a point.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:If the pressure field needs to be pinned at a point.

  • pressure_pin_pointThe point where the pressure needs to be pinned.

    C++ Type:libMesh::Point

    Controllable:No

    Description:The point where the pressure needs to be pinned.

  • pressure_pin_value0The value which needs to be enforced for the pressure.

    Default:0

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The value which needs to be enforced for the pressure.

Pressure Pin Parameters

  • pressure_absolute_tolerance1e-05The absolute tolerance on the normalized residual of the pressure equation.

    Default:1e-05

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The absolute tolerance on the normalized residual of the pressure equation.

  • pressure_l_abs_tol1e-10The absolute tolerance on the normalized residual in the linear solver of the pressure equation.

    Default:1e-10

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The absolute tolerance on the normalized residual in the linear solver of the pressure equation.

  • pressure_l_max_its10000The maximum allowed iterations in the linear solver of the pressure equation.

    Default:10000

    C++ Type:unsigned int

    Controllable:No

    Description:The maximum allowed iterations in the linear solver of the pressure equation.

  • pressure_l_tol1e-05The relative tolerance on the normalized residual in the linear solver of the pressure equation.

    Default:1e-05

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The relative tolerance on the normalized residual in the linear solver of the pressure equation.

  • pressure_petsc_optionsSingleton PETSc options for the pressure equation

    C++ Type:MultiMooseEnum

    Options:-dm_moose_print_embedding, -dm_view, -KSP_CONVERGED_REASON, -KSP_GMRES_MODIFIEDGRAMSCHMIDT, -KSP_MONITOR, -KSP_MONITOR_SNES_LG, -SNES_KSP_EW, -KSP_SNES_EW, -SNES_CONVERGED_REASON, -SNES_KSP, -SNES_LINESEARCH_MONITOR, -SNES_MF, -SNES_MF_OPERATOR, -SNES_MONITOR, -SNES_TEST_DISPLAY, -SNES_VIEW, -SNES_MONITOR_CANCEL

    Controllable:No

    Description:Singleton PETSc options for the pressure equation

  • pressure_petsc_options_inameNames of PETSc name/value pairs for the pressure equation

    C++ Type:MultiMooseEnum

    Options:-mat_fd_coloring_err, -mat_fd_type, -mat_mffd_type, -pc_asm_overlap, -pc_factor_levels, -pc_factor_mat_ordering_type, -pc_hypre_boomeramg_grid_sweeps_all, -pc_hypre_boomeramg_max_iter, -pc_hypre_boomeramg_strong_threshold, -pc_hypre_type, -pc_type, -sub_pc_type, -KSP_ATOL, -KSP_GMRES_RESTART, -KSP_MAX_IT, -KSP_PC_SIDE, -KSP_RTOL, -KSP_TYPE, -SUB_KSP_TYPE, -SNES_ATOL, -SNES_LINESEARCH_TYPE, -SNES_LS, -SNES_MAX_IT, -SNES_RTOL, -SNES_DIVERGENCE_TOLERANCE, -SNES_TYPE

    Controllable:No

    Description:Names of PETSc name/value pairs for the pressure equation

  • pressure_petsc_options_valueValues of PETSc name/value pairs (must correspond with "petsc_options_iname" for the pressure equation

    C++ Type:std::vector<std::string>

    Controllable:No

    Description:Values of PETSc name/value pairs (must correspond with "petsc_options_iname" for the pressure equation

  • pressure_systemThe solver system for the pressure equation.

    C++ Type:SolverSystemName

    Controllable:No

    Description:The solver system for the pressure equation.

  • pressure_variable_relaxation1The relaxation which should be used for the pressure variable (=1 for no relaxation).

    Default:1

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The relaxation which should be used for the pressure variable (=1 for no relaxation).

Pressure Equation Parameters

    Restart Parameters

    • solid_energy_absolute_tolerance1e-05The absolute tolerance on the normalized residual of the solid energy equation.

      Default:1e-05

      C++ Type:double

      Unit:(no unit assumed)

      Controllable:No

      Description:The absolute tolerance on the normalized residual of the solid energy equation.

    • solid_energy_l_abs_tol1e-10The absolute tolerance on the normalized residual in the linear solver of the solid energy equation.

      Default:1e-10

      C++ Type:double

      Unit:(no unit assumed)

      Controllable:No

      Description:The absolute tolerance on the normalized residual in the linear solver of the solid energy equation.

    • solid_energy_l_max_its10000The maximum allowed iterations in the linear solver of the solid energy equation.

      Default:10000

      C++ Type:unsigned int

      Controllable:No

      Description:The maximum allowed iterations in the linear solver of the solid energy equation.

    • solid_energy_l_tol1e-05The relative tolerance on the normalized residual in the linear solver of the solid energy equation.

      Default:1e-05

      C++ Type:double

      Unit:(no unit assumed)

      Controllable:No

      Description:The relative tolerance on the normalized residual in the linear solver of the solid energy equation.

    • solid_energy_petsc_optionsSingleton PETSc options for the solid energy equation

      C++ Type:MultiMooseEnum

      Options:-dm_moose_print_embedding, -dm_view, -KSP_CONVERGED_REASON, -KSP_GMRES_MODIFIEDGRAMSCHMIDT, -KSP_MONITOR, -KSP_MONITOR_SNES_LG, -SNES_KSP_EW, -KSP_SNES_EW, -SNES_CONVERGED_REASON, -SNES_KSP, -SNES_LINESEARCH_MONITOR, -SNES_MF, -SNES_MF_OPERATOR, -SNES_MONITOR, -SNES_TEST_DISPLAY, -SNES_VIEW, -SNES_MONITOR_CANCEL

      Controllable:No

      Description:Singleton PETSc options for the solid energy equation

    • solid_energy_petsc_options_inameNames of PETSc name/value pairs for the solid energy equation

      C++ Type:MultiMooseEnum

      Options:-mat_fd_coloring_err, -mat_fd_type, -mat_mffd_type, -pc_asm_overlap, -pc_factor_levels, -pc_factor_mat_ordering_type, -pc_hypre_boomeramg_grid_sweeps_all, -pc_hypre_boomeramg_max_iter, -pc_hypre_boomeramg_strong_threshold, -pc_hypre_type, -pc_type, -sub_pc_type, -KSP_ATOL, -KSP_GMRES_RESTART, -KSP_MAX_IT, -KSP_PC_SIDE, -KSP_RTOL, -KSP_TYPE, -SUB_KSP_TYPE, -SNES_ATOL, -SNES_LINESEARCH_TYPE, -SNES_LS, -SNES_MAX_IT, -SNES_RTOL, -SNES_DIVERGENCE_TOLERANCE, -SNES_TYPE

      Controllable:No

      Description:Names of PETSc name/value pairs for the solid energy equation

    • solid_energy_petsc_options_valueValues of PETSc name/value pairs (must correspond with "petsc_options_iname" for the solid energy equation

      C++ Type:std::vector<std::string>

      Controllable:No

      Description:Values of PETSc name/value pairs (must correspond with "petsc_options_iname" for the solid energy equation

    Solid Energy Equation Parameters

    • time_period_endsThe end times of time periods

      C++ Type:std::vector<double>

      Unit:(no unit assumed)

      Controllable:No

      Description:The end times of time periods

    • time_period_startsThe start times of time periods

      C++ Type:std::vector<double>

      Unit:(no unit assumed)

      Controllable:No

      Description:The start times of time periods

    • time_periodsThe names of periods

      C++ Type:std::vector<std::string>

      Controllable:No

      Description:The names of periods

    Time Periods Parameters

    Input Files

    References

    1. Christopher Greenshields and Henry Weller. Notes on Computational Fluid Dynamics: General Principles. CFD Direct Ltd, Reading, UK, 2022.[BibTeX]
    2. Raad I Issa. Solution of the implicitly discretised fluid flow equations by operator-splitting. Journal of computational physics, 62(1):40–65, 1986.[BibTeX]