Navier Stokes Solid Heat Transfer / PNSFVSolidHeatTransfer

Define the Navier Stokes porous media solid energy equation

Equation

This Physics object creates the kernels and boundary conditions to solve the heat transfer equation for the solid phase in a porous media flow problem.

where:

  • is the solid specific enthalpy, computed from the specific heat and the solid temperature

  • is the solid density

  • is the porosity

  • is the solid temperature

  • is the fluid temperature

  • the solid thermal conductivity

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

  • is the ambient convection volumetric heat transfer coefficient

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

The kernels potentially created for this equation are:

commentnote

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

Automatically created variables

The PNSFVSolidHeatTransferPhysics by default will automatically create a nonlinear variable for the solid phase temperature. It should be named as follows:

static const std::string T_solidus = "T_solidus";
(moose/modules/navier_stokes/include/base/NS.h)

Coupling with other Physics

In the presence of fluid flow, a Navier Stokes Fluid Heat Transfer / WCNSFVFluidHeatTransferPhysics should be created using the [Physics/NavierStokes/FluidHeatTransfer/<name>] syntax. The following input performs the coupling between the fluid equations and the solid temperature equations. The coupling between the fluid and solid domain is performed through a volumetric ambient convection term.

[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."}>>> = 'weakly-compressible'
        porous_medium_treatment<<<{"description": "Whether to use porous medium kernels or not."}>>> = true
        define_variables<<<{"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."}>>> = true

        pressure_variable<<<{"description": "If supplied, the system checks for available pressure variable. Otherwise, it is created within the action."}>>> = 'pressure'

        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-6 0'
        initial_pressure<<<{"description": "The initial pressure, assumed constant everywhere"}>>> = '${p_outlet}'

        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"}>>> = 'top bottom'
        momentum_wall_types<<<{"description": "Types of wall boundaries for the momentum equation"}>>> = 'noslip symmetry'

        outlet_boundaries<<<{"description": "Names of outlet boundaries"}>>> = 'right'
        momentum_outlet_types<<<{"description": "Types of outlet boundaries for the momentum equation"}>>> = 'fixed-pressure'
        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"}>>>]
      [fluid]
        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'
        effective_conductivity<<<{"description": "Whether the conductivity should be multiplied by porosity, or whether the provided conductivity is an effective conductivity taking porosity effects into account"}>>> = true
        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'

        initial_temperature<<<{"description": "The initial temperature, assumed constant everywhere"}>>> = '${T_inlet}'

        # See 'flow' for inlet boundaries
        energy_inlet_types<<<{"description": "Types for the inlet boundaries for the energy equation."}>>> = 'fixed-temperature'
        energy_inlet_function = '${T_inlet}'

        # See 'flow' for wall boundaries
        energy_wall_types<<<{"description": "Types for the wall boundaries for the energy equation."}>>> = 'heatflux heatflux'
        energy_wall_function = '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'
      []
    []

    [SolidHeatTransfer<<<{"href": "../../syntax/Physics/NavierStokes/SolidHeatTransfer/index.html"}>>>]
      [solid]
        block<<<{"description": "Blocks (subdomains) that this Physics is active on."}>>> = 0
        initial_temperature<<<{"description": "Initial value of the temperature variable"}>>> = 100
        transient<<<{"description": "Whether the physics is to be solved as a transient"}>>> = true
        # To match the previous test results
        solid_temperature_two_term_bc_expansion<<<{"description": "If a two-term Taylor expansion is needed for the determination of the boundary valuesof the temperature/energy."}>>> = true

        thermal_conductivity_solid<<<{"description": "Thermal conductivity, which may have different names depending on the subdomain. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '${k_s}'
        cp_solid<<<{"description": "Specific heat functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${cp_s}
        rho_solid<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho_s}

        fixed_temperature_boundaries<<<{"description": "Boundaries on which to apply a fixed temperature"}>>> = 'top'
        boundary_temperatures<<<{"description": "Functors to compute the heat flux on each boundary in 'fixed_temperature_boundaries'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '${top_side_temperature}'

        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_convection_temperature<<<{"description": "The fluid 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_fluid'
        verbose<<<{"description": "Flag to facilitate debugging a Physics"}>>> = true
      []
    []
  []
[]
(moose/modules/navier_stokes/test/tests/finite_volume/pwcns/channel-flow/2d-transient-physics.i)
warningwarning

Conjugate heat transfer on a surface on the boundary of the fluid domain is not currently implemented with the Physics syntax. Please use a FVConvectionCorrelationInterface for that purpose.

Input Parameters

  • blockBlocks (subdomains) that this Physics is active on.

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

    Controllable:No

    Description:Blocks (subdomains) that this Physics is active on.

  • 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

  • heat_source_functorFunctor providing the heat source. 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:Functor providing the heat source. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.

  • initial_temperature300Initial value of the temperature variable

    Default:300

    C++ Type:FunctionName

    Unit:(no unit assumed)

    Controllable:No

    Description:Initial value of the temperature variable

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

  • solid_temperature_variableT_solidName of the solid phase temperature variable

    Default:T_solid

    C++ Type:VariableName

    Unit:(no unit assumed)

    Controllable:No

    Description:Name of the solid phase temperature variable

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

  • temperature_scaling1Scaling factor for the heat conduction equation

    Default:1

    C++ Type:double

    Unit:(no unit assumed)

    Controllable:No

    Description:Scaling factor for the heat conduction equation

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

  • 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_convection_temperatureT_fluid The fluid 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.

    Default:T_fluid

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

    Unit:(no unit assumed)

    Controllable:No

    Description:The fluid 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.

Ambient Convection Parameters

  • boundary_heat_fluxesFunctors to compute the heat flux on each boundary in 'heat_flux_boundaries'. 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:Functors to compute the heat flux on each boundary in 'heat_flux_boundaries'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.

  • boundary_temperaturesFunctors to compute the heat flux on each boundary in 'fixed_temperature_boundaries'. 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:Functors to compute the heat flux on each boundary in 'fixed_temperature_boundaries'. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.

  • fixed_temperature_boundariesBoundaries on which to apply a fixed temperature

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

    Controllable:No

    Description:Boundaries on which to apply a fixed temperature

  • heat_flux_boundariesBoundaries on which to apply a heat flux

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

    Controllable:No

    Description:Boundaries on which to apply a heat flux

  • insulated_boundariesBoundaries on which to apply a zero heat flux

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

    Controllable:No

    Description:Boundaries on which to apply a zero heat flux

Thermal Boundaries Parameters

  • cp_solidcp_solidSpecific heat functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.

    Default:cp_solid

    C++ Type:MooseFunctorName

    Unit:(no unit assumed)

    Controllable:No

    Description:Specific heat functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.

  • rho_solidrho_solidDensity functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.

    Default:rho_solid

    C++ Type:MooseFunctorName

    Unit:(no unit assumed)

    Controllable:No

    Description:Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.

  • thermal_conductivity_blocksBlocks which each thermal conductivity is defined

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

    Controllable:No

    Description:Blocks which each thermal conductivity is defined

  • thermal_conductivity_solidThermal conductivity, which may have different names depending on the subdomain. 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:Thermal conductivity, which may have different names depending on the subdomain. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.

  • 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

  • density_functorDensity functor. 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:Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.

  • specific_heat_functorSpecific heat functor. 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:Specific heat functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.

Thermal 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

  • 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_blocksThe blocks where the heat source is present.

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

    Controllable:No

    Description:The blocks where the heat source is present.

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

Solid Porous Medium 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

  • solid_temperature_face_interpolationaverageThe numerical scheme to interpolate the temperature/energy to the face for conduction (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 for conduction (separate from the advected quantity interpolation).

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

Numerical Scheme Parameters