Navier Stokes Fluid Heat Transfer / WCNSFVFluidHeatTransferPhysics

Define the Navier Stokes weakly-compressible energy equation

Equation

This Physics object creates the kernels and boundary conditions to solve the advection-diffusion equation for the fluid temperature. For free flow in a non-porous media:

For flow in a porous medium:

where:

  • is the fluid enthalpy, computed from the specific heat

  • is the fluid density

  • is the porosity

  • is the fluid temperature

  • is the advecting velocity (clean flow)

  • is the advecting superficial velocity (porous media flow)

  • the fluid effective thermal conductivity

  • is the source term, corresponding to energy deposited directly in the fluid

  • is the ambient convection volumetric heat transfer coefficient

  • is the ambient temperature

The enthalpy is used in lieu of to be able to model gases with temperature dependent specific heat.

The kernels created for flow in a non-porous medium are:

For flow in a porous medium:

commentnote

Additional details on porous media flow equations can be found on this page.

Automatically defined variables

The WCNSFVFluidHeatTransferPhysics automatically sets up the variables which are necessary for solving the energy transport equation:

For the default names of other variables used in this action, visit this site.

Coupling with other Physics

The heat advection equation can be solved concurrently with the flow equations by combining both the Navier Stokes Fluid Heat Transfer / WCNSFVFluidHeatTransferPhysics and the Navier Stokes Flow / WCNSFVFlowPhysics. The following input performs this coupling for incompressible flow in a 2D 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'
      []
    []
  []
[]
(moose/modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-rc-transient-physics.i)

Input Parameters

  • add_energy_equationFalseWhether to add the energy equation. This parameter is not necessary if using the Physics syntax

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Whether to add the energy equation. This parameter is not necessary if 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.

  • effective_conductivityTrueWhether the conductivity should be multiplied by porosity, or whether the provided conductivity is an effective conductivity taking porosity effects into account

    Default:True

    C++ Type:bool

    Controllable:No

    Description:Whether the conductivity should be multiplied by porosity, or whether the provided conductivity is an effective conductivity taking porosity effects into account

  • energy_wall_boundariesWall boundaries to apply energy boundary conditions on. If not specified, the flow equation Physics wall boundaries will be used

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

    Controllable:No

    Description:Wall boundaries to apply energy boundary conditions on. If not specified, the flow equation Physics wall boundaries will be used

  • fluid_temperature_variableT_fluidName of the fluid temperature variable

    Default:T_fluid

    C++ Type:NonlinearVariableName

    Unit:(no unit assumed)

    Controllable:No

    Description:Name of the fluid temperature variable

  • initial_enthalpyInitial value of the enthalpy variable, only to be used when solving for enthalpy

    C++ Type:FunctionName

    Unit:(no unit assumed)

    Controllable:No

    Description:Initial value of the enthalpy variable, only to be used when solving for enthalpy

  • initial_temperature300The initial temperature, assumed constant everywhere

    Default:300

    C++ Type:FunctionName

    Unit:(no unit assumed)

    Controllable:No

    Description:The initial temperature, assumed constant everywhere

  • solve_for_enthalpyFalseWhether to solve for the enthalpy or the temperature of the fluid

    Default:False

    C++ Type:bool

    Controllable:No

    Description:Whether to solve for the enthalpy or the temperature of the fluid

  • 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

    Options:true, false, same_as_problem

    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

Optional 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

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

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

    Unit:(no unit assumed)

    Controllable:No

    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.

  • ambient_convection_blocksThe blocks where the ambient convection is present.

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

    Controllable:No

    Description:The blocks where the ambient convection is present.

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

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

    Unit:(no unit assumed)

    Controllable:No

    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.

Volumetric Heat Convection Parameters

  • coupled_flow_physicsWCNSFVFlowPhysics generating the velocities

    C++ Type:PhysicsName

    Controllable:No

    Description:WCNSFVFlowPhysics generating the velocities

  • coupled_turbulence_physicsTurbulence Physics coupled with this Physics

    C++ Type:PhysicsName

    Controllable:No

    Description:Turbulence Physics coupled with this Physics

Coupled Physics 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

  • energy_advection_interpolationupwindThe numerical scheme to use for interpolating energy/temperature, as an advected quantity, to the face.

    Default:upwind

    C++ Type:MooseEnum

    Options:average, upwind, sou, min_mod, vanLeer, quick, venkatakrishnan, skewness-corrected

    Controllable:No

    Description:The numerical scheme to use for interpolating energy/temperature, as an advected quantity, to the face.

  • energy_face_interpolationaverageThe numerical scheme to interpolate the temperature/energy to the face (separate from the advected quantity interpolation).

    Default:average

    C++ Type:MooseEnum

    Options:average, skewness-corrected

    Controllable:No

    Description:The numerical scheme to interpolate the temperature/energy to the face (separate from the advected quantity interpolation).

  • energy_scaling1The scaling factor for the energy variable.

    Default:1

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:The scaling factor for the energy variable.

  • energy_two_term_bc_expansionTrueIf a two-term Taylor expansion is needed for the determination of the boundary valuesof the temperature/energy.

    Default:True

    C++ Type:bool

    Controllable:No

    Description:If a two-term Taylor expansion is needed for the determination of the boundary valuesof the temperature/energy.

  • preconditioningdeferWhich preconditioning to use/add for this Physics, or whether to defer to the Preconditioning block, or another Physics

    Default:defer

    C++ Type:MooseEnum

    Options:default, defer

    Controllable:No

    Description:Which preconditioning to use/add for this Physics, or whether to defer to the Preconditioning block, or another Physics

Numerical Scheme Parameters

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

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

    Unit:(no unit assumed)

    Controllable:No

    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.

  • energy_inlet_typesTypes for the inlet boundaries for the energy equation.

    C++ Type:MultiMooseEnum

    Options:fixed-temperature, flux-mass, flux-velocity, heatflux

    Controllable:No

    Description:Types for the inlet boundaries for the energy equation.

Inlet Boundary Conditions Parameters

  • energy_wall_functorsFunctions for Dirichlet/Neumann boundaries in the energy equation. For wall types requiring multiple functions, the syntax is ::... So, 'convection' types are ':'. 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 Dirichlet/Neumann boundaries in the energy equation. For wall types requiring multiple functions, the syntax is ::... So, 'convection' types are ':'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.

  • energy_wall_typesTypes for the wall boundaries for the energy equation.

    C++ Type:MultiMooseEnum

    Options:fixed-temperature, heatflux, wallfunction, convection

    Controllable:No

    Description:Types for the wall boundaries for the energy equation.

Wall Boundary Conditions Parameters

  • external_heat_sourceThe name of a functor which contains the external heat source for the energy equation. 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:The name of a functor which contains the external heat source for the energy equation. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.

  • external_heat_source_coeff1Multiplier for the coupled heat source term.

    Default:1

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:Multiplier for the coupled heat source term.

Heat Source 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

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

    Default:cp

    C++ Type:MooseFunctorName

    Unit:(no unit assumed)

    Controllable:No

    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.

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

    Default:k

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

    Unit:(no unit assumed)

    Controllable:No

    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.

  • thermal_conductivity_blocksThe blocks where the user wants define different thermal conductivities.

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

    Controllable:No

    Description:The blocks where the user wants define different thermal conductivities.

  • use_external_enthalpy_materialFalseTo indicate if the enthalpy material is set up outside of the action.

    Default:False

    C++ Type:bool

    Controllable:No

    Description:To indicate if the enthalpy material is set up outside of the action.

Material Properties Parameters