- add_flow_equationsTrueWhether to add the flow equations. This parameter is not necessary when using the Physics syntax
Default:True
C++ Type:bool
Controllable:No
Description:Whether to add the flow equations. This parameter is not necessary when using the Physics syntax
- blockBlocks (subdomains) that this Physics is active on.
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:Blocks (subdomains) that this Physics is active on.
- compressibilityincompressibleCompressibility constraint for the Navier-Stokes equations.
Default:incompressible
C++ Type:MooseEnum
Controllable:No
Description:Compressibility constraint for the Navier-Stokes equations.
- density_for_gravity_termsIf specified, replaces the 'density' for the Boussinesq and gravity momentum kernels. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:If specified, replaces the 'density' for the Boussinesq and gravity momentum kernels. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
- fluid_temperature_variableIf supplied, the system checks for available fluid temperature variable. Otherwise, it is created within the action.
C++ Type:NonlinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:If supplied, the system checks for available fluid temperature variable. Otherwise, it is created within the action.
- hydraulic_separator_sidesetsSidesets which serve as hydraulic separators.
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:Sidesets which serve as hydraulic separators.
- porosityporosityThe name of the auxiliary variable for the porosity field. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
Default:porosity
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:The name of the auxiliary variable for the porosity field. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
- porous_medium_treatmentFalseWhether to use porous medium kernels or not.
Default:False
C++ Type:bool
Controllable:No
Description:Whether to use porous medium kernels or not.
- solve_for_dynamic_pressureFalseWhether to solve for the dynamic pressure instead of the total pressure
Default:False
C++ Type:bool
Controllable:No
Description:Whether to solve for the dynamic pressure instead of the total pressure
- system_namesnl0 Name of the solver system(s) for the variables. If a single name is specified, that system is used for all solver variables.
Default:nl0
C++ Type:std::vector<SolverSystemName>
Controllable:No
Description:Name of the solver system(s) for the variables. If a single name is specified, that system is used for all solver variables.
- transientsame_as_problemWhether the physics is to be solved as a transient
Default:same_as_problem
C++ Type:MooseEnum
Controllable:No
Description:Whether the physics is to be solved as a transient
- verboseFalseFlag to facilitate debugging a Physics
Default:False
C++ Type:bool
Controllable:No
Description:Flag to facilitate debugging a Physics
Navier Stokes Flow / WCNSFVFlowPhysics
Define the Navier Stokes weakly-compressible mass and momentum equations
Equations
This Physics object creates the kernels and boundary conditions to solve the mass and momentum Navier Stokes equations for the flow. The time derivatives in the mass equations are omitted for incompressible flow. For regular flow in a non-porous medium:
For porous media flow:
where:
is the density
is the dynamic viscosity
is the porosity
is the velocity (non-porous flow)
is the superficial velocity (porous flow)
is the pressure
is the gravity term
is the friction / inter-phase friction term
Additional details on porous media flow equations can be found on this page.
The kernels created for free flow for the mass equation:
WCNSFVMassTimeDerivative for weakly-compressible flow in a transient case
INSFVMassAdvection for the mass advection term
for porous media flow:
PWCNSFVMassTimeDerivative for weakly-compressible flow
PINSFVMassAdvection for mass advection
The kernels created for the momentum equation for free flow:
WCNSFVMomentumTimeDerivative for the time derivative for weakly-compressible flow in a transient case
INSFVMomentumTimeDerivative for the time derivative incompressible flow in a transient case
INSFVMomentumAdvection for the momentum advection term
INSFVMomentumDiffusion for the momentum diffusion term
INSFVMomentumPressure for the pressure gradient term
PINSFVMomentumFriction for the friction term if specified
INSFVMomentumGravity for the gravity term if specified
for porous media flow:
PINSFVMomentumTimeDerivative for the time derivative
PINSFVMomentumAdvection for the momentum advection term
PINSFVMomentumDiffusion for the momentum diffusion term
PINSFVMomentumPressure for the pressure gradient term
PINSFVMomentumFriction for the friction term if specified
PINSFVMomentumGravity for the gravity term if specified
Automatically defined variables
The WCNSFVFlowPhysics
automatically sets up the variables which are necessary for the solution of a given problem. These variables can then be used to couple fluid flow simulations with other physics. The list of variable names commonly used in the action syntax is presented below:
Velocities for non-porous-medium simulations:
(modules/navier_stokes/include/base/NS.h)static const std::string velocity_x = "vel_x"; static const std::string velocity_y = "vel_y"; static const std::string velocity_z = "vel_z";
Velocities for porous medium simulations:
(modules/navier_stokes/include/base/NS.h)static const std::string superficial_velocity_x = "superficial_vel_x"; static const std::string superficial_velocity_y = "superficial_vel_y"; static const std::string superficial_velocity_z = "superficial_vel_z";
Pressure:
(modules/navier_stokes/include/base/NS.h)static const std::string total_pressure = "pressure_total";
For the default names of other variables used in this action, visit this site.
Coupling with other Physics
The energy advection equation can be solved concurrently with the flow equations using an additional Navier Stokes Fluid Heat Transfer / WCNSFVFluidHeatTransferPhysics. The following input performs this coupling for incompressible flow in a 2D flow channel. No system parameters are passed, so the equations are solved in a fully coupled manner in the same nonlinear system.
[Physics<<<{"href": "../../syntax/Physics/index.html"}>>>]
[NavierStokes<<<{"href": "../../syntax/Physics/NavierStokes/index.html"}>>>]
[Flow<<<{"href": "../../syntax/Physics/NavierStokes/Flow/index.html"}>>>]
[flow]
compressibility<<<{"description": "Compressibility constraint for the Navier-Stokes equations."}>>> = 'incompressible'
density<<<{"description": "The name of the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'rho'
dynamic_viscosity<<<{"description": "The name of the dynamic viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'mu'
initial_velocity<<<{"description": "The initial velocity, assumed constant everywhere"}>>> = '${u_inlet} 1e-12 0'
initial_pressure<<<{"description": "The initial pressure, assumed constant everywhere"}>>> = 0.0
inlet_boundaries<<<{"description": "Names of inlet boundaries"}>>> = 'left'
momentum_inlet_types<<<{"description": "Types of inlet boundaries for the momentum equation."}>>> = 'fixed-velocity'
momentum_inlet_function = '${u_inlet} 0'
wall_boundaries<<<{"description": "Names of wall boundaries"}>>> = 'bottom top'
momentum_wall_types<<<{"description": "Types of wall boundaries for the momentum equation"}>>> = 'symmetry noslip'
outlet_boundaries<<<{"description": "Names of outlet boundaries"}>>> = 'right'
momentum_outlet_types<<<{"description": "Types of outlet boundaries for the momentum equation"}>>> = 'fixed-pressure-zero-gradient'
pressure_function = '${p_outlet}'
mass_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating density, as an advected quantity, to the face."}>>> = 'average'
momentum_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating momentum/velocity, as an advected quantity, to the face."}>>> = 'average'
[]
[]
[FluidHeatTransfer<<<{"href": "../../syntax/Physics/NavierStokes/FluidHeatTransfer/index.html"}>>>]
[heat]
thermal_conductivity<<<{"description": "The name of the fluid thermal conductivity for each block. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'k'
specific_heat<<<{"description": "The name of the specific heat. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'cp'
fluid_temperature_variable<<<{"description": "Name of the fluid temperature variable"}>>> = 'T_fluid'
initial_temperature<<<{"description": "The initial temperature, assumed constant everywhere"}>>> = '${T_inlet}'
energy_inlet_types<<<{"description": "Types for the inlet boundaries for the energy equation."}>>> = 'heatflux'
energy_inlet_functors<<<{"description": "Functions for fixed-value boundaries in the energy equation. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '${fparse u_inlet * rho * cp * T_inlet}'
energy_wall_types<<<{"description": "Types for the wall boundaries for the energy equation."}>>> = 'heatflux heatflux'
energy_wall_functors<<<{"description": "Functions for Dirichlet/Neumann boundaries in the energy equation. For wall types requiring multiple functions, the syntax is <function_1>:<function_2>:... So, 'convection' types are '<Tinf_function>:<htc_function>'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '0 0'
ambient_convection_alpha<<<{"description": "The heat exchange coefficients for each block in 'ambient_convection_blocks'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'h_cv'
ambient_temperature<<<{"description": "The ambient temperature for each block in 'ambient_convection_blocks'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'T_solid'
energy_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating energy/temperature, as an advected quantity, to the face."}>>> = 'average'
[]
[]
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-rc-transient-physics.i)Other advected scalar equations can be solved concurrently with the flow equations using an additional Navier Stokes Scalar Transport / WCNSFVScalarTransportPhysics. The following input performs this coupling for incompressible flow in a 2D flow channel. No system parameters are passed, so the equations are solved in a fully coupled manner in the same nonlinear system.
[Physics<<<{"href": "../../syntax/Physics/index.html"}>>>]
[NavierStokes<<<{"href": "../../syntax/Physics/NavierStokes/index.html"}>>>]
[Flow<<<{"href": "../../syntax/Physics/NavierStokes/Flow/index.html"}>>>]
[flow]
compressibility<<<{"description": "Compressibility constraint for the Navier-Stokes equations."}>>> = 'incompressible'
density<<<{"description": "The name of the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
dynamic_viscosity<<<{"description": "The name of the dynamic viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
inlet_boundaries<<<{"description": "Names of inlet boundaries"}>>> = 'left'
momentum_inlet_types<<<{"description": "Types of inlet boundaries for the momentum equation."}>>> = 'fixed-velocity'
momentum_inlet_function = '1 0'
wall_boundaries<<<{"description": "Names of wall boundaries"}>>> = 'top bottom'
momentum_wall_types<<<{"description": "Types of wall boundaries for the momentum equation"}>>> = 'noslip noslip'
outlet_boundaries<<<{"description": "Names of outlet boundaries"}>>> = 'right'
momentum_outlet_types<<<{"description": "Types of outlet boundaries for the momentum equation"}>>> = 'fixed-pressure'
pressure_function = '0'
mass_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating density, as an advected quantity, to the face."}>>> = 'average'
momentum_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating momentum/velocity, as an advected quantity, to the face."}>>> = 'average'
[]
[]
[FluidHeatTransfer<<<{"href": "../../syntax/Physics/NavierStokes/FluidHeatTransfer/index.html"}>>>]
[heat]
thermal_conductivity<<<{"description": "The name of the fluid thermal conductivity for each block. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${k}
specific_heat<<<{"description": "The name of the specific heat. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${cp}
energy_inlet_types<<<{"description": "Types for the inlet boundaries for the energy equation."}>>> = 'fixed-temperature'
energy_inlet_function = '1'
energy_wall_types<<<{"description": "Types for the wall boundaries for the energy equation."}>>> = 'heatflux heatflux'
energy_wall_function = '0 0'
energy_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating energy/temperature, as an advected quantity, to the face."}>>> = 'average'
[]
[]
[ScalarTransport<<<{"href": "../../syntax/Physics/NavierStokes/ScalarTransport/index.html"}>>>]
[heat]
passive_scalar_names<<<{"description": "Vector containing the names of the advected scalar variables."}>>> = 'scalar'
passive_scalar_diffusivity<<<{"description": "Functor names for the diffusivities used for the passive scalar fields. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${diff}
passive_scalar_source<<<{"description": "Passive scalar sources. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 0.1
passive_scalar_coupled_source<<<{"description": "Coupled variable names for the sources used for the passive scalar fields. If multiple sources for each equation are specified, major (outer) ordering by equation. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = U
passive_scalar_coupled_source_coeff<<<{"description": "Coupled variable multipliers for the sources used for the passive scalar fields. If multiple sources for each equation are specified, major (outer) ordering by equation."}>>> = 0.1
passive_scalar_inlet_types<<<{"description": "Types for the inlet boundaries for the passive scalar equation."}>>> = 'fixed-value'
passive_scalar_inlet_function = '1'
passive_scalar_advection_interpolation<<<{"description": "The numerical scheme to use for interpolating passive scalar field, as an advected quantity, to the face."}>>> = 'average'
[]
[]
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-scalar-transport-physics.i)Modeling options and implementation details
This physics only supports Rhie-Chow interpolation for the determination of face velocities in the advection terms. The face interpolation of the advected quantities (e.g. upwind, average) can be controlled through the *_advection_interpolation
physics parameters.
Bernoulli pressure jump treatment
Please see the Bernoulli pressure variable documentation for more information.
Input Parameters
- active__all__ If specified only the blocks named will be visited and made active
Default:__all__
C++ Type:std::vector<std::string>
Controllable:No
Description:If specified only the blocks named will be visited and made active
- 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.
- define_variablesTrueWhether to define variables if the variables with the specified names do not exist. Note that if the variables are defined externally from the Physics, the initial conditions will not be created in the Physics either.
Default:True
C++ Type:bool
Controllable:No
Description:Whether to define variables if the variables with the specified names do not exist. Note that if the variables are defined externally from the Physics, the initial conditions will not be created in the Physics either.
- ghost_layers2Number of layers of elements to ghost near process domain boundaries
Default:2
C++ Type:unsigned short
Controllable:No
Description:Number of layers of elements to ghost near process domain boundaries
- inactiveIf specified blocks matching these identifiers will be skipped.
C++ Type:std::vector<std::string>
Controllable:No
Description:If specified blocks matching these identifiers will be skipped.
Advanced Parameters
- boussinesq_approximationFalseTrue to have Boussinesq approximation
Default:False
C++ Type:bool
Controllable:No
Description:True to have Boussinesq approximation
- gravity0 0 0The gravitational acceleration vector.
Default:0 0 0
C++ Type:libMesh::VectorValue<double>
Unit:(no unit assumed)
Controllable:No
Description:The gravitational acceleration vector.
- ref_temperature273.15Value for reference temperature in case of Boussinesq approximation
Default:273.15
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:Value for reference temperature in case of Boussinesq approximation
- thermal_expansionalphaThe name of the thermal expansion coefficient in the Boussinesq approximation. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
Default:alpha
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:The name of the thermal expansion coefficient in the Boussinesq approximation. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
Gravity Treatment Parameters
- characteristic_speedThe characteristic speed. For porous medium simulations, this characteristic speed should correspond to the superficial velocity, not the interstitial.
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:The characteristic speed. For porous medium simulations, this characteristic speed should correspond to the superficial velocity, not the interstitial.
- include_deviatoric_stressFalseWhether to include the full expansion (the transposed term as well) of the stress tensor
Default:False
C++ Type:bool
Controllable:No
Description:Whether to include the full expansion (the transposed term as well) of the stress tensor
- mass_advection_interpolationupwindThe numerical scheme to use for interpolating density, as an advected quantity, to the face.
Default:upwind
C++ Type:MooseEnum
Controllable:No
Description:The numerical scheme to use for interpolating density, as an advected quantity, to the face.
- mass_scaling1The scaling factor for the mass variables (for incompressible simulation this is pressure scaling).
Default:1
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:The scaling factor for the mass variables (for incompressible simulation this is pressure scaling).
- momentum_advection_interpolationupwindThe numerical scheme to use for interpolating momentum/velocity, as an advected quantity, to the face.
Default:upwind
C++ Type:MooseEnum
Controllable:No
Description:The numerical scheme to use for interpolating momentum/velocity, as an advected quantity, to the face.
- momentum_face_interpolationaverageThe numerical scheme to interpolate the velocity/momentum to the face (separate from the advected quantity interpolation).
Default:average
C++ Type:MooseEnum
Controllable:No
Description:The numerical scheme to interpolate the velocity/momentum to the face (separate from the advected quantity interpolation).
- momentum_scaling1The scaling factor for the momentum variables.
Default:1
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:The scaling factor for the momentum variables.
- momentum_two_term_bc_expansionTrueIf a two-term Taylor expansion is needed for the determination of the boundary valuesof the velocity/momentum.
Default:True
C++ Type:bool
Controllable:No
Description:If a two-term Taylor expansion is needed for the determination of the boundary valuesof the velocity/momentum.
- mu_interp_methodharmonicSwitch that can select face interpolation method for the viscosity.
Default:harmonic
C++ Type:MooseEnum
Controllable:No
Description:Switch that can select face interpolation method for the viscosity.
- preconditioningdeferWhich preconditioning to use/add for this Physics, or whether to defer to the Preconditioning block, or another Physics
Default:defer
C++ Type:MooseEnum
Controllable:No
Description:Which preconditioning to use/add for this Physics, or whether to defer to the Preconditioning block, or another Physics
- pressure_face_interpolationaverageThe numerical scheme to interpolate the pressure to the face (separate from the advected quantity interpolation).
Default:average
C++ Type:MooseEnum
Controllable:No
Description:The numerical scheme to interpolate the pressure to the face (separate from the advected quantity interpolation).
- pressure_two_term_bc_expansionTrueIf a two-term Taylor expansion is needed for the determination of the boundary valuesof the pressure.
Default:True
C++ Type:bool
Controllable:No
Description:If a two-term Taylor expansion is needed for the determination of the boundary valuesof the pressure.
- time_derivative_contributes_to_RC_coefficientsTrueWhether the time derivative term should contribute to the Rhie Chow coefficients. This adds stabilization, but makes the solution dependent on the time step size
Default:True
C++ Type:bool
Controllable:No
Description:Whether the time derivative term should contribute to the Rhie Chow coefficients. This adds stabilization, but makes the solution dependent on the time step size
- velocity_interpolationrcThe interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow.
Default:rc
C++ Type:MooseEnum
Controllable:No
Description:The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow.
Numerical Scheme Parameters
- consistent_scalingScaling parameter for the friction correction in the momentum equation (if requested).
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:Scaling parameter for the friction correction in the momentum equation (if requested).
- porosity_interface_pressure_treatmentautomaticHow to treat pressure at a porosity interface
Default:automatic
C++ Type:MooseEnum
Controllable:No
Description:How to treat pressure at a porosity interface
- porosity_smoothing_layersThe number of interpolation-reconstruction operations to perform on the porosity.
C++ Type:unsigned short
Controllable:No
Description:The number of interpolation-reconstruction operations to perform on the porosity.
- pressure_allow_expansion_on_bernoulli_facesFalseSwitch to enable the two-term extrapolation on porosity jump faces. WARNING: Depending on the mesh, enabling this parameter may lead to termination in parallel runs due to insufficient ghosting between processors. An example can be the presence of multiple porosity jumps separated by only one cell while using the Bernoulli pressure treatment. In such cases adjust the `ghost_layers` parameter.
Default:False
C++ Type:bool
Controllable:No
Description:Switch to enable the two-term extrapolation on porosity jump faces. WARNING: Depending on the mesh, enabling this parameter may lead to termination in parallel runs due to insufficient ghosting between processors. An example can be the presence of multiple porosity jumps separated by only one cell while using the Bernoulli pressure treatment. In such cases adjust the `ghost_layers` parameter.
- pressure_drop_form_factorsUser-supplied form loss coefficients to be applied over the sidesets listed above
C++ Type:std::vector<double>
Unit:(no unit assumed)
Controllable:No
Description:User-supplied form loss coefficients to be applied over the sidesets listed above
- pressure_drop_sidesetsSidesets over which form loss coefficients are to be applied
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:Sidesets over which form loss coefficients are to be applied
- use_friction_correctionFalseIf friction correction should be applied in the momentum equation.
Default:False
C++ Type:bool
Controllable:No
Description:If friction correction should be applied in the momentum equation.
Flow Medium Discontinuity Treatment Parameters
- coupled_turbulence_physicsTurbulence Physics coupled with the flow
C++ Type:PhysicsName
Controllable:No
Description:Turbulence Physics coupled with the flow
Coupled Physics Parameters
- densityrhoThe name of the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
Default:rho
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:The name of the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
- dynamic_viscositymuThe name of the dynamic viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
Default:mu
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:The name of the dynamic viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
Material Properties Parameters
- dont_create_aux_kernelsFalseWhether to skip the 'add_aux_kernel' task
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip the 'add_aux_kernel' task
- dont_create_aux_variablesFalseWhether to skip the 'add_aux_variable' task
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip the 'add_aux_variable' task
- dont_create_bcsFalseWhether to skip the 'add_bc' task for each boundary condition type
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip the 'add_bc' task for each boundary condition type
- dont_create_correctorsFalseWhether to skip the 'add_correctors' task
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip the 'add_correctors' task
- dont_create_functionsFalseWhether to skip the 'add_function' task
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip the 'add_function' task
- dont_create_icsFalseWhether to skip the 'add_ic' task
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip the 'add_ic' task
- dont_create_kernelsFalseWhether to skip the 'add_kernel' task for each kernel type
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip the 'add_kernel' task for each kernel type
- dont_create_materialsFalseWhether to skip the 'add_material' task for each material type
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip the 'add_material' task for each material type
- dont_create_postprocessorsFalseWhether to skip the 'add_postprocessors' task
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip the 'add_postprocessors' task
- dont_create_solver_variablesFalseWhether to skip the 'add_variable' task
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip the 'add_variable' task
- dont_create_user_objectsFalseWhether to skip the 'add_user_object' task. This does not apply to UserObject derived classes being created on a different task (for example: postprocessors, VPPs, correctors)
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip the 'add_user_object' task. This does not apply to UserObject derived classes being created on a different task (for example: postprocessors, VPPs, correctors)
- dont_create_vectorpostprocessorsFalseWhether to skip the 'add_vectorpostprocessors' task
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip the 'add_vectorpostprocessors' task
Reduce Physics Object Creation Parameters
- flux_inlet_directionsThe directions which can be used to define the orientation of the flux with respect to the mesh. This can be used to define a flux which is incoming with an angle or to adjust the flux direction with respect to the normal. If the inlet surface is defined on an internal face, this is necessary to ensure the arbitrary orientation of the normal does not result in non-physical results.
C++ Type:std::vector<libMesh::Point>
Controllable:No
Description:The directions which can be used to define the orientation of the flux with respect to the mesh. This can be used to define a flux which is incoming with an angle or to adjust the flux direction with respect to the normal. If the inlet surface is defined on an internal face, this is necessary to ensure the arbitrary orientation of the normal does not result in non-physical results.
- flux_inlet_ppsThe name of the postprocessors which compute the mass flow/ velocity magnitude. Mainly used for coupling between different applications.
C++ Type:std::vector<PostprocessorName>
Unit:(no unit assumed)
Controllable:No
Description:The name of the postprocessors which compute the mass flow/ velocity magnitude. Mainly used for coupling between different applications.
Boundary Condition Parameters
- friction_blocksThe blocks where the friction factors are applied to emulate flow resistances.
C++ Type:std::vector<std::vector<SubdomainName>>
Controllable:No
Description:The blocks where the friction factors are applied to emulate flow resistances.
- friction_coeffsThe friction coefficients for every item in 'friction_types'. Note that if 'porous_medium_treatment' is enabled, the coefficients already contain a velocity multiplier but they are not multiplied with density yet!
C++ Type:std::vector<std::vector<std::string>>
Controllable:No
Description:The friction coefficients for every item in 'friction_types'. Note that if 'porous_medium_treatment' is enabled, the coefficients already contain a velocity multiplier but they are not multiplied with density yet!
- friction_typesThe types of friction forces for every block in 'friction_blocks'.
C++ Type:std::vector<std::vector<std::string>>
Controllable:No
Description:The types of friction forces for every block in 'friction_blocks'.
- standard_friction_formulationTrueFlag to enable the standard friction formulation or its alternative, which is a simplified version (see user documentation for PINSFVMomentumFriction).
Default:True
C++ Type:bool
Controllable:No
Description:Flag to enable the standard friction formulation or its alternative, which is a simplified version (see user documentation for PINSFVMomentumFriction).
Friction Control Parameters
- initial_from_file_timestepLATESTGives the time step number (or "LATEST") for which to read the Exodus solution
Default:LATEST
C++ Type:std::string
Controllable:No
Description:Gives the time step number (or "LATEST") for which to read the Exodus solution
- initialize_variables_from_mesh_fileFalseDetermines if the variables that are added by the action are initializedfrom the mesh file (only for Exodus format)
Default:False
C++ Type:bool
Controllable:No
Description:Determines if the variables that are added by the action are initializedfrom the mesh file (only for Exodus format)
Restart From Exodus Parameters
- initial_pressure1e5The initial pressure, assumed constant everywhere
Default:1e5
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:No
Description:The initial pressure, assumed constant everywhere
- initial_velocity1e-15 1e-15 1e-15 The initial velocity, assumed constant everywhere
Default:1e-15 1e-15 1e-15
C++ Type:std::vector<FunctionName>
Unit:(no unit assumed)
Controllable:No
Description:The initial velocity, assumed constant everywhere
- pressure_variableIf supplied, the system checks for available pressure variable. Otherwise, it is created within the action.
C++ Type:NonlinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:If supplied, the system checks for available pressure variable. Otherwise, it is created within the action.
- velocity_variableIf supplied, the system checks for available velocity variables. Otherwise, they are created within the action.
C++ Type:std::vector<std::string>
Controllable:No
Description:If supplied, the system checks for available velocity variables. Otherwise, they are created within the action.
Variables Parameters
- inlet_boundariesNames of inlet boundaries
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:Names of inlet boundaries
- momentum_inlet_functorsFunctions for inlet boundary velocities or pressures (for fixed-pressure option). Provide a double vector where the leading dimension corresponds to the number of fixed-velocity and fixed-pressure entries in momentum_inlet_types and the second index runs either over dimensions for fixed-velocity boundaries or is a single function name for pressure inlets. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
C++ Type:std::vector<std::vector<MooseFunctorName>>
Unit:(no unit assumed)
Controllable:No
Description:Functions for inlet boundary velocities or pressures (for fixed-pressure option). Provide a double vector where the leading dimension corresponds to the number of fixed-velocity and fixed-pressure entries in momentum_inlet_types and the second index runs either over dimensions for fixed-velocity boundaries or is a single function name for pressure inlets. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
- momentum_inlet_typesTypes of inlet boundaries for the momentum equation.
C++ Type:MultiMooseEnum
Controllable:No
Description:Types of inlet boundaries for the momentum equation.
Inlet Boundary Conditions Parameters
- momentum_outlet_typesTypes of outlet boundaries for the momentum equation
C++ Type:MultiMooseEnum
Controllable:No
Description:Types of outlet boundaries for the momentum equation
- outlet_boundariesNames of outlet boundaries
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:Names of outlet boundaries
- pressure_functorsFunctions for boundary pressures at outlets. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
C++ Type:std::vector<MooseFunctorName>
Unit:(no unit assumed)
Controllable:No
Description:Functions for boundary pressures at outlets. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
Outlet Boundary Conditions Parameters
- momentum_wall_functorsFunctors for each component of the velocity value on walls. This is only necessary for the fixed-velocity momentum wall types. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
C++ Type:std::vector<std::vector<MooseFunctorName>>
Unit:(no unit assumed)
Controllable:No
Description:Functors for each component of the velocity value on walls. This is only necessary for the fixed-velocity momentum wall types. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
- momentum_wall_typesTypes of wall boundaries for the momentum equation
C++ Type:MultiMooseEnum
Controllable:No
Description:Types of wall boundaries for the momentum equation
- wall_boundariesNames of wall boundaries
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:Names of wall boundaries
Wall Boundary Conditions Parameters
- pin_pressureFalseSwitch to enable pressure shifting for incompressible simulations.
Default:False
C++ Type:bool
Controllable:No
Description:Switch to enable pressure shifting for incompressible simulations.
- pinned_pressure_point0 0 0The XYZ coordinates where pressure needs to be pinned for incompressible simulations.
Default:0 0 0
C++ Type:libMesh::Point
Controllable:No
Description:The XYZ coordinates where pressure needs to be pinned for incompressible simulations.
- pinned_pressure_typeaverage-uoTypes for shifting (pinning) the pressure in case of incompressible simulations.
Default:average-uo
C++ Type:MooseEnum
Controllable:No
Description:Types for shifting (pinning) the pressure in case of incompressible simulations.
- pinned_pressure_value1e5The value used for pinning the pressure (point value/domain average).
Default:1e5
C++ Type:PostprocessorName
Unit:(no unit assumed)
Controllable:No
Description:The value used for pinning the pressure (point value/domain average).