SIMPLE

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

Overview

This executioner is based on the algorithm proposed by Patankar and Spalding (1983). The algorithm is based on the splitting of operators and successive correction for the momentum and pressure fields. The formulation implemented in MOOSE has been presented in Jasak (1996) and Juretic (2005). See also the examples and derivations in Moukalled et al. (2016). The concept relies on deriving a pressure equation using the discretized form of the momentum equations together with the continuity constraint. Let's take the steady-state incompressible Navier-Stokes equations in the following form:

(1)(2)

Where denotes the velocity, the pressure, the density, and the effective dynamic viscosity which potentially includes the contributions of eddy viscosity derived from turbulence models. Term expresses a volumetric source term which can be potentially velocity-dependent. As a first step, we assume that we have a guess for the pressure field, therefore the gradient is known. Furthermore, we assume that the advecting velocity field is known from the previous iteration. By explicitly showing the iteration index, Eq. (1) and Eq. (2) become:

(3)(4)

At this point, we should note that the finite volume discretization in MOOSE uses a collocated formulation which has an advantage of being flexible for unstructured meshes. However, in certain scenarios it can exhibit numerical pressure checker-boarding due to the discretization of the pressure gradient and continuity terms. A common approach for tackling this issue is the utilization of the Rhie-Chow interpolation method (See Rhie and Chow (1983) and Moukalled et al. (2016) for a detailed explanation). This means that the face velocities (or face fluxes) are determined using pressure corrections. As we will see later, due to this behavior, the iteration between pressure and velocity will in fact be an iteration between pressure and face velocity. Nevertheless, to keep this in mind we add a subscript to the advecting velocity in our formulation:

(5)(6)

Next, we split the operator acting on in the momentum equation into two components: a component that incorporates effects that result in contributions to the diagonal of a soon-to-be-generated system matrix and another component that contains everything else. With this in mind, we can rewrite the equation the following, semi-discretized way:

(7)

where is the diagonal contribution, and includes the off-diagonal contributions multiplied by the solution together with any additional volumetric source and sink terms (i.e. the discretized forms of ). One can solve this equation to obtain a new guess for the velocity field. This guess, however, will not respect the continuity equation, therefore we need to correct it. For this, a pressure equation is derived from the following formulation:

(8)

By applying the inverse of the diagonal operator (a very cheap process computationally), we arrive to the following expression:

(9)

By applying the continuity equation onto (which is a constraint) and assuming that the Rhie-Chow interpolation is used for the velocity, we arrive to a Poisson equation for pressure:

(10)

This equation is solved for a pressure which can be used to correct the face velocities in a sense that they respect the continuity equation. This correction already involves a Rhie-Chow interpolation, considering that the and fields are interpolated to the faces in a discretized form:

(11)

This correction applies the continuity constraint in an iterative manner, while ensuring the lack of numerical pressure checker-boarding phenomena.

The next guess for the velocity, however, does not necessarily respect the momentum equation. Therefore, the momentum prediction and pressure correction steps need to be repeated until both the momentum and continuity equations are satisfied.

commentnote

The iterative process above is not stable if the full update is applied every time. This means that the variables need to be relaxed. Specifically, it is a common practice to relax the pressure when plugging it back to the gradient term in the momentum predictor:

(12)

where is the relaxed field and is the corresponding relaxation parameter.

commentnote

To help the solution process of the linear solver, we add options to ensure diagonal dominance through the relaxation of equations. This is done using the method mentioned in Juretic (2005), meaning that a numerical correction is added to the diagonal of the system matrix and the right hand side. This is especially useful for advection-dominated systems.

commentnote

Currently, this solver only respects the following execute_on flags: INITAL, TIMESTEP_BEGIN, and FINAL, other flags are ignored. MultiApps and the corresponding MultiappTransfers are executed at FINAL only.

Example Input Syntax

The setup of a problem with the segregated solver in MOOSE is slightly different compared to conventional monolithic solvers. In this section, we highlight the main differences. For setting up a 2D simulation with the SIMPLE algorithm, we need three linear systems in MOOSE: one for each momentum component and another for the pressure. The different systems can be created within the Problem block:

[Problem<<<{"href": "../../syntax/Problem/index.html"}>>>]
  linear_sys_names = 'u_system v_system pressure_system'
  previous_nl_solution_required = true
[]
(moose/modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-velocity-pressure.i)

