- flow_channelName of flow channel component to connect to
C++ Type:std::string
Controllable:No
Description:Name of flow channel component to connect to
- q_wallSpecified wall heat flux [W/m^2]
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:No
Description:Specified wall heat flux [W/m^2]
HeatTransferFromHeatFlux1Phase
This component is a single-phase heat transfer component that uses a wall heat flux from a user-provided function.
Usage
The user must supply the name of the connected flow channel via the parameter "flow_channel".
The parameter "P_hf" is optional and specifies the heated perimeter ; if unspecified, this is computed from the cross-sectional area assuming a circular cross section.
The parameter "Hw" is optional and specifies the heat transfer coefficient ; if unspecified, it is computed using the selected closures. Note that depending on the type of heat transfer and the chosen closures, it may not be relevant.
The parameter "q_wall" specifies the wall heat flux function .
Input Parameters
- HwConvective heat transfer coefficient [W/(m^2-K)]
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:Yes
Description:Convective heat transfer coefficient [W/(m^2-K)]
- P_hfHeat flux perimeter [m]
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:Yes
Description:Heat flux perimeter [m]
- P_hf_transferredFalseIs heat flux perimeter transferred from an external source?
Default:False
C++ Type:bool
Controllable:No
Description:Is heat flux perimeter transferred from an external source?
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:No
Description:Set the enabled status of the MooseObject.
Advanced Parameters
Formulation
In general, a single-phase heat transfer adds some heat flux term to the energy equation:
where is a heat flux at the flow channel wall, and is the heated perimeter.
Input Files
- (modules/thermal_hydraulics/test/tests/problems/natural_circulation/base.i)
- (modules/thermal_hydraulics/test/tests/components/heat_transfer_from_heat_flux_1phase/phy.q_wall_multiple_3eqn.i)
- (modules/thermal_hydraulics/test/tests/components/heat_transfer_base/err.mixed_heat_modes.i)
- (modules/combined/test/tests/subchannel_thm_coupling/THM_SCM_coupling_pump.i)
- (modules/combined/test/tests/subchannel_thm_coupling/THM_SCM_coupling.i)
- (modules/thermal_hydraulics/test/tests/components/heat_transfer_from_heat_flux_1phase/phy.energy_heatflux_ss_1phase.i)
flow_channel
C++ Type:std::string
Controllable:No
Description:Name of flow channel component to connect to
P_hf
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:Yes
Description:Heat flux perimeter [m]
Hw
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:Yes
Description:Convective heat transfer coefficient [W/(m^2-K)]
q_wall
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:No
Description:Specified wall heat flux [W/m^2]
(modules/thermal_hydraulics/test/tests/problems/natural_circulation/base.i)
# Natural circulation loop
#
# The setup consists of 4 connected 1-m pipes, forming a square:
#
# top_pipe
# *--------------* (1,1)
# | |
# | <- <- | | g
# heated_pipe | <- <- | cooled_pipe V
# | <- <- |
# | |
# (0,0) *--------------*
# bottom_pipe
#
# Heating and cooling occurs in the range z = (0.2 m, 0.8 m) with uniform heat fluxes.
[GlobalParams]
gravity_vector = '0 0 -9.81'
length = ${length}
n_elems = ${n_elems}
A = ${area}
initial_T = ${T_ambient}
initial_p = ${p_initial}
initial_vel = 0
fp = fp
closures = closures
f = 0
Hw = ${htc}
rdg_slope_reconstruction = full
scaling_factor_1phase = '1 1 1e-5'
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
emit_on_nan = none
[]
[]
[Closures]
[closures]
type = Closures1PhaseSimple
[]
[]
[Functions]
[heating_flux_fn]
type = PiecewiseConstant
axis = z
x = '0 0.2 0.8'
y = '0 ${fparse power / (S_heated)} 0'
[]
[cooling_flux_fn]
type = PiecewiseConstant
axis = z
x = '0 0.2 0.8'
y = '0 ${fparse -power / (S_cooled)} 0'
[]
[]
[Components]
[heated_pipe]
type = FlowChannel1Phase
position = '0 0 0'
orientation = '0 0 1'
[]
[top_pipe]
type = FlowChannel1Phase
position = '0 0 1'
orientation = '1 0 0'
[]
[cooled_pipe]
type = FlowChannel1Phase
position = '1 0 1'
orientation = '0 0 -1'
[]
[bottom_pipe]
type = FlowChannel1Phase
position = '1 0 0'
orientation = '-1 0 0'
[]
[heating]
type = HeatTransferFromHeatFlux1Phase
flow_channel = 'heated_pipe'
q_wall = heating_flux_fn
[]
[cooling]
type = HeatTransferFromHeatFlux1Phase
flow_channel = 'cooled_pipe'
q_wall = cooling_flux_fn
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
scheme = bdf2
start_time = 0
end_time = 50
[TimeStepper]
type = IterationAdaptiveDT
dt = 0.01
optimal_iterations = 6
iteration_window = 0
growth_factor = 1.2
cutback_factor = 0.8
[]
steady_state_detection = true
petsc_options_iname = '-pc_type'
petsc_options_value = ' lu '
nl_rel_tol = 1e-10
nl_abs_tol = 1e-10
nl_max_its = 15
l_tol = 1e-4
l_max_its = 10
[]
[VectorPostprocessors]
[heated_pipe_vpp]
type = ElementValueSampler
block = 'heated_pipe'
variable = ${output_variables}
sort_by = z
execute_on = 'FINAL'
[]
[top_pipe_vpp]
type = ElementValueSampler
block = 'top_pipe'
variable = ${output_variables}
sort_by = x
execute_on = 'FINAL'
[]
[cooled_pipe_vpp]
type = ElementValueSampler
block = 'cooled_pipe'
variable = ${output_variables}
sort_by = z
execute_on = 'FINAL'
[]
[bottom_pipe_vpp]
type = ElementValueSampler
block = 'bottom_pipe'
variable = ${output_variables}
sort_by = x
execute_on = 'FINAL'
[]
[]
[Outputs]
xml = true
velocity_as_vector = false
execute_on = 'FINAL'
[]
(modules/thermal_hydraulics/test/tests/components/heat_transfer_from_heat_flux_1phase/phy.q_wall_multiple_3eqn.i)
# Tests that energy conservation is satisfied in 1-phase flow when there are
# multiple heat transfer components connected to the same pipe, using specified
# wall heat flux.
#
# This problem has 2 wall heat flux sources, each with differing parameters.
# Solid wall boundary conditions are imposed such that there should be no flow,
# and the solution should be spatially uniform. With no other sources, the
# energy balance is
# (rho*e*A)^{n+1} = (rho*e*A)^n + dt * [(q1*P1) + (q2*P2)]
# Note that spatial integration is dropped here due to spatial uniformity, and
# E has been replaced with e since velocity should be zero.
#
# For the initial conditions
# p = 100 kPa
# T = 300 K
# the density and specific internal energy should be
# rho = 1359.792245 kg/m^3
# e = 1.1320645935e+05 J/kg
#
# With the following heat source parameters:
# q1 = 10 MW/m^2 P1 = 0.2 m
# q2 = 20 MW/m^2 P2 = 0.4 m
# and A = 1 m^2 and dt = 2 s, the new energy solution value should be
# (rho*e*A)^{n+1} = 1359.792245 * 1.1320645935e+05 * 1 + 2 * (10e6 * 0.2 + 20e6 * 0.4)
# = 173937265.50803775 J/m
#
[GlobalParams]
gravity_vector = '0 0 0'
initial_T = 300
initial_p = 100e3
initial_vel = 0
closures = simple_closures
[]
[FluidProperties]
[fp]
type = StiffenedGasFluidProperties
gamma = 2.35
q = -1167e3
q_prime = 0
p_inf = 1.e9
cv = 1816
[]
[]
[Closures]
[simple_closures]
type = Closures1PhaseSimple
[]
[]
[Components]
[pipe]
type = FlowChannel1Phase
fp = fp
position = '0 0 0'
orientation = '1 0 0'
A = 1
f = 0
# length and number of elements should be arbitrary for the test
length = 10
n_elems = 1
[]
[ht1]
type = HeatTransferFromHeatFlux1Phase
flow_channel = pipe
q_wall = 10e6
P_hf = 0.2
Hw = 1
[]
[ht2]
type = HeatTransferFromHeatFlux1Phase
flow_channel = pipe
q_wall = 20e6
P_hf = 0.4
Hw = 1
[]
[left]
type = SolidWall1Phase
input = 'pipe:in'
[]
[right]
type = SolidWall1Phase
input = 'pipe:out'
[]
[]
[Preconditioning]
[preconditioner]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
scheme = 'bdf2'
start_time = 0
dt = 2
num_steps = 1
abort_on_solve_fail = true
solve_type = 'NEWTON'
line_search = 'basic'
nl_rel_tol = 0
nl_abs_tol = 1e-6
nl_max_its = 5
l_tol = 1e-10
l_max_its = 10
[]
[Postprocessors]
[rhoEA_predicted]
type = ElementAverageValue
variable = rhoEA
block = pipe
[]
# This is included to test the naming of heat transfer quantities in the case
# of multiple heat transfers connected to a flow channel. This PP is not used
# in output but just included to ensure that an error does not occur (which is
# the case if the expected material property name does not exist).
# See https://github.com/idaholab/moose/issues/26286.
[q_wall_name_check]
type = ADElementAverageMaterialProperty
mat_prop = 'q_wall:2'
[]
[]
[Outputs]
[out]
type = CSV
show = 'rhoEA_predicted'
execute_on = 'final'
[]
[]
(modules/thermal_hydraulics/test/tests/components/heat_transfer_base/err.mixed_heat_modes.i)
# Tests that an error is thrown if the user specifies a mixture of heat source
# types (temperature and heat flux).
[GlobalParams]
initial_T = 300
initial_p = 100e3
initial_vel = 0
closures = simple_closures
[]
[FluidProperties]
[fp_water]
type = StiffenedGasFluidProperties
gamma = 2.35
cv = 1816.0
q = -1.167e6
p_inf = 1.0e9
q_prime = 0
[]
[]
[Closures]
[simple_closures]
type = Closures1PhaseSimple
[]
[]
[Components]
[pipe]
type = FlowChannel1Phase
fp = fp_water
position = '0 0 0'
orientation = '1 0 0'
A = 1
f = 0
length = 1
n_elems = 1
[]
[ht1]
type = HeatTransferFromHeatFlux1Phase
flow_channel = pipe
q_wall = 1
P_hf = 1
Hw = 1
[]
[ht2]
type = HeatTransferFromSpecifiedTemperature1Phase
flow_channel = pipe
T_wall = 300
P_hf = 1
Hw = 1
[]
[left]
type = SolidWall
input = 'pipe:in'
[]
[right]
type = SolidWall
input = 'pipe:out'
[]
[]
[Preconditioning]
[preconditioner]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu mumps'
[]
[]
[Executioner]
type = Transient
scheme = 'bdf2'
start_time = 0
dt = 1
num_steps = 1
solve_type = 'NEWTON'
line_search = 'basic'
nl_rel_tol = 0
nl_abs_tol = 1e-6
nl_max_its = 5
l_tol = 1e-10
l_max_its = 10
[Quadrature]
type = GAUSS
order = SECOND
[]
[]
(modules/combined/test/tests/subchannel_thm_coupling/THM_SCM_coupling_pump.i)
# THM file based on https://mooseframework.inl.gov/modules/thermal_hydraulics/tutorials/single_phase_flow/step05.html
# Used to loosely couple THM with SCM
# This is a simple closed loop with a pump providing pressure head, core, pressurizer and HX.
# THM sends massflux and temperature at the inlet of the core, and pressure at the outlet of the core
# to subchannel. Subchannel returns total pressure drop of the assembly and total power to THM and THM calculates an
# average friction factor for the core region.
T_in = 583.0 # K
press = 2e5 # Pa
SC_core = 0.0004980799633447909 #m2
# core parameters
core_length = 1. # m
core_n_elems = 1
A_core = 0.005 #dummy
# pipe parameters
pipe_dia = '${units 10. cm -> m}'
A_pipe = '${fparse 0.25 * pi * pipe_dia^2}'
# heat exchanger parameters
hx_dia_inner = '${units 12. cm -> m}'
hx_wall_thickness = '${units 5. mm -> m}'
hx_dia_outer = '${units 50. cm -> m}'
hx_radius_wall = '${fparse hx_dia_inner / 2. + hx_wall_thickness}'
hx_length = 1.5 # m
hx_n_elems = 25
m_dot_sec_in = 1. # kg/s
[GlobalParams]
initial_p = ${press}
initial_vel = 0.0001
initial_T = ${T_in}
initial_vel_x = 0
initial_vel_y = 0
initial_vel_z = 0
gravity_vector = '0 0 0'
rdg_slope_reconstruction = minmod
scaling_factor_1phase = '1 1e-2 1e-4'
scaling_factor_rhoV = 1
scaling_factor_rhouV = 1e-2
scaling_factor_rhovV = 1e-2
scaling_factor_rhowV = 1e-2
scaling_factor_rhoEV = 1e-4
closures = thm_closures
fp = sodium_eos
[]
[Functions]
[q_wall_fn]
type = ParsedFunction
symbol_names = 'core_power length'
symbol_values = 'core_power ${core_length}'
expression = 'core_power/length'
[]
[]
[FluidProperties]
[water]
type = StiffenedGasFluidProperties
gamma = 2.35
cv = 1816.0
q = -1.167e6
p_inf = 1.0e9
q_prime = 0
[]
[sodium_eos]
type = StiffenedGasFluidProperties
gamma = 1.24
cv = 1052.8
q = -2.6292e+05
p_inf = 1.1564e+08
q_prime = 0
mu = 3.222e-04
k = 73.82
[]
[]
[Closures]
[thm_closures]
type = Closures1PhaseTHM
[]
[none_closures]
type = Closures1PhaseNone
[]
[]
[Materials]
[f_mat]
type = ADParsedMaterial
property_name = f_D
postprocessor_names = 'core_f'
expression = 'core_f'
block = 'core_chan'
[]
[]
[HeatStructureMaterials]
[steel]
type = SolidMaterialProperties
rho = 8050
k = 45
cp = 466
[]
[]
[Components]
[up_pipe_1]
type = FlowChannel1Phase
position = '0 0 -0.5'
orientation = '0 0 1'
length = 0.5
n_elems = 15
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[jct1]
type = JunctionParallelChannels1Phase
position = '0 0 0'
connections = 'up_pipe_1:out core_chan:in'
volume = 1e-5
[]
[core_chan]
type = FlowChannel1Phase
position = '0 0 0'
orientation = '0 0 1'
length = ${core_length}
n_elems = ${core_n_elems}
A = ${A_core}
closures = none_closures
[]
[core_ht]
type = HeatTransferFromHeatFlux1Phase
flow_channel = core_chan
q_wall = q_wall_fn
P_hf = 1
[]
[jct2]
type = JunctionParallelChannels1Phase
position = '0 0 1'
connections = 'core_chan:out up_pipe_2:in'
volume = 1e-5
[]
[up_pipe_2]
type = FlowChannel1Phase
position = '0 0 1'
orientation = '0 0 1'
length = 0.5
n_elems = 10
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[jct3]
type = JunctionOneToOne1Phase
connections = 'up_pipe_2:out top_pipe_1:in'
[]
[top_pipe_1]
type = FlowChannel1Phase
position = '0 0 1.5'
orientation = '1 0 0'
length = 0.5
n_elems = 10
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[top_pipe_2]
type = FlowChannel1Phase
position = '0.5 0 1.5'
orientation = '1 0 0'
length = 0.5
n_elems = 10
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[jct4]
type = VolumeJunction1Phase
position = '0.5 0 1.5'
volume = 1e-5
connections = 'top_pipe_1:out top_pipe_2:in press_pipe:in'
[]
[press_pipe]
type = FlowChannel1Phase
position = '0.5 0 1.5'
orientation = '0 1 0'
length = 0.2
n_elems = 5
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[pressurizer]
type = InletStagnationPressureTemperature1Phase
p0 = ${press}
T0 = 580
input = press_pipe:out
[]
[jct5]
type = JunctionOneToOne1Phase
connections = 'top_pipe_2:out down_pipe_1:in'
[]
[down_pipe_1]
type = FlowChannel1Phase
position = '1 0 1.5'
orientation = '0 0 -1'
length = 0.25
A = ${A_pipe}
n_elems = 5
[]
[jct6]
type = JunctionParallelChannels1Phase
position = '1 0 1.25'
connections = 'down_pipe_1:out hx/pri:in'
volume = 1e-5
[]
[hx]
[pri]
type = FlowChannel1Phase
position = '1 0 1.25'
orientation = '0 0 -1'
length = ${hx_length}
n_elems = ${hx_n_elems}
roughness = 1e-5
A = '${fparse pi * hx_dia_inner * hx_dia_inner / 4.}'
D_h = ${hx_dia_inner}
[]
[ht_pri]
type = HeatTransferFromHeatStructure1Phase
hs = hx/wall
hs_side = inner
flow_channel = hx/pri
P_hf = '${fparse pi * hx_dia_inner}'
[]
[wall]
type = HeatStructureCylindrical
position = '1 0 1.25'
orientation = '0 0 -1'
length = ${hx_length}
n_elems = ${hx_n_elems}
widths = '${hx_wall_thickness}'
n_part_elems = '3'
materials = 'steel'
names = '0'
inner_radius = '${fparse hx_dia_inner / 2.}'
[]
[ht_sec]
type = HeatTransferFromHeatStructure1Phase
hs = hx/wall
hs_side = outer
flow_channel = hx/sec
P_hf = '${fparse 2 * pi * hx_radius_wall}'
[]
[sec]
type = FlowChannel1Phase
position = '${fparse 1 + hx_wall_thickness} 0 -0.25'
orientation = '0 0 1'
length = ${hx_length}
n_elems = ${hx_n_elems}
A = '${fparse pi * (hx_dia_outer * hx_dia_outer / 4. - hx_radius_wall * hx_radius_wall)}'
D_h = '${fparse hx_dia_outer - (2 * hx_radius_wall)}'
fp = water
initial_T = 300
[]
[]
[jct7]
type = JunctionParallelChannels1Phase
position = '1 0 -0.25'
connections = 'hx/pri:out down_pipe_2:in'
volume = 1e-5
[]
[down_pipe_2]
type = FlowChannel1Phase
position = '1 0 -0.25'
orientation = '0 0 -1'
length = 0.25
n_elems = 10
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[jct8]
type = JunctionOneToOne1Phase
connections = 'down_pipe_2:out bottom_1:in'
[]
[bottom_1]
type = FlowChannel1Phase
position = '1 0 -0.5'
orientation = '-1 0 0'
length = 0.5
n_elems = 5
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[pump]
type = Pump1Phase
position = '0.5 0 -0.5'
connections = 'bottom_1:out bottom_2:in'
volume = 1e-4
A_ref = ${A_pipe}
head = 3.56
[]
[bottom_2]
type = FlowChannel1Phase
position = '0.5 0 -0.5'
orientation = '-1 0 0'
length = 0.5
n_elems = 5
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[jct9]
type = JunctionOneToOne1Phase
connections = 'bottom_2:out up_pipe_1:in'
[]
[inlet_sec]
type = InletMassFlowRateTemperature1Phase
input = 'hx/sec:in'
m_dot = ${m_dot_sec_in}
T = 300
[]
[outlet_sec]
type = Outlet1Phase
input = 'hx/sec:out'
p = 1e5
[]
[]
[Postprocessors]
[power_to_coolant]
type = ADHeatRateDirectFlowChannel
q_wall_prop = q_wall
block = core_chan
P_hf = 1
[]
[core_T_out]
type = SideAverageValue
boundary = core_chan:out
variable = T
[]
[T_out]
type = SideAverageValue
boundary = bottom_1:out
variable = T
[]
[core_p_in]
type = SideAverageValue
boundary = up_pipe_1:out
variable = p
[]
[core_p_out]
type = SideAverageValue
boundary = up_pipe_2:in
variable = p
[]
[core_delta_p]
type = ParsedPostprocessor
pp_names = 'core_p_in core_p_out'
function = 'core_p_in - core_p_out'
[]
[hx_pri_T_out]
type = SideAverageValue
boundary = hx/pri:out
variable = T
[]
[hx_sec_T_in]
type = SideAverageValue
boundary = inlet_sec
variable = T
[]
[hx_sec_T_out]
type = SideAverageValue
boundary = outlet_sec
variable = T
[]
[m_dot_sec]
type = ADFlowBoundaryFlux1Phase
boundary = inlet_sec
equation = mass
[]
############## Friction Factor Calculation #############
[av_rhouA]
type = ElementAverageValue
variable = 'rhouA'
block = 'core_chan'
[]
[av_rho]
type = ElementAverageValue
variable = 'rho'
block = 'core_chan'
[]
[Kloss]
type = ParsedPostprocessor
pp_names = 'core_delta_p_tgt av_rhouA av_rho'
function = '2.0 * core_delta_p_tgt * av_rho * ${A_core} * ${A_core} / (av_rhouA * av_rhouA)'
[]
[Dh]
type = ADElementAverageMaterialProperty
mat_prop = D_h
block = core_chan
[]
[core_f]
type = ParsedPostprocessor
pp_names = 'Kloss Dh'
function = 'Kloss * Dh / ${core_length}'
[]
### INFO to send to SC
[outlet_pressure]
type = SideAverageValue
boundary = up_pipe_2:in
variable = p
[]
[inlet_mass_flow_rate]
type = ADFlowJunctionFlux1Phase
boundary = up_pipe_1:out
connection_index = 0
equation = mass
junction = jct1
[]
[inlet_temperature]
type = SideAverageValue
boundary = up_pipe_1:out
variable = T
[]
[inlet_mass_flux]
type = ParsedPostprocessor
pp_names = 'inlet_mass_flow_rate'
function = 'abs(inlet_mass_flow_rate/${SC_core})'
[]
#####
##### Info received from subchannel
[core_delta_p_tgt]
type = Receiver
default = 100
[]
[core_power]
type = Receiver
default = 100
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
start_time = 0
[TimeStepper]
type = IterationAdaptiveDT
dt = 2
[]
dtmax = 50
end_time = 10
line_search = basic
solve_type = NEWTON
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_rel_tol = 1e-7
nl_abs_tol = 1e-7
nl_max_its = 25
fixed_point_min_its = 1
fixed_point_max_its = 5
accept_on_max_fixed_point_iteration = true
auto_advance = true
relaxation_factor = 0.5
[]
[Outputs]
csv = true
[console]
type = Console
max_rows = 1
outlier_variable_norms = false
[]
print_linear_residuals = false
[]
################################################################################
# A multiapp that couples THM to subchannel
################################################################################
[MultiApps]
# active = ''
[subchannel]
type = FullSolveMultiApp
input_files = 'subchannel.i'
execute_on = 'timestep_end'
positions = '0 0 0'
max_procs_per_app = 1
output_in_position = true
bounding_box_padding = '0 0 0.1'
[]
[]
[Transfers]
# active = ''
[pressure_drop_transfer] # Get pressure drop to THM from subchannel
type = MultiAppPostprocessorTransfer
from_multi_app = subchannel
from_postprocessor = total_pressure_drop_SC
to_postprocessor = core_delta_p_tgt
reduction_type = average
execute_on = 'timestep_end'
[]
[power_transfer] # Get Total power to THM from subchannel
type = MultiAppPostprocessorTransfer
from_multi_app = subchannel
from_postprocessor = Total_power
to_postprocessor = core_power
reduction_type = average
execute_on = 'timestep_end'
[]
[mass_flux_tranfer] # Send mass_flux at the inlet of THM core to subchannel
type = MultiAppPostprocessorTransfer
to_multi_app = subchannel
from_postprocessor = inlet_mass_flux
to_postprocessor = report_mass_flux_inlet
execute_on = 'timestep_end'
[]
[outlet_pressure_tranfer] # Send pressure at the outlet of THM core to subchannel
type = MultiAppPostprocessorTransfer
to_multi_app = subchannel
from_postprocessor = outlet_pressure
to_postprocessor = report_pressure_outlet
execute_on = 'timestep_end'
[]
[inlet_temperature_transfer]
type = MultiAppPostprocessorTransfer
to_multi_app = subchannel
from_postprocessor = inlet_temperature
to_postprocessor = report_temperature_inlet
execute_on = 'timestep_end'
[]
[]
(modules/combined/test/tests/subchannel_thm_coupling/THM_SCM_coupling.i)
# THM file based on https://mooseframework.inl.gov/modules/thermal_hydraulics/tutorials/single_phase_flow/step05.html
# Used to loosely couple THM with SCM
# This is a simple open loop with fixed massflow at the inlet and pressure at the outlet.
# THM sends massflux and temperature at the inlet of the core, and pressure at the outlet of the core
# to subchannel. Subchannel returns total pressure drop of the assembly and total power to THM and THM calculates an
# average friction factor for the core region.
T_in = 583.0 # K
m_dot_in = 1 # kg/s
press = 2e5 # Pa
SC_core = 0.0004980799633447909 #m2
# core parameters
core_length = 1. # m
core_n_elems = 1
A_core = 0.005 #dummy
# pipe parameters
pipe_dia = '${units 10. cm -> m}'
A_pipe = '${fparse 0.25 * pi * pipe_dia^2}'
# heat exchanger parameters
hx_dia_inner = '${units 12. cm -> m}'
hx_wall_thickness = '${units 5. mm -> m}'
hx_dia_outer = '${units 50. cm -> m}'
hx_radius_wall = '${fparse hx_dia_inner / 2. + hx_wall_thickness}'
hx_length = 1.5 # m
hx_n_elems = 25
m_dot_sec_in = 1. # kg/s
[GlobalParams]
initial_p = ${press}
initial_vel = 0.0001
initial_T = ${T_in}
initial_vel_x = 0
initial_vel_y = 0
initial_vel_z = 0
gravity_vector = '0 0 0'
rdg_slope_reconstruction = minmod
scaling_factor_1phase = '1 1e-2 1e-4'
scaling_factor_rhoV = 1
scaling_factor_rhouV = 1e-2
scaling_factor_rhovV = 1e-2
scaling_factor_rhowV = 1e-2
scaling_factor_rhoEV = 1e-4
closures = thm_closures
fp = sodium_eos
[]
[Functions]
[q_wall_fn]
type = ParsedFunction
symbol_names = 'core_power length'
symbol_values = 'core_power ${core_length}'
expression = 'core_power/length'
[]
[]
[FluidProperties]
[water]
type = StiffenedGasFluidProperties
gamma = 2.35
cv = 1816.0
q = -1.167e6
p_inf = 1.0e9
q_prime = 0
[]
[sodium_eos]
type = StiffenedGasFluidProperties
gamma = 1.24
cv = 1052.8
q = -2.6292e+05
p_inf = 1.1564e+08
q_prime = 0
mu = 3.222e-04
k = 73.82
[]
[]
[Closures]
[thm_closures]
type = Closures1PhaseTHM
[]
[none_closures]
type = Closures1PhaseNone
[]
[]
[Materials]
[f_mat]
type = ADParsedMaterial
property_name = f_D
postprocessor_names = 'core_f'
expression = 'core_f'
block = 'core_chan'
[]
[]
[HeatStructureMaterials]
[steel]
type = SolidMaterialProperties
rho = 8050
k = 45
cp = 466
[]
[]
[Components]
[inlet]
type = InletMassFlowRateTemperature1Phase
input = 'bottom_2:in'
m_dot = ${m_dot_in}
T = ${T_in}
[]
[outlet]
type = Outlet1Phase
input = 'bottom_1:out'
p = ${press}
[]
[up_pipe_1]
type = FlowChannel1Phase
position = '0 0 -0.5'
orientation = '0 0 1'
length = 0.5
n_elems = 15
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[jct1]
type = JunctionParallelChannels1Phase
position = '0 0 0'
connections = 'up_pipe_1:out core_chan:in'
volume = 1e-5
[]
[core_chan]
type = FlowChannel1Phase
position = '0 0 0'
orientation = '0 0 1'
length = ${core_length}
n_elems = ${core_n_elems}
A = ${A_core}
closures = none_closures
[]
[core_ht]
type = HeatTransferFromHeatFlux1Phase
flow_channel = core_chan
q_wall = q_wall_fn
P_hf = 1
[]
[jct2]
type = JunctionParallelChannels1Phase
position = '0 0 1'
connections = 'core_chan:out up_pipe_2:in'
volume = 1e-5
[]
[up_pipe_2]
type = FlowChannel1Phase
position = '0 0 1'
orientation = '0 0 1'
length = 0.5
n_elems = 10
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[jct3]
type = JunctionOneToOne1Phase
connections = 'up_pipe_2:out top_pipe_1:in'
[]
[top_pipe_1]
type = FlowChannel1Phase
position = '0 0 1.5'
orientation = '1 0 0'
length = 0.5
n_elems = 10
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[top_pipe_2]
type = FlowChannel1Phase
position = '0.5 0 1.5'
orientation = '1 0 0'
length = 0.5
n_elems = 10
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[jct4]
type = VolumeJunction1Phase
position = '0.5 0 1.5'
volume = 1e-5
connections = 'top_pipe_1:out top_pipe_2:in'
[]
[jct5]
type = JunctionOneToOne1Phase
connections = 'top_pipe_2:out down_pipe_1:in'
[]
[down_pipe_1]
type = FlowChannel1Phase
position = '1 0 1.5'
orientation = '0 0 -1'
length = 0.25
A = ${A_pipe}
n_elems = 5
[]
[jct6]
type = JunctionParallelChannels1Phase
position = '1 0 1.25'
connections = 'down_pipe_1:out hx/pri:in'
volume = 1e-5
[]
[hx]
[pri]
type = FlowChannel1Phase
position = '1 0 1.25'
orientation = '0 0 -1'
length = ${hx_length}
n_elems = ${hx_n_elems}
roughness = 1e-5
A = '${fparse pi * hx_dia_inner * hx_dia_inner / 4.}'
D_h = ${hx_dia_inner}
[]
[ht_pri]
type = HeatTransferFromHeatStructure1Phase
hs = hx/wall
hs_side = inner
flow_channel = hx/pri
P_hf = '${fparse pi * hx_dia_inner}'
[]
[wall]
type = HeatStructureCylindrical
position = '1 0 1.25'
orientation = '0 0 -1'
length = ${hx_length}
n_elems = ${hx_n_elems}
widths = '${hx_wall_thickness}'
n_part_elems = '3'
materials = 'steel'
names = '0'
inner_radius = '${fparse hx_dia_inner / 2.}'
[]
[ht_sec]
type = HeatTransferFromHeatStructure1Phase
hs = hx/wall
hs_side = outer
flow_channel = hx/sec
P_hf = '${fparse 2 * pi * hx_radius_wall}'
[]
[sec]
type = FlowChannel1Phase
position = '${fparse 1 + hx_wall_thickness} 0 -0.25'
orientation = '0 0 1'
length = ${hx_length}
n_elems = ${hx_n_elems}
A = '${fparse pi * (hx_dia_outer * hx_dia_outer / 4. - hx_radius_wall * hx_radius_wall)}'
D_h = '${fparse hx_dia_outer - (2 * hx_radius_wall)}'
fp = water
initial_T = 300
[]
[]
[jct7]
type = JunctionParallelChannels1Phase
position = '1 0 -0.25'
connections = 'hx/pri:out down_pipe_2:in'
volume = 1e-5
[]
[down_pipe_2]
type = FlowChannel1Phase
position = '1 0 -0.25'
orientation = '0 0 -1'
length = 0.25
n_elems = 10
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[jct8]
type = JunctionOneToOne1Phase
connections = 'down_pipe_2:out bottom_1:in'
[]
[bottom_1]
type = FlowChannel1Phase
position = '1 0 -0.5'
orientation = '-1 0 0'
length = 0.5
n_elems = 5
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[bottom_2]
type = FlowChannel1Phase
position = '0.5 0 -0.5'
orientation = '-1 0 0'
length = 0.5
n_elems = 5
A = ${A_pipe}
D_h = ${pipe_dia}
[]
[jct9]
type = JunctionOneToOne1Phase
connections = 'bottom_2:out up_pipe_1:in'
[]
[inlet_sec]
type = InletMassFlowRateTemperature1Phase
input = 'hx/sec:in'
m_dot = ${m_dot_sec_in}
T = 300
[]
[outlet_sec]
type = Outlet1Phase
input = 'hx/sec:out'
p = 1e5
[]
[]
[Postprocessors]
[power_to_coolant]
type = ADHeatRateDirectFlowChannel
q_wall_prop = q_wall
block = core_chan
P_hf = 1
[]
[core_T_out]
type = SideAverageValue
boundary = core_chan:out
variable = T
[]
[T_out]
type = SideAverageValue
boundary = bottom_1:out
variable = T
[]
[core_p_in]
type = SideAverageValue
boundary = up_pipe_1:out
variable = p
[]
[core_p_out]
type = SideAverageValue
boundary = up_pipe_2:in
variable = p
[]
[core_delta_p]
type = ParsedPostprocessor
pp_names = 'core_p_in core_p_out'
function = 'core_p_in - core_p_out'
[]
[hx_pri_T_out]
type = SideAverageValue
boundary = hx/pri:out
variable = T
[]
[hx_sec_T_in]
type = SideAverageValue
boundary = inlet_sec
variable = T
[]
[hx_sec_T_out]
type = SideAverageValue
boundary = outlet_sec
variable = T
[]
[m_dot_sec]
type = ADFlowBoundaryFlux1Phase
boundary = inlet_sec
equation = mass
[]
############## Friction Factor Calculation #############
[av_rhouA]
type = ElementAverageValue
variable = 'rhouA'
block = 'core_chan'
[]
[av_rho]
type = ElementAverageValue
variable = 'rho'
block = 'core_chan'
[]
[Kloss]
type = ParsedPostprocessor
pp_names = 'core_delta_p_tgt av_rhouA av_rho'
function = '2.0 * core_delta_p_tgt * av_rho * ${A_core} * ${A_core} / (av_rhouA * av_rhouA)'
[]
[Dh]
type = ADElementAverageMaterialProperty
mat_prop = D_h
block = core_chan
[]
[core_f]
type = ParsedPostprocessor
pp_names = 'Kloss Dh'
function = 'Kloss * Dh / ${core_length}'
[]
### INFO to send to SC
[outlet_pressure]
type = SideAverageValue
boundary = up_pipe_2:in
variable = p
[]
[inlet_mass_flow_rate]
type = ADFlowJunctionFlux1Phase
boundary = up_pipe_1:out
connection_index = 0
equation = mass
junction = jct1
[]
[inlet_temperature]
type = SideAverageValue
boundary = up_pipe_1:out
variable = T
[]
[inlet_mass_flux]
type = ParsedPostprocessor
pp_names = 'inlet_mass_flow_rate'
function = 'abs(inlet_mass_flow_rate/${SC_core})'
[]
#####
##### Info received from subchannel
[core_delta_p_tgt]
type = Receiver
default = 100
[]
[core_power]
type = Receiver
default = 100
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
start_time = 0
[TimeStepper]
type = IterationAdaptiveDT
dt = 0.01
[]
# dtmax = 5
end_time = 5
line_search = basic
solve_type = NEWTON
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_rel_tol = 1e-8
nl_abs_tol = 1e-8
nl_max_its = 25
[]
[Outputs]
csv = true
[console]
type = Console
max_rows = 1
outlier_variable_norms = false
[]
print_linear_residuals = false
[]
################################################################################
# A multiapp that couples THM to subchannel
################################################################################
[MultiApps]
[subchannel]
type = FullSolveMultiApp
input_files = 'subchannel.i'
execute_on = 'timestep_end'
positions = '0 0 0'
max_procs_per_app = 1
output_in_position = true
bounding_box_padding = '0 0 0.1'
[]
[]
[Transfers]
[pressure_drop_transfer] # Get pressure drop to THM from subchannel
type = MultiAppPostprocessorTransfer
from_multi_app = subchannel
from_postprocessor = total_pressure_drop_SC
to_postprocessor = core_delta_p_tgt
reduction_type = average
execute_on = 'timestep_end'
[]
[power_transfer] # Get Total power to THM from subchannel
type = MultiAppPostprocessorTransfer
from_multi_app = subchannel
from_postprocessor = Total_power
to_postprocessor = core_power
reduction_type = average
execute_on = 'timestep_end'
[]
[mass_flux_tranfer] # Send mass_flux at the inlet of THM core to subchannel
type = MultiAppPostprocessorTransfer
to_multi_app = subchannel
from_postprocessor = inlet_mass_flux
to_postprocessor = report_mass_flux_inlet
execute_on = 'timestep_end'
[]
[outlet_pressure_tranfer] # Send pressure at the outlet of THM core to subchannel
type = MultiAppPostprocessorTransfer
to_multi_app = subchannel
from_postprocessor = outlet_pressure
to_postprocessor = report_pressure_outlet
execute_on = 'timestep_end'
[]
[inlet_temperature_transfer]
type = MultiAppPostprocessorTransfer
to_multi_app = subchannel
from_postprocessor = inlet_temperature
to_postprocessor = report_temperature_inlet
execute_on = 'timestep_end'
[]
[]
(modules/thermal_hydraulics/test/tests/components/heat_transfer_from_heat_flux_1phase/phy.energy_heatflux_ss_1phase.i)
# This test tests conservation of energy at steady state for 1-phase flow when a
# heat flux is specified. Conservation is checked by comparing the integral of
# the heat flux against the difference of the boundary fluxes.
[GlobalParams]
initial_p = 7.0e6
initial_vel = 0
initial_T = 513
gravity_vector = '0.0 0.0 0.0'
closures = simple_closures
[]
[FluidProperties]
[eos]
type = StiffenedGasFluidProperties
gamma = 2.35
q = -1167e3
q_prime = 0
p_inf = 1.e9
cv = 1816
[]
[]
[Closures]
[simple_closures]
type = Closures1PhaseSimple
[]
[]
[Components]
[pipe]
type = FlowChannel1Phase
position = '0 0 0'
orientation = '0 0 1'
length = 3.66
n_elems = 10
A = 1.907720E-04
D_h = 1.698566E-02
f = 0.0
fp = eos
[]
[ht_pipe]
type = HeatTransferFromHeatFlux1Phase
flow_channel = pipe
q_wall = 1.0e5
Hw = 1.0e4
P_hf = 4.4925e-2
[]
[inlet]
type = SolidWall1Phase
input = 'pipe:in'
[]
[outlet]
type = SolidWall1Phase
input = 'pipe:out'
[]
[]
[Postprocessors]
[E]
type = ElementIntegralVariablePostprocessor
variable = rhoEA
execute_on = 'initial timestep_end'
[]
[E_change]
type = ChangeOverTimePostprocessor
postprocessor = E
execute_on = 'initial timestep_end'
[]
[]
[Preconditioning]
[SMP_PJFNK]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
abort_on_solve_fail = true
dt = 1
solve_type = 'NEWTON'
line_search = 'basic'
nl_rel_tol = 1e-10
nl_abs_tol = 1e-7
nl_max_its = 50
l_tol = 1e-3
l_max_its = 60
start_time = 0
num_steps = 10
[]
[Outputs]
[out]
type = CSV
show = 'E_change'
[]
[console]
type = Console
show = 'E_change'
[]
[]