- boundaryThe list of boundary IDs from the mesh where this object applies
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:The list of boundary IDs from the mesh where this object applies
- variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable that this residual object operates on
ADConvectiveHeatFluxBC
Convective heat transfer boundary condition with temperature and heat transfer coefficient given by material properties.
This boundary condition computes convective heat flux , where is convective heat transfer coefficient, is the temperature, and is far field temperature. Both and are coupled as material properties.
[./right]
type = ADConvectiveHeatFluxBC
variable = temp
boundary = 'right'
T_infinity = 100.0
heat_transfer_coefficient = 1
[../]
(modules/heat_transfer/test/tests/ad_convective_heat_flux/equilibrium.i)Input Parameters
- T_infinityMaterial property for far-field temperature
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:Material property for far-field temperature
- T_infinity_functorFunctor for far-field temperature. 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 for far-field temperature. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
- displacementsThe displacements
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The displacements
- heat_transfer_coefficientMaterial property for heat transfer coefficient
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:Material property for heat transfer coefficient
- heat_transfer_coefficient_functorFunctor for heat transfer coefficient. 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 for heat transfer coefficient. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
- matrix_onlyFalseWhether this object is only doing assembly to matrices (no vectors)
Default:False
C++ Type:bool
Controllable:No
Description:Whether this object is only doing assembly to matrices (no vectors)
Optional Parameters
- absolute_value_vector_tagsThe tags for the vectors this residual object should fill with the absolute value of the residual contribution
C++ Type:std::vector<TagName>
Controllable:No
Description:The tags for the vectors this residual object should fill with the absolute value of the residual contribution
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
Description:The extra tags for the matrices this Kernel should fill
- extra_vector_tagsThe extra tags for the vectors this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
Description:The extra tags for the vectors this Kernel should fill
- matrix_tagssystemThe tag for the matrices this Kernel should fill
Default:system
C++ Type:MultiMooseEnum
Options:nontime, system
Controllable:No
Description:The tag for the matrices this Kernel should fill
- vector_tagsnontimeThe tag for the vectors this Kernel should fill
Default:nontime
C++ Type:MultiMooseEnum
Options:nontime, time
Controllable:No
Description:The tag for the vectors this Kernel should fill
Contribution To Tagged Field Data 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.
- diag_save_inThe name of auxiliary variables to save this BC's diagonal jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The name of auxiliary variables to save this BC's diagonal jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:Yes
Description:Set the enabled status of the MooseObject.
- implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Controllable:No
Description:Determines whether this object is calculated using an implicit or explicit form
- save_inThe name of auxiliary variables to save this BC's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The name of auxiliary variables to save this BC's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Controllable:No
Description:The seed for the master random number generator
- skip_execution_outside_variable_domainFalseWhether to skip execution of this boundary condition when the variable it applies to is not defined on the boundary. This can facilitate setups with moving variable domains and fixed boundaries. Note that the FEProblem boundary-restricted integrity checks will also need to be turned off if using this option
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip execution of this boundary condition when the variable it applies to is not defined on the boundary. This can facilitate setups with moving variable domains and fixed boundaries. Note that the FEProblem boundary-restricted integrity checks will also need to be turned off if using this option
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
- use_interpolated_stateFalseFor the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Default:False
C++ Type:bool
Controllable:No
Description:For the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Material Property Retrieval Parameters
Input Files
- (tutorials/shield_multiphysics/inputs/step06_transient_heat_conduction/step6_transient.i)
- (modules/combined/test/tests/electromagnetic_joule_heating/fusing_current_through_copper_wire.i)
- (tutorials/shield_multiphysics/inputs/step08_adaptivity/step8_adapt.i)
- (modules/heat_transfer/test/tests/ad_convective_heat_flux/fe_fv_coupled_error.i)
- (tutorials/shield_multiphysics/inputs/step13_restart/step13b_initialization_from_exodus.i)
- (modules/heat_transfer/test/tests/radiative_bcs/ad_function_radiative_bc.i)
- (modules/heat_transfer/test/tests/ad_convective_heat_flux/coupled.i)
- (tutorials/shield_multiphysics/inputs/step11_multiapps/step11_2d_heat_conduction.i)
- (tutorials/shield_multiphysics/inputs/step13_restart/step13c_restart_from_checkpoint.i)
- (tutorials/shield_multiphysics/inputs/step08_adaptivity/step8_uniform.i)
- (tutorials/shield_multiphysics/inputs/step07_mechanics/step7.i)
- (tutorials/shield_multiphysics/inputs/step05_auxiliary_variables/step5.i)
- (tutorials/shield_multiphysics/inputs/step04_heat_conduction/step4.i)
- (modules/heat_transfer/test/tests/ad_convective_heat_flux/fe_fv_coupled.i)
- (modules/heat_transfer/test/tests/ad_convective_heat_flux/equilibrium.i)
- (modules/heat_transfer/test/tests/ad_convective_heat_flux/flux.i)
- (tutorials/shield_multiphysics/inputs/step06_transient_heat_conduction/step6_pseudo_transient.i)
- (tutorials/shield_multiphysics/inputs/step13_restart/step13a_base_calc.i)
(modules/heat_transfer/test/tests/ad_convective_heat_flux/equilibrium.i)
[Mesh]
type = GeneratedMesh
dim = 3
nx = 10
[]
[Variables]
[./temp]
initial_condition = 200.0
[../]
[]
[Kernels]
[./heat_dt]
type = ADTimeDerivative
variable = temp
[../]
[./heat_conduction]
type = ADDiffusion
variable = temp
[../]
[]
[BCs]
[./right]
type = ADConvectiveHeatFluxBC
variable = temp
boundary = 'right'
T_infinity = 100.0
heat_transfer_coefficient = 1
[../]
[]
[Postprocessors]
[./left_temp]
type = SideAverageValue
variable = temp
boundary = left
execute_on = 'TIMESTEP_END initial'
[../]
[./right_temp]
type = SideAverageValue
variable = temp
boundary = right
[../]
[./right_flux]
type = SideDiffusiveFluxAverage
variable = temp
boundary = right
diffusivity = 1
[../]
[]
[Executioner]
type = Transient
num_steps = 10
dt = 1e1
nl_abs_tol = 1e-12
[]
[Outputs]
[./out]
type = CSV
time_step_interval = 10
[../]
[]
(tutorials/shield_multiphysics/inputs/step06_transient_heat_conduction/step6_transient.i)
[Mesh]
[fmg]
type = FileMeshGenerator
file = '../step03_boundary_conditions/mesh_in.e'
[]
[]
[Variables]
[T]
# Adds a Linear Lagrange variable by default
block = 'concrete_hd concrete Al'
initial_condition = 300
[]
[]
[Kernels]
[diffusion_concrete]
type = ADHeatConduction
variable = T
[]
[time_derivative]
type = ADHeatConductionTimeDerivative
variable = T
[]
[]
[Materials]
[concrete_hd]
type = ADHeatConductionMaterial
block = concrete_hd
temp = 'T'
# we specify a function of time, temperature is passed as the time argument
# in the material
thermal_conductivity_temperature_function = '5.0 + 0.001 * t'
specific_heat = 1050
[]
[concrete]
type = ADHeatConductionMaterial
block = concrete
temp = 'T'
thermal_conductivity_temperature_function = '2.25 + 0.001 * t'
specific_heat = 1050
[]
[Al]
type = ADHeatConductionMaterial
block = Al
temp = T
thermal_conductivity_temperature_function = '175'
specific_heat = 875
[]
[density_concrete_hd]
type = ADGenericConstantMaterial
block = 'concrete_hd'
prop_names = 'density'
prop_values = '3524' # kg / m3
[]
[density_concrete]
type = ADGenericConstantMaterial
block = 'concrete'
prop_names = 'density'
prop_values = '2403' # kg / m3
[]
[density_Al]
type = ADGenericConstantMaterial
block = 'Al'
prop_names = 'density'
prop_values = '2270' # kg / m3
[]
[]
[BCs]
[from_reactor]
type = NeumannBC
variable = T
boundary = inner_cavity_solid
# 5 MW reactor, only 50 kW removed from radiation, 144 m2 cavity area
value = '${fparse 5e4 / 144}'
[]
[air_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'air_boundary'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 10
[]
[ground]
type = DirichletBC
variable = T
value = 300
boundary = 'ground'
[]
[water_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'water_boundary_inwards'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 600
[]
[]
[Problem]
# No kernels on the water domain
kernel_coverage_check = false
# No materials on the water domain
material_coverage_check = false
[]
[Executioner]
type = Transient
num_steps = 10
dt = ${units 12 h -> s}
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(modules/combined/test/tests/electromagnetic_joule_heating/fusing_current_through_copper_wire.i)
# This test is a simpified coupled case between the electromagnetic and
# heat transfer modules. While the file microwave_heating.i is a test
# utilizing the method of manufactured solutions, where both real and
# complex components of the electromagnetic properties are provided
# (such that no term is zeroed out), this test involves only the
# real components of the electromagnetic properties. In particular,
# this test supplies the fusing current to a copper wire and simulations
# the spatial and temporal heating profile until the wire reaches its
# melting point. The PDE's of this test file are as follows:
#
# curl(curl(A)) + j*mu*omega*(sigma*A) = J
# mag(E) = mag(-j*omega*A) + mag(J/sigma)
# rho*C*dT/dt - div(k*grad(T)) = Q
# Q = 0.5*sigma*mag(E)^2
#
# Where:
# - A is the magnetic vector potential
# - j is the sqrt(-1)
# - mu is the permeability of free space
# - omega is the angular frequency of the system
# - sigma is the electric conductivity of the wire
# - J is the supplied DC current
# - E is the electric field
# - rho is the density of copper
# - C is the heat capacity of copper
# - T is the temperature
# - k is the thermal conductivity of the wire
# - Q is the Joule heating
#
# The BCs are as follows:
#
# curl(n) x curl(A) = 0, where n is the normal vector
# q * n = h (T - T_infty), where q is the heat flux,
# h is the convective heat transfer coefficient,
# and T_infty is the far-field temperature.
[Mesh]
# Mesh of the copper wire
[fmg]
type = FileMeshGenerator
file = copper_wire.msh
[]
[]
[Variables]
# The real and complex components of the magnetic vector
# potential in the frequency domain
[A_real]
family = NEDELEC_ONE
order = FIRST
[]
[A_imag]
family = NEDELEC_ONE
order = FIRST
[]
# The temperature of the air in the copper wire
[T]
initial_condition = 293.0 #in K
[]
[]
[Kernels]
### Physics to determine the magnetic vector potential propagation ###
# The propagation of the real component
[curl_curl_real]
type = CurlCurlField
variable = A_real
[]
# Current induced by the electrical conductivity
# of the copper wire
[conduction_real]
type = ADConductionCurrent
variable = A_real
field_imag = A_imag
field_real = A_real
conductivity_real = electrical_conductivity
conductivity_imag = 0.0
ang_freq_real = omega_real
ang_freq_imag = 0.0
permeability_real = mu_real
permeability_imag = 0.0
component = real
[]
# Current supplied to the wire
[source_real]
type = VectorBodyForce
variable = A_real
function = mu_curr_real
[]
# The propagation of the complex component
[curl_curl_imag]
type = CurlCurlField
variable = A_imag
[]
# Current induced by the electrical conductivity
# of the copper wire
[conduction_imag]
type = ADConductionCurrent
variable = A_imag
field_imag = A_imag
field_real = A_real
conductivity_real = electrical_conductivity
conductivity_imag = 0.0
ang_freq_real = omega_real
ang_freq_imag = 0.0
permeability_real = mu_real
permeability_imag = 0.0
component = imaginary
[]
### Physics to determine the heat transfer ###
# Heat transfer in the copper wire
[HeatTdot_in_copper]
type = ADHeatConductionTimeDerivative
variable = T
specific_heat = specific_heat_copper
density_name = density_copper
[]
[HeatDiff_in_copper]
type = ADHeatConduction
variable = T
thermal_conductivity = thermal_conductivity_copper
[]
# Heating due the total current
[HeatSrc]
type = ADJouleHeatingSource
variable = T
heating_term = 'electric_field_heating'
[]
[]
[AuxVariables]
# Decomposing the magnetic vector potential
# for the electric field calculations
[A_x_real]
family = MONOMIAL
order = FIRST
[]
[A_y_real]
family = MONOMIAL
order = FIRST
[]
[A_x_imag]
family = MONOMIAL
order = FIRST
[]
[A_y_imag]
family = MONOMIAL
order = FIRST
[]
# The electrical conductivity for the electric
# field calculations
[elec_cond]
family = MONOMIAL
order = FIRST
[]
# The electric field profile determined from
# the magnetic vector potential
[E_real]
family = NEDELEC_ONE
order = FIRST
[]
[E_imag]
family = NEDELEC_ONE
order = FIRST
[]
[]
[AuxKernels]
# Decomposing the magnetic vector potential
# for the electric field calculations
[A_x_real]
type = VectorVariableComponentAux
variable = A_x_real
vector_variable = A_real
component = X
[]
[A_y_real]
type = VectorVariableComponentAux
variable = A_y_real
vector_variable = A_real
component = Y
[]
[A_x_imag]
type = VectorVariableComponentAux
variable = A_x_imag
vector_variable = A_imag
component = X
[]
[A_y_imag]
type = VectorVariableComponentAux
variable = A_y_imag
vector_variable = A_imag
component = Y
[]
# The electrical conductivity for the electric
# field calculations
[cond]
type = ADMaterialRealAux
property = electrical_conductivity
variable = elec_cond
execute_on = 'INITIAL LINEAR TIMESTEP_END'
[]
# The magnitude of electric field profile determined
# from the magnetic vector potential using:
# abs(E) = abs(-j*omega*A) + abs(supplied current / elec_cond)
# NOTE: The reason for calculating the magnitude of the electric
# field is the heating term is defined as:
# Q = 1/2 abs(E)^2 for frequency domain field formulations
[E_real]
type = ParsedVectorAux
coupled_variables = 'A_x_imag A_y_imag elec_cond'
expression_x = 'abs(2*3.14*60*A_x_imag) + abs(60e6/elec_cond)'
expression_y = 'abs(2*3.14*60*A_y_imag)'
variable = E_real
[]
[E_imag]
type = ParsedVectorAux
coupled_variables = 'A_x_real A_y_real'
expression_x = 'abs(-2*3.14*60*A_x_real)'
expression_y = 'abs(-2*3.14*60*A_y_real)'
variable = E_imag
[]
[]
[Functions]
# The supplied current density to the wire
# where only the real x-component is considered
[curr_real_x]
type = ParsedFunction
expression = '60e6' # Units in A/m^2, equivalent to 1178 A in a 5mm diameter wire
[]
# Permeability of free space
[mu_real_func]
type = ParsedFunction
expression = '4*pi*1e-7' # Units in N/A^2
[]
# The angular drive frequency of the system
[omega_real_func]
type = ParsedFunction
expression = '2*pi*60' # Units in rad/s
[]
# The angular frequency time permeability of free space
[omegaMu]
type = ParsedFunction
symbol_names = 'omega mu'
symbol_values = 'omega_real_func mu_real_func'
expression = 'omega*mu'
[]
# The supplied current density time permeability of free space
[mu_curr_real]
type = ParsedVectorFunction
symbol_names = 'current_mag mu'
symbol_values = 'curr_real_x mu_real_func'
expression_x = 'mu * current_mag'
[]
[]
[BCs]
### Temperature boundary conditions ###
# Convective heat flux BC with copper wire
# exposed to air
[surface]
type = ADConvectiveHeatFluxBC
variable = T
boundary = walls
T_infinity = 293
heat_transfer_coefficient = 10
[]
### Magnetic vector potential boundary conditions ###
# No defined boundary conditions represents
# zero curl conditions at the boundaries, such that:
# A x n = 0
[]
[Materials]
[k]
type = ADGenericConstantMaterial
prop_names = 'thermal_conductivity_copper'
prop_values = '397.48' #in W/(m K)
[]
[cp]
type = ADGenericConstantMaterial
prop_names = 'specific_heat_copper'
prop_values = '385.0' #in J/(kg K)
[]
[rho]
type = ADGenericConstantMaterial
prop_names = 'density_copper'
prop_values = '8920.0' #in kg/(m^3)
[]
# Electrical conductivity (copper is default material)
[sigma]
type = ADElectricalConductivity
temperature = T
block = copper
[]
# Material that supplies the correct Joule heating formulation
[ElectromagneticMaterial]
type = ElectromagneticHeatingMaterial
electric_field = E_real
complex_electric_field = E_imag
electric_field_heating_name = electric_field_heating
electrical_conductivity = electrical_conductivity
formulation = FREQUENCY
solver = ELECTROMAGNETIC
block = copper
[]
# Coefficient for wave propagation
[mu_real]
type = ADGenericFunctionMaterial
prop_names = mu_real
prop_values = mu_real_func
[]
[omega_real]
type = ADGenericFunctionMaterial
prop_names = omega_real
prop_values = omega_real_func
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
scheme = bdf2
solve_type = NEWTON
line_search = NONE
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
dt = 1.0
# NOTE: Change 'end_time' to 10s to accurately simulate the fusing current
# end_time = 10
end_time = 5
automatic_scaling = true
[]
[Outputs]
exodus = true
perf_graph = true
[]
(tutorials/shield_multiphysics/inputs/step08_adaptivity/step8_adapt.i)
[Mesh]
[fmg]
type = FileMeshGenerator
file = '../step03_boundary_conditions/mesh_in.e'
[]
[]
[Adaptivity]
marker = jump_threshold
max_h_level = 2
[Indicators]
[temperature_jump]
type = GradientJumpIndicator
variable = T
scale_by_flux_faces = true
[]
[]
[Markers]
[jump_threshold]
type = ValueThresholdMarker
coarsen = 0.3
variable = temperature_jump
refine = 2
block = 'concrete_hd concrete Al'
[]
[]
[]
[Variables]
[T]
# Adds a Linear Lagrange variable by default
block = 'concrete_hd concrete Al'
initial_condition = 300
[]
[]
[Kernels]
[diffusion_concrete]
type = ADHeatConduction
variable = T
[]
[time_derivative]
type = ADHeatConductionTimeDerivative
variable = T
[]
[]
[Materials]
[concrete_hd]
type = ADHeatConductionMaterial
block = concrete_hd
temp = 'T'
# we specify a function of time, temperature is passed as the time argument
# in the material
thermal_conductivity_temperature_function = '5.0 + 0.001 * t'
specific_heat = 1050
[]
[concrete]
type = ADHeatConductionMaterial
block = concrete
temp = 'T'
thermal_conductivity_temperature_function = '2.25 + 0.001 * t'
specific_heat = 1050
[]
[Al]
type = ADHeatConductionMaterial
block = Al
temp = T
thermal_conductivity_temperature_function = '175'
specific_heat = 875
[]
[density_concrete_hd]
type = ADGenericConstantMaterial
block = 'concrete_hd'
prop_names = 'density'
prop_values = '3524' # kg / m3
[]
[density_concrete]
type = ADGenericConstantMaterial
block = 'concrete'
prop_names = 'density'
prop_values = '2403' # kg / m3
[]
[density_Al]
type = ADGenericConstantMaterial
block = 'Al'
prop_names = 'density'
prop_values = '2270' # kg / m3
[]
[]
[BCs]
[from_reactor]
type = NeumannBC
variable = T
boundary = inner_cavity_solid
# 5 MW reactor, only 50 kW removed from radiation, 144 m2 cavity area
value = '${fparse 5e4 / 144}'
[]
[air_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'air_boundary'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 10
[]
[ground]
type = DirichletBC
variable = T
value = 300
boundary = 'ground'
[]
[water_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'water_boundary_inwards'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 600
[]
[]
[Problem]
# No kernels on the water domain
kernel_coverage_check = false
# No materials on the water domain
material_coverage_check = false
[]
[Executioner]
type = Transient
num_steps = 5
dt = ${units 12 h -> s}
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(modules/heat_transfer/test/tests/ad_convective_heat_flux/fe_fv_coupled_error.i)
[Mesh]
[cmg]
type = CartesianMeshGenerator
dim = 2
dx = '0.5 0.5'
dy = '0.5 0.5'
ix = '5 5'
iy = '5 5'
subdomain_id = '0 1
2 3'
[]
[add_sideset]
type = SideSetsBetweenSubdomainsGenerator
input = cmg
new_boundary = 'middle middle'
primary_block = '0 2'
paired_block = '1 3'
[]
[]
[Variables]
[u_fe]
block = '0 2 3'
[]
[]
[AuxVariables]
[u_fv]
type = MooseVariableFVReal
block = '0 3'
initial_condition = 1
[]
[]
[Kernels]
[u_fe_diff]
type = ADDiffusion
variable = u_fe
[]
[]
[BCs]
[u_fe_left]
type = ADDirichletBC
boundary = left
variable = u_fe
value = 0
[]
[u_fe_middle]
type = ADConvectiveHeatFluxBC
boundary = middle
variable = u_fe
T_infinity_functor = u_fv
heat_transfer_coefficient_functor = 1.0
[]
[]
[Executioner]
type = Steady
[]
[Problem]
kernel_coverage_check = false
[]
(tutorials/shield_multiphysics/inputs/step13_restart/step13b_initialization_from_exodus.i)
[Mesh]
[fmg]
type = FileMeshGenerator
file = 'step13a_base_calc_out.e'
use_for_exodus_restart = true
[]
[]
[Variables]
[T]
# Adds a Linear Lagrange variable by default
block = 'concrete_hd concrete Al'
initial_from_file_var = 'T'
[]
[]
[Kernels]
[diffusion_concrete]
type = ADHeatConduction
variable = T
[]
[time_derivative]
type = ADHeatConductionTimeDerivative
variable = T
[]
[]
[Materials]
[concrete_hd]
type = ADHeatConductionMaterial
block = concrete_hd
temp = 'T'
# we specify a function of time, temperature is passed as the time argument
# in the material
thermal_conductivity_temperature_function = '5.0 + 0.001 * t'
specific_heat = 1050
[]
[concrete]
type = ADHeatConductionMaterial
block = concrete
temp = 'T'
thermal_conductivity_temperature_function = '2.25 + 0.001 * t'
specific_heat = 1050
[]
[Al]
type = ADHeatConductionMaterial
block = Al
temp = T
thermal_conductivity_temperature_function = '175'
specific_heat = 875
[]
[density_concrete_hd]
type = ADGenericConstantMaterial
block = 'concrete_hd'
prop_names = 'density'
prop_values = '3524' # kg / m3
[]
[density_concrete]
type = ADGenericConstantMaterial
block = 'concrete'
prop_names = 'density'
prop_values = '2403' # kg / m3
[]
[density_Al]
type = ADGenericConstantMaterial
block = 'Al'
prop_names = 'density'
prop_values = '2270' # kg / m3
[]
[]
[BCs]
[from_reactor]
type = NeumannBC
variable = T
boundary = inner_cavity_solid
# 5 MW reactor, only 50 kW removed from radiation, 144 m2 cavity area
value = '${fparse 5e4 / 144}'
[]
[air_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'air_boundary'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 10
[]
[ground]
type = DirichletBC
variable = T
value = 300
boundary = 'ground'
[]
[water_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'water_boundary_inwards'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 600
[]
[]
[Problem]
# No kernels on the water domain
kernel_coverage_check = false
# No materials on the water domain
material_coverage_check = false
[]
[Executioner]
type = Transient
start_time = '${units 2 day -> s}'
num_steps = 6
dt = '${units 12 h -> s}'
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(modules/heat_transfer/test/tests/radiative_bcs/ad_function_radiative_bc.i)
#
# If we assume that epsilon*sigma*(T_inf^4-T_s^4) is approximately equal to
# epsilon*sigma*4*T_inf^3*(T_inf-T_s), that form is equivalent to
# h*(T_inf-T_s), the convective flux bc. So, the radiative and convective
# flux bcs should give nearly the same answer if the leading terms are equal.
#
[Mesh]
[top]
type = GeneratedMeshGenerator
dim = 3
nx = 10
bias_x = 0.8
ymin = 1.2
ymax = 2.2
boundary_name_prefix = top
[]
[bottom]
type = GeneratedMeshGenerator
dim = 3
nx = 10
bias_x = 0.8
boundary_name_prefix = bot
boundary_id_offset = 6
[]
[two_blocks]
type = MeshCollectionGenerator
inputs = 'top bottom'
[]
[]
[Variables]
[temp]
initial_condition = 600.0
[]
[]
[Kernels]
[heat_dt]
type = ADTimeDerivative
variable = temp
[]
[heat_conduction]
type = ADHeatConduction
variable = temp
[]
[]
[BCs]
[top_right]
type = ADConvectiveHeatFluxBC
variable = temp
boundary = top_right
T_infinity = 300.0
heat_transfer_coefficient = 3.0
[]
[bot_right]
type = ADFunctionRadiativeBC
variable = temp
boundary = bot_right
# htc/(stefan-boltzmann*4*T_inf^3)
emissivity_function = '3/(5.670367e-8*4*300*300*300)'
# Using previous default
Tinfinity = 0
[]
[]
[Materials]
[thermal]
type = ADGenericConstantMaterial
prop_names = 'density thermal_conductivity specific_heat'
prop_values = '1 10 100'
[]
[]
[Postprocessors]
[top_left_temp]
type = SideAverageValue
variable = temp
boundary = top_left
execute_on = 'TIMESTEP_END initial'
[]
[bot_left_temp]
type = SideAverageValue
variable = temp
boundary = bot_left
execute_on = 'TIMESTEP_END initial'
[]
[top_right_temp]
type = SideAverageValue
variable = temp
boundary = top_right
[]
[bot_right_temp]
type = SideAverageValue
variable = temp
boundary = bot_right
[]
[]
[Executioner]
type = Transient
num_steps = 10
dt = 1e1
nl_abs_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(modules/heat_transfer/test/tests/ad_convective_heat_flux/coupled.i)
[Mesh]
type = GeneratedMesh
dim = 3
nx = 10
[]
[Variables]
[./temp]
initial_condition = 200.0
[../]
[]
[Kernels]
[./heat_dt]
type = ADTimeDerivative
variable = temp
[../]
[./heat_conduction]
type = Diffusion
variable = temp
[../]
[./heat]
type = ADBodyForce
variable = temp
value = 0
[../]
[]
[BCs]
[./right]
type = ADConvectiveHeatFluxBC
variable = temp
boundary = 'right'
T_infinity = T_inf
heat_transfer_coefficient = htc
[../]
[]
[Materials]
[chf_mat]
type = ADConvectiveHeatFluxTest
temperature = temp
boundary = 'right'
[]
[]
[Postprocessors]
[./left_temp]
type = SideAverageValue
variable = temp
boundary = left
execute_on = 'TIMESTEP_END initial'
[../]
[./right_temp]
type = SideAverageValue
variable = temp
boundary = right
[../]
[./right_flux]
type = SideDiffusiveFluxAverage
variable = temp
boundary = right
diffusivity = 1
[../]
[]
[Executioner]
type = Transient
num_steps = 10
dt = 1
nl_abs_tol = 1e-12
[]
[Outputs]
[./out]
type = CSV
time_step_interval = 10
[../]
[]
(tutorials/shield_multiphysics/inputs/step11_multiapps/step11_2d_heat_conduction.i)
# Real facility uses forced convection to cool the water tank at full power
# Need to lower power for natural convection so concrete doesn't get too hot.
power = '${fparse 5e4 / 144 * 0.5}'
[Mesh]
[fmg]
type = FileMeshGenerator
file = 'mesh2d_coarse_in.e'
[]
[]
[Variables]
[T]
# Adds a Linear Lagrange variable by default
block = 'concrete_hd concrete Al'
[]
[]
[Kernels]
[diffusion_concrete]
type = ADHeatConduction
variable = T
[]
[]
[Materials]
[concrete_hd]
type = ADHeatConductionMaterial
block = concrete_hd
temp = 'T'
# we specify a function of time, temperature is passed as the time argument
# in the material
thermal_conductivity_temperature_function = '5.0 + 0.001 * t'
[]
[concrete]
type = ADHeatConductionMaterial
block = concrete
temp = 'T'
thermal_conductivity_temperature_function = '2.25 + 0.001 * t'
[]
[Al]
type = ADHeatConductionMaterial
block = Al
temp = T
thermal_conductivity_temperature_function = '175'
[]
[]
[BCs]
[from_reactor]
type = NeumannBC
variable = T
boundary = inner_cavity_solid
# 5 MW reactor, only 50 kW removed from radiation, 144 m2 cavity area
value = '${power}'
[]
[air_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'air_boundary'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 10
[]
[ground]
type = DirichletBC
variable = T
value = 300
boundary = 'ground'
[]
[water_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'water_boundary_inwards'
T_infinity_functor = T_fluid
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient_functor = 600
[]
[]
[Problem]
# No kernels on the water domain
kernel_coverage_check = false
# No materials on the water domain
material_coverage_check = false
[]
[Executioner]
# For pseudo-transient
type = Transient
start_time = -1
end_time = ${units 4 h -> s}
dtmax = 100
[TimeStepper]
type = FunctionDT
function = 'if(t<0.1, 0.1, t)'
[]
# For steady-state fixed-point iteration
# type = Steady
# fixed_point_max_its = 20
# accept_on_max_fixed_point_iteration = true
solve_type = NEWTON # Perform a Newton solve, uses AD to compute Jacobian terms
petsc_options_iname = '-pc_type -pc_hypre_type' # PETSc option pairs with values below
petsc_options_value = 'hypre boomeramg'
nl_abs_tol = 1e-8
[]
[Positions]
[detector_positions]
type = FilePositions
files = detector_positions_2d.txt
[]
[]
[MultiApps]
[fluid]
# For pseudo-transient
type = TransientMultiApp
# For steady-state fixed-point iteration
# type = FullSolveMultiApp
input_files = step11_2d_fluid.i
execute_on = 'TIMESTEP_END'
# Pass in parameter values as if from command line
cli_args = 'power=${power}'
[]
[detectors]
type = FullSolveMultiApp
input_files = 'step11_local.i'
# Create one app at each position
positions_objects = 'detector_positions'
# displace the subapp output to their position in the parent app frame
output_in_position = true
# compute the global temperature first
execute_on = 'TIMESTEP_END'
# Pass in parameter values as if from command line
cli_args = 'Outputs/console=false'
[]
[]
[Transfers]
# transfers solid temperature to nearest node on fluid mesh
[send_T_solid]
type = MultiAppCopyTransfer
to_multi_app = fluid
source_variable = T
variable = T_solid
[]
# Receive fluid temperature
[recv_T_fluid]
type = MultiAppCopyTransfer
from_multi_app = fluid
source_variable = T_fluid
variable = T_fluid
to_blocks = 'water'
from_blocks = 'water'
[]
# transfers local boundary temperature to the each child app
[send_exterior_temperature]
type = MultiAppVariableValueSamplePostprocessorTransfer
to_multi_app = detectors
source_variable = T
postprocessor = T_boundary
[]
# transfers local flux conditions to each child app
[send_local_flux]
type = MultiAppVariableValueSampleTransfer
to_multi_app = detectors
source_variable = flux
variable = flux
[]
# retrieve outputs from the child apps
[hdpe_temperature]
type = MultiAppPostprocessorInterpolationTransfer
from_multi_app = detectors
postprocessor = T_hdpe_inner
variable = T_hdpe_inner
[]
[boron_temperature]
type = MultiAppPostprocessorInterpolationTransfer
from_multi_app = detectors
postprocessor = T_boron_inner
variable = T_boron_inner
[]
[]
[AuxVariables]
[T_fluid]
type = INSFVEnergyVariable
initial_condition = 300
block = 'water'
[]
[flux]
[InitialCondition]
type = FunctionIC
function = '1e4 * exp(-((x-3.25)^2 + (y-2.225)^2))'
[]
[]
# We only output two fields as an example
[T_hdpe_inner]
family = MONOMIAL
order = CONSTANT
block = Al
[]
[T_boron_inner]
family = MONOMIAL
order = CONSTANT
block = Al
[]
[]
[Outputs]
exodus = true
[]
(tutorials/shield_multiphysics/inputs/step13_restart/step13c_restart_from_checkpoint.i)
[Mesh]
[fmg]
type = FileMeshGenerator
file = 'step13a_base_calc_out_cp/LATEST'
[]
[]
[Problem]
# all variables, both nonlinear and auxiliary, are 'restarted'
restart_file_base = 'step13a_base_calc_out_cp/LATEST'
# No kernels on the water domain
kernel_coverage_check = false
# No materials on the water domain
material_coverage_check = false
[]
[Variables]
[T]
# Adds a Linear Lagrange variable by default
block = 'concrete_hd concrete Al'
[]
[]
[Kernels]
[diffusion_concrete]
type = ADHeatConduction
variable = T
[]
[time_derivative]
type = ADHeatConductionTimeDerivative
variable = T
[]
[]
[Materials]
[concrete_hd]
type = ADHeatConductionMaterial
block = concrete_hd
temp = 'T'
# we specify a function of time, temperature is passed as the time argument
# in the material
thermal_conductivity_temperature_function = '5.0 + 0.001 * t'
specific_heat = 1050
[]
[concrete]
type = ADHeatConductionMaterial
block = concrete
temp = 'T'
thermal_conductivity_temperature_function = '2.25 + 0.001 * t'
specific_heat = 1050
[]
[Al]
type = ADHeatConductionMaterial
block = Al
temp = T
thermal_conductivity_temperature_function = '175'
specific_heat = 875
[]
[density_concrete_hd]
type = ADGenericConstantMaterial
block = 'concrete_hd'
prop_names = 'density'
prop_values = '3524' # kg / m3
[]
[density_concrete]
type = ADGenericConstantMaterial
block = 'concrete'
prop_names = 'density'
prop_values = '2403' # kg / m3
[]
[density_Al]
type = ADGenericConstantMaterial
block = 'Al'
prop_names = 'density'
prop_values = '2270' # kg / m3
[]
[]
[BCs]
[from_reactor]
type = NeumannBC
variable = T
boundary = inner_cavity_solid
# 5 MW reactor, only 50 kW removed from radiation, 144 m2 cavity area
value = '${fparse 5e4 / 144}'
[]
[air_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'air_boundary'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 10
[]
[ground]
type = DirichletBC
variable = T
value = 300
boundary = 'ground'
[]
[water_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'water_boundary_inwards'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 600
[]
[]
[Executioner]
type = Transient
num_steps = 6
dt = '${units 12 h -> s}'
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(tutorials/shield_multiphysics/inputs/step08_adaptivity/step8_uniform.i)
[Mesh]
[fmg]
type = FileMeshGenerator
file = '../step03_boundary_conditions/mesh_in.e'
[]
uniform_refine = 0
[]
[Variables]
[T]
# Adds a Linear Lagrange variable by default
block = 'concrete_hd concrete Al'
[]
[]
[Kernels]
[diffusion_concrete]
type = ADHeatConduction
variable = T
[]
[]
[Materials]
[concrete_hd]
type = ADHeatConductionMaterial
block = concrete_hd
temp = 'T'
# we specify a function of time, temperature is passed as the time argument
# in the material
thermal_conductivity_temperature_function = '5.0 + 0.001 * t'
[]
[concrete]
type = ADHeatConductionMaterial
block = concrete
temp = 'T'
thermal_conductivity_temperature_function = '2.25 + 0.001 * t'
[]
[Al]
type = ADHeatConductionMaterial
block = Al
temp = T
thermal_conductivity_temperature_function = '175'
[]
[]
[BCs]
[from_reactor]
type = NeumannBC
variable = T
boundary = inner_cavity_solid
# 5 MW reactor, only 50 kW removed from radiation, 144 m2 cavity area
value = '${fparse 5e4 / 144}'
[]
[air_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'air_boundary'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 10
[]
[ground]
type = DirichletBC
variable = T
value = 300
boundary = 'ground'
[]
[water_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'water_boundary_inwards'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 600
[]
[]
[Problem]
# No kernels on the water domain
kernel_coverage_check = false
# No materials on the water domain
material_coverage_check = false
[]
[Executioner]
type = Steady # Steady state problem
solve_type = NEWTON # Perform a Newton solve, uses AD to compute Jacobian terms
petsc_options_iname = '-pc_type -pc_hypre_type' # PETSc option pairs with values below
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true # Output Exodus format
[]
(tutorials/shield_multiphysics/inputs/step07_mechanics/step7.i)
[Mesh]
[fmg]
type = FileMeshGenerator
file = '../step03_boundary_conditions/mesh_in.e'
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y disp_z'
[]
[Variables]
[T]
# Adds a Linear Lagrange variable by default
block = 'concrete_hd concrete Al'
[]
[]
[Kernels]
[diffusion_concrete]
type = ADHeatConduction
variable = T
[]
[gravity]
type = Gravity
variable = 'disp_z'
value = '-9.81'
block = 'concrete_hd concrete Al'
[]
[]
[Physics/SolidMechanics/QuasiStatic]
[all]
# This block adds all of the proper Kernels, strain calculators, and Variables
# for Solid Mechanics equations in the correct coordinate system (autodetected)
add_variables = true
strain = FINITE
eigenstrain_names = eigenstrain
use_automatic_differentiation = true
generate_output = 'vonmises_stress elastic_strain_xx elastic_strain_yy strain_xx strain_yy'
block = 'concrete_hd concrete Al'
[]
[]
[BCs]
[from_reactor]
type = NeumannBC
variable = T
boundary = inner_cavity_solid
# 5 MW reactor, only 50 kW removed from radiation, 144 m2 cavity area
value = '${fparse 5e4 / 144}'
[]
[air_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'air_boundary'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 10
[]
[ground]
type = DirichletBC
variable = T
value = 300
boundary = 'ground'
[]
[water_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'water_boundary_inwards'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 600
[]
[hold_ground_x]
type = DirichletBC
variable = disp_x
boundary = ground
value = 0
[]
[hold_ground_y]
type = DirichletBC
variable = disp_y
boundary = ground
value = 0
[]
[hold_ground_z]
type = DirichletBC
variable = disp_z
boundary = ground
value = 0
[]
[]
[Materials]
[concrete_hd]
type = ADHeatConductionMaterial
block = concrete_hd
temp = 'T'
# we specify a function of time, temperature is passed as the time argument
# in the material
thermal_conductivity_temperature_function = '5.0 + 0.001 * t'
[]
[concrete]
type = ADHeatConductionMaterial
block = concrete
temp = 'T'
thermal_conductivity_temperature_function = '2.25 + 0.001 * t'
[]
[Al]
type = ADHeatConductionMaterial
block = Al
temp = T
thermal_conductivity_temperature_function = '175'
[]
[elasticity_tensor_concrete_hd]
type = ADComputeIsotropicElasticityTensor
youngs_modulus = 2.75e9 # (Pa)
poissons_ratio = 0.15
block = 'concrete_hd'
[]
[elasticity_tensor_concrete]
type = ADComputeIsotropicElasticityTensor
youngs_modulus = 30e9 # (Pa)
poissons_ratio = 0.2
block = 'concrete'
[]
[elasticity_tensor_Al]
type = ADComputeIsotropicElasticityTensor
youngs_modulus = 68e9 # (Pa)
poissons_ratio = 0.36
block = 'Al'
[]
[elastic_stress]
type = ADComputeFiniteStrainElasticStress
block = 'concrete_hd concrete Al'
[]
[thermal_strain_concrete_hd]
type = ADComputeThermalExpansionEigenstrain
stress_free_temperature = 300
eigenstrain_name = eigenstrain
temperature = T
thermal_expansion_coeff = 1e-5 # 1/K
block = 'concrete_hd'
[]
[thermal_strain_concrete]
type = ADComputeThermalExpansionEigenstrain
stress_free_temperature = 300
eigenstrain_name = eigenstrain
temperature = T
thermal_expansion_coeff = 1e-5 # 1/K
block = 'concrete'
[]
[thermal_strain_Al]
type = ADComputeThermalExpansionEigenstrain
stress_free_temperature = 300 # arbitrary value
eigenstrain_name = eigenstrain
temperature = T
thermal_expansion_coeff = 2.4e-5 # 1/K
block = 'Al'
[]
# NOTE: This handles thermal expansion by coupling to the displacements
[density_concrete_hd]
type = Density
block = 'concrete_hd'
density = '3524' # kg / m3
[]
[density_concrete]
type = Density
block = 'concrete'
density = '2403' # kg / m3
[]
[density_Al]
type = Density
block = 'Al'
density = '2270' # kg / m3
[]
[]
[Problem]
# No kernels on the water domain
kernel_coverage_check = false
# No materials defined on the water domain
material_coverage_check = false
[]
[Executioner]
type = Steady
solve_type = NEWTON
automatic_scaling = true
petsc_options_iname = '-pc_type -pc_hypre_type -ksp_gmres_restart'
petsc_options_value = 'hypre boomeramg 500'
line_search = none
[]
[Outputs]
exodus = true
[]
(tutorials/shield_multiphysics/inputs/step05_auxiliary_variables/step5.i)
[Mesh]
[fmg]
type = FileMeshGenerator
file = '../step03_boundary_conditions/mesh_in.e'
[]
[]
[Variables]
[T]
# Adds a Linear Lagrange variable by default
block = 'concrete_hd concrete Al'
[]
[]
[Kernels]
[diffusion_concrete]
type = ADHeatConduction
variable = T
[]
[]
[Materials]
[concrete_hd]
type = ADHeatConductionMaterial
block = concrete_hd
temp = 'T'
# we specify a function of time, temperature is passed as the time argument
# in the material
thermal_conductivity_temperature_function = '5.0 + 0.001 * t'
[]
[concrete]
type = ADHeatConductionMaterial
block = concrete
temp = 'T'
thermal_conductivity_temperature_function = '2.25 + 0.001 * t'
[]
[Al]
type = ADHeatConductionMaterial
block = Al
temp = T
thermal_conductivity_temperature_function = '175'
[]
[]
[AuxVariables]
[T_fluid]
family = MONOMIAL
order = CONSTANT
block = 'water'
initial_condition = 300
[]
[heat_flux_x]
family = MONOMIAL
order = CONSTANT
block = 'concrete_hd concrete'
[]
[heat_flux_y]
family = MONOMIAL
order = CONSTANT
block = 'concrete_hd concrete'
[]
[heat_flux_z]
family = MONOMIAL
order = CONSTANT
block = 'concrete_hd concrete'
[]
[]
[AuxKernels]
[diff_flux_x]
type = DiffusionFluxAux
variable = heat_flux_x
diffusion_variable = T
diffusivity = 'thermal_conductivity'
component = 'x'
[]
[diff_flux_y]
type = DiffusionFluxAux
variable = heat_flux_y
diffusion_variable = T
diffusivity = 'thermal_conductivity'
component = 'y'
[]
[diff_flux_z]
type = DiffusionFluxAux
variable = heat_flux_z
diffusion_variable = T
diffusivity = 'thermal_conductivity'
component = 'z'
[]
[]
[BCs]
[from_reactor]
type = NeumannBC
variable = T
boundary = inner_cavity_solid
# 5 MW reactor, only 50 kW removed from radiation, 144 m2 cavity area
value = '${fparse 5e4 / 144}'
[]
[air_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'air_boundary'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 10
[]
[ground]
type = DirichletBC
variable = T
value = 300
boundary = 'ground'
[]
[water_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'water_boundary_inwards'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 600
[]
[]
[Problem]
# No kernels on the water domain
kernel_coverage_check = false
# No materials on the water domain
material_coverage_check = false
[]
[Executioner]
type = Steady # Steady state problem
solve_type = NEWTON # Perform a Newton solve, uses AD to compute Jacobian terms
petsc_options_iname = '-pc_type -pc_hypre_type' # PETSc option pairs with values below
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true # Output Exodus format
[]
(tutorials/shield_multiphysics/inputs/step04_heat_conduction/step4.i)
[Mesh]
[fmg]
type = FileMeshGenerator
file = '../step03_boundary_conditions/mesh_in.e'
[]
[]
[Variables]
[T]
# Adds a Linear Lagrange variable by default
block = 'concrete_hd concrete Al'
[]
[]
[Kernels]
[diffusion_concrete]
type = ADHeatConduction
variable = T
[]
[]
[Materials]
[concrete_hd]
type = ADHeatConductionMaterial
block = concrete_hd
temp = 'T'
# we specify a function of time, temperature is passed as the time argument
# in the material
thermal_conductivity_temperature_function = '5.0 + 0.001 * t'
[]
[concrete]
type = ADHeatConductionMaterial
block = concrete
temp = 'T'
thermal_conductivity_temperature_function = '2.25 + 0.001 * t'
[]
[Al]
type = ADHeatConductionMaterial
block = Al
temp = T
thermal_conductivity_temperature_function = '175'
[]
[]
[BCs]
[from_reactor]
type = NeumannBC
variable = T
boundary = inner_cavity_solid
# 5 MW reactor, only 50 kW removed from radiation, 144 m2 cavity area
value = '${fparse 5e4 / 144}'
[]
[air_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'air_boundary'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 10
[]
[ground]
type = DirichletBC
variable = T
value = 300
boundary = 'ground'
[]
[water_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'water_boundary_inwards'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 600
[]
[]
[Problem]
# No kernels on the water domain
kernel_coverage_check = false
# No materials on the water domain
material_coverage_check = false
[]
[Executioner]
type = Steady # Steady state problem
solve_type = NEWTON # Perform a Newton solve, uses AD to compute Jacobian terms
petsc_options_iname = '-pc_type -pc_hypre_type' # PETSc option pairs with values below
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true # Output Exodus format
[]
(modules/heat_transfer/test/tests/ad_convective_heat_flux/fe_fv_coupled.i)
[Mesh]
[cmg]
type = CartesianMeshGenerator
dim = 2
dx = '0.5 0.5'
dy = '0.5 0.5'
ix = '5 5'
iy = '5 5'
subdomain_id = '0 1
0 1'
[]
[add_sideset0]
type = SideSetsBetweenSubdomainsGenerator
input = cmg
new_boundary = middle01
primary_block = 0
paired_block = 1
[]
[add_sideset1]
type = SideSetsBetweenSubdomainsGenerator
input = add_sideset0
new_boundary = middle10
primary_block = 1
paired_block = 0
[]
[]
[Variables]
[u_fe]
block = 0
[]
[u_fv]
type = MooseVariableFVReal
block = 1
[]
[]
[Kernels]
[u_fe_diff]
type = ADDiffusion
variable = u_fe
[]
[]
[BCs]
[u_fe_left]
type = ADDirichletBC
boundary = left
variable = u_fe
value = 0
[]
[u_fe_middle]
type = ADConvectiveHeatFluxBC
boundary = middle01
variable = u_fe
T_infinity_functor = u_fv
heat_transfer_coefficient_functor = 1.0
[]
[]
[FVKernels]
[u_fv_diff]
type = FVDiffusion
variable = u_fv
coeff = 1.0
[]
[]
[FVBCs]
[u_fv_right]
type = FVDirichletBC
boundary = right
variable = u_fv
value = 1.0
[]
[u_fv_middle]
type = FVFunctorConvectiveHeatFluxBC
boundary = middle10
variable = u_fv
T_bulk = u_fv
T_solid = u_fe
heat_transfer_coefficient = 1.0
is_solid = false
[]
[]
[Executioner]
type = Steady
[]
[Outputs]
exodus = true
[]
(modules/heat_transfer/test/tests/ad_convective_heat_flux/equilibrium.i)
[Mesh]
type = GeneratedMesh
dim = 3
nx = 10
[]
[Variables]
[./temp]
initial_condition = 200.0
[../]
[]
[Kernels]
[./heat_dt]
type = ADTimeDerivative
variable = temp
[../]
[./heat_conduction]
type = ADDiffusion
variable = temp
[../]
[]
[BCs]
[./right]
type = ADConvectiveHeatFluxBC
variable = temp
boundary = 'right'
T_infinity = 100.0
heat_transfer_coefficient = 1
[../]
[]
[Postprocessors]
[./left_temp]
type = SideAverageValue
variable = temp
boundary = left
execute_on = 'TIMESTEP_END initial'
[../]
[./right_temp]
type = SideAverageValue
variable = temp
boundary = right
[../]
[./right_flux]
type = SideDiffusiveFluxAverage
variable = temp
boundary = right
diffusivity = 1
[../]
[]
[Executioner]
type = Transient
num_steps = 10
dt = 1e1
nl_abs_tol = 1e-12
[]
[Outputs]
[./out]
type = CSV
time_step_interval = 10
[../]
[]
(modules/heat_transfer/test/tests/ad_convective_heat_flux/flux.i)
# This is a test of the ConvectiveHeatFluxBC.
# There is a single 1x1 element with a prescribed temperature
# on the left side and a convective flux BC on the right side.
# The temperature on the left is 100, and the far-field temp is 200.
# The conductance of the body (conductivity * length) is 10
#
# If the conductance in the BC is also 10, the temperature on the
# right side of the solid element should be 150 because half of the
# temperature drop should occur over the body and half in the BC.
#
# The integrated flux is deltaT * conductance, or -50 * 10 = -500.
# The negative sign indicates that heat is going into the body.
[Mesh]
type = GeneratedMesh
dim = 2
[]
[Problem]
extra_tag_vectors = 'bcs'
[]
[Variables]
[./temp]
initial_condition = 100.0
[../]
[]
[Kernels]
[./heat_conduction]
type = ADHeatConduction
variable = temp
thermal_conductivity = 10
[../]
[]
[BCs]
[./left]
type = ADDirichletBC
variable = temp
boundary = left
value = 100.0
[../]
[./right]
type = ADConvectiveHeatFluxBC
variable = temp
boundary = right
T_infinity = 200.0
heat_transfer_coefficient = 10
[../]
[]
[Postprocessors]
[./right_flux]
type = SideDiffusiveFluxAverage
variable = temp
boundary = right
diffusivity = 10
[../]
[]
[Executioner]
type = Transient
num_steps = 1.0
nl_rel_tol = 1e-12
[]
[Outputs]
csv = true
[]
(tutorials/shield_multiphysics/inputs/step06_transient_heat_conduction/step6_pseudo_transient.i)
cp_multiplier = 1e-6
[Mesh]
[fmg]
type = FileMeshGenerator
file = '../step03_boundary_conditions/mesh_in.e'
[]
[]
[Variables]
[T]
# Adds a Linear Lagrange variable by default
block = 'concrete_hd concrete Al'
initial_condition = 300
[]
[]
[Kernels]
[diffusion_concrete]
type = ADHeatConduction
variable = T
[]
[time_derivative]
type = ADHeatConductionTimeDerivative
variable = T
[]
[]
[Materials]
[concrete_hd]
type = ADHeatConductionMaterial
block = concrete_hd
temp = 'T'
# we specify a function of time, temperature is passed as the time argument
# in the material
thermal_conductivity_temperature_function = '5.0 + 0.001 * t'
specific_heat = ${fparse cp_multiplier * 1050}
[]
[concrete]
type = ADHeatConductionMaterial
block = concrete
temp = 'T'
thermal_conductivity_temperature_function = '2.25 + 0.001 * t'
specific_heat = ${fparse cp_multiplier * 1050}
[]
[Al]
type = ADHeatConductionMaterial
block = Al
temp = T
thermal_conductivity_temperature_function = '175'
specific_heat = ${fparse cp_multiplier * 875}
[]
[density_concrete_hd]
type = ADGenericConstantMaterial
block = 'concrete_hd'
prop_names = 'density'
prop_values = '3524' # kg / m3
[]
[density_concrete]
type = ADGenericConstantMaterial
block = 'concrete'
prop_names = 'density'
prop_values = '2403' # kg / m3
[]
[density_Al]
type = ADGenericConstantMaterial
block = 'Al'
prop_names = 'density'
prop_values = '2270' # kg / m3
[]
[]
[BCs]
[from_reactor]
type = NeumannBC
variable = T
boundary = inner_cavity_solid
# 5 MW reactor, only 50 kW removed from radiation, 144 m2 cavity area
value = '${fparse 5e4 / 144}'
[]
[air_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'air_boundary'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 10
[]
[ground]
type = DirichletBC
variable = T
value = 300
boundary = 'ground'
[]
[water_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'water_boundary_inwards'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 600
[]
[]
[Problem]
# No kernels on the water domain
kernel_coverage_check = false
# No materials on the water domain
material_coverage_check = false
[]
[Executioner]
type = Transient
steady_state_detection = true
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
# Difficult to converge with relative tolerance close to steady-state
nl_abs_tol = 1e-8
[]
[Outputs]
exodus = true
[]
(tutorials/shield_multiphysics/inputs/step13_restart/step13a_base_calc.i)
[Mesh]
[fmg]
type = FileMeshGenerator
file = '../step03_boundary_conditions/mesh_in.e'
[]
[]
[Variables]
[T]
# Adds a Linear Lagrange variable by default
block = 'concrete_hd concrete Al'
initial_condition = 300
[]
[]
[Kernels]
[diffusion_concrete]
type = ADHeatConduction
variable = T
[]
[time_derivative]
type = ADHeatConductionTimeDerivative
variable = T
[]
[]
[Materials]
[concrete_hd]
type = ADHeatConductionMaterial
block = concrete_hd
temp = 'T'
# we specify a function of time, temperature is passed as the time argument
# in the material
thermal_conductivity_temperature_function = '5.0 + 0.001 * t'
specific_heat = 1050
[]
[concrete]
type = ADHeatConductionMaterial
block = concrete
temp = 'T'
thermal_conductivity_temperature_function = '2.25 + 0.001 * t'
specific_heat = 1050
[]
[Al]
type = ADHeatConductionMaterial
block = Al
temp = T
thermal_conductivity_temperature_function = '175'
specific_heat = 875
[]
[density_concrete_hd]
type = ADGenericConstantMaterial
block = 'concrete_hd'
prop_names = 'density'
prop_values = '3524' # kg / m3
[]
[density_concrete]
type = ADGenericConstantMaterial
block = 'concrete'
prop_names = 'density'
prop_values = '2403' # kg / m3
[]
[density_Al]
type = ADGenericConstantMaterial
block = 'Al'
prop_names = 'density'
prop_values = '2270' # kg / m3
[]
[]
[BCs]
[from_reactor]
type = NeumannBC
variable = T
boundary = inner_cavity_solid
# 5 MW reactor, only 50 kW removed from radiation, 144 m2 cavity area
value = '${fparse 5e4 / 144}'
[]
[air_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'air_boundary'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 10
[]
[ground]
type = DirichletBC
variable = T
value = 300
boundary = 'ground'
[]
[water_convection]
type = ADConvectiveHeatFluxBC
variable = T
boundary = 'water_boundary_inwards'
T_infinity = 300.0
# The heat transfer coefficient should be obtained from a correlation
heat_transfer_coefficient = 600
[]
[]
[Problem]
# No kernels on the water domain
kernel_coverage_check = false
# No materials on the water domain
material_coverage_check = false
[]
[Executioner]
type = Transient
num_steps = 4
dt = ${units 12 h -> s}
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
checkpoint = true
[]