It is visible that we requested that MOOSE keeps previous solution iterates as well. This is necessary to facilitate the relaxation processes mentioned in the overview. Next, we create linear FV variables and assign them to the given systems.

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [vel_x]
    type = MooseLinearVariableFVReal<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../variables/MooseLinearVariableFV.html"}>>>
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.5
    solver_sys<<<{"description": "If this variable is a solver variable, this is the solver system to which it should be added."}>>> = u_system
  []
  [vel_y]
    type = MooseLinearVariableFVReal<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../variables/MooseLinearVariableFV.html"}>>>
    solver_sys<<<{"description": "If this variable is a solver variable, this is the solver system to which it should be added."}>>> = v_system
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.0
  []
  [pressure]
    type = MooseLinearVariableFVReal<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../variables/MooseLinearVariableFV.html"}>>>
    solver_sys<<<{"description": "If this variable is a solver variable, this is the solver system to which it should be added."}>>> = pressure_system
    initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.2
  []
[]
(moose/modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-velocity-pressure.i)

The kernels are then created within the LinearFVKernels block. The fundamental terms that contribute to the face fluxes in the momentum equation (stress and advection terms) are lumped into one kernel. Furthermore, instead of adding contribution from the continuity equation, we build an anisotropic diffusion (Poisson) equation for pressure:

[LinearFVKernels<<<{"href": "../../syntax/LinearFVKernels/index.html"}>>>]
  [u_advection_stress]
    type = LinearWCNSFVMomentumFlux<<<{"description": "Represents the matrix and right hand side contributions of the stress and advection terms of the momentum equation.", "href": "../linearfvkernels/LinearWCNSFVMomentumFlux.html"}>>>
    variable<<<{"description": "The name of the variable whose linear system this object contributes to"}>>> = vel_x
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    mu<<<{"description": "The diffusion coefficient. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    u<<<{"description": "The velocity in the x direction."}>>> = vel_x
    v<<<{"description": "The velocity in the y direction."}>>> = vel_y
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
    rhie_chow_user_object<<<{"description": "The rhie-chow user-object which is used to determine the face velocity."}>>> = 'rc'
    use_nonorthogonal_correction<<<{"description": "If the nonorthogonal correction should be used when computing the normal gradient."}>>> = false
  []
  [v_advection_stress]
    type = LinearWCNSFVMomentumFlux<<<{"description": "Represents the matrix and right hand side contributions of the stress and advection terms of the momentum equation.", "href": "../linearfvkernels/LinearWCNSFVMomentumFlux.html"}>>>
    variable<<<{"description": "The name of the variable whose linear system this object contributes to"}>>> = vel_y
    advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
    mu<<<{"description": "The diffusion coefficient. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
    u<<<{"description": "The velocity in the x direction."}>>> = vel_x
    v<<<{"description": "The velocity in the y direction."}>>> = vel_y
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
    rhie_chow_user_object<<<{"description": "The rhie-chow user-object which is used to determine the face velocity."}>>> = 'rc'
    use_nonorthogonal_correction<<<{"description": "If the nonorthogonal correction should be used when computing the normal gradient."}>>> = false
  []
  [u_pressure]
    type = LinearFVMomentumPressure<<<{"description": "Represents the pressure gradient term in the Navier Stokes momentum equations, added to the right hand side.", "href": "../linearfvkernels/LinearFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable whose linear system this object contributes to"}>>> = vel_x
    pressure<<<{"description": "The pressure variable whose gradient should be used."}>>> = pressure
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
  []
  [v_pressure]
    type = LinearFVMomentumPressure<<<{"description": "Represents the pressure gradient term in the Navier Stokes momentum equations, added to the right hand side.", "href": "../linearfvkernels/LinearFVMomentumPressure.html"}>>>
    variable<<<{"description": "The name of the variable whose linear system this object contributes to"}>>> = vel_y
    pressure<<<{"description": "The pressure variable whose gradient should be used."}>>> = pressure
    momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
  []
  [p_diffusion]
    type = LinearFVAnisotropicDiffusion<<<{"description": "Represents the matrix and right hand side contributions of a diffusion term in a partial differential equation.", "href": "../linearfvkernels/LinearFVAnisotropicDiffusion.html"}>>>
    variable<<<{"description": "The name of the variable whose linear system this object contributes to"}>>> = pressure
    diffusion_tensor<<<{"description": "Functor describing a diagonal diffusion tensor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = Ainv
    use_nonorthogonal_correction<<<{"description": "If the nonorthogonal correction should be used when computing the normal gradient."}>>> = false
  []
  [HbyA_divergence]
    type = LinearFVDivergence<<<{"description": "Represents a divergence term. Note, this term does not contribute to the system matrix, only takes the divergence of a face flux field and adds it to the right hand side of the linear system.", "href": "../linearfvkernels/LinearFVDivergence.html"}>>>
    variable<<<{"description": "The name of the variable whose linear system this object contributes to"}>>> = pressure
    face_flux<<<{"description": "Functor for the face flux. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = HbyA
    force_boundary_execution<<<{"description": "Whether to force execution of this object on all external boundaries."}>>> = true
  []
[]
(moose/modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-velocity-pressure.i)

By default, the coupling fields corresponding to and are handled by functor called HbyA and Ainv, respectively. These fields are generated by RhieChowMassFlux under the hood. This means that we need to add the user object responsible for generating these fields:

[UserObjects<<<{"href": "../../syntax/UserObjects/index.html"}>>>]
  [rc]
    type = RhieChowMassFlux<<<{"description": "Computes H/A and 1/A together with face mass fluxes for segregated momentum-pressure equations using linear systems.", "href": "../userobjects/RhieChowMassFlux.html"}>>>
    u<<<{"description": "The x-component of velocity"}>>> = vel_x
    v<<<{"description": "The y-component of velocity"}>>> = vel_y
    pressure<<<{"description": "The pressure variable."}>>> = pressure
    rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
    p_diffusion_kernel<<<{"description": "The diffusion kernel acting on the pressure."}>>> = p_diffusion
  []
[]
(moose/modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-velocity-pressure.i)

As a last step, we add the SIMPLE executioner:

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  type = SIMPLE
  momentum_l_abs_tol = 1e-10
  pressure_l_abs_tol = 1e-10
  momentum_l_tol = 0
  pressure_l_tol = 0
  rhie_chow_user_object = 'rc'
  momentum_systems = 'u_system v_system'
  pressure_system = 'pressure_system'
  momentum_equation_relaxation = 0.8
  pressure_variable_relaxation = 0.3
  num_iterations = 100
  pressure_absolute_tolerance = 1e-10
  momentum_absolute_tolerance = 1e-10
  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'
  print_fields = false
[]
(moose/modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-velocity-pressure.i)

Passive scalar advection

The SIMPLE executioner can be used to solve coupled problems involving both flow and passive scalar advection. Advected passive scalars do not affect the flow distribution, and therefore can be solved after the velocity and pressure fields have been computed using the SIMPLE algorithm. Several systems may be used, for each passive scalar.

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.

  • energy_systemThe solver system for the energy equation.

    C++ Type:SolverSystemName

    Controllable:No

    Description:The solver system for the energy equation.

  • 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.

  • 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.

  • solid_energy_systemThe solver system for the solid energy equation.

    C++ Type:SolverSystemName

    Controllable:No

    Description:The solver system for the solid energy equation.

  • time0System time

    Default:0

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:System time

  • verboseFalseSet to true to print additional information

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Set to true to print additional information

Optional 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

  • 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.

  • enableTrueSet the enabled status of the MooseObject.

    Default:True

    C++ Type:bool

    Controllable:No

    Description:Set the enabled status of the MooseObject.

  • 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

Advanced 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

    References

    1. Hrvoje Jasak. Error analysis and estimation for the finite volume method with applications to fluid flows. PhD thesis, Imperial College London (University of London), 1996.[BibTeX]
    2. Franjo Juretic. Error analysis in finite volume CFD. PhD thesis, Imperial College London (University of London), 2005.[BibTeX]
    3. Fadl Moukalled, L Mangani, Marwan Darwish, and others. The finite volume method in computational fluid dynamics. Volume 6. Springer, 2016.[BibTeX]
    4. Suhas V Patankar and D Brian Spalding. A calculation procedure for heat, mass and momentum transfer in three-dimensional parabolic flows. In Numerical prediction of flow, heat transfer, turbulence and combustion, pages 54–73. Elsevier, 1983.[BibTeX]
    5. Chae M Rhie and Wei-Liang Chow. Numerical study of the turbulent flow past an airfoil with trailing edge separation. AIAA journal, 21(11):1525–1532, 1983.[BibTeX]