- inputName of the input
C++ Type:BoundaryName
Unit:(no unit assumed)
Controllable:No
Description:Name of the input
FreeBoundary1Phase
This is a single-phase 1-D flow boundary component in which no boundary data is supplied by the user. This component should be used only if the boundary is a supersonic exit, since this is the only situation in which no external characteristics enter the domain.
Usage
This component must be connected to a FlowChannel1Phase. See how to connect a flow boundary component.
Input Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
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
Unit:(no unit assumed)
Controllable:No
Description:Set the enabled status of the MooseObject.
Advanced Parameters
Input Files
- (modules/thermal_hydraulics/test/tests/problems/abrupt_area_change_liquid/without_junction.i)
- (modules/thermal_hydraulics/test/tests/problems/abrupt_area_change_liquid/with_junction.i)
- (modules/thermal_hydraulics/test/tests/components/free_boundary_1phase/phy.conservation_free_boundary_1phase.i)
- (modules/thermal_hydraulics/test/tests/components/junction_one_to_one_1phase/junction_one_to_one_1phase.i)
- (modules/thermal_hydraulics/test/tests/problems/lax_shock_tube/lax_shock_tube.i)
- (modules/thermal_hydraulics/test/tests/components/junction_one_to_one_1phase/no_junction_1phase.i)
- (modules/thermal_hydraulics/test/tests/problems/sod_shock_tube/sod_shock_tube.i)
- (modules/thermal_hydraulics/test/tests/problems/square_wave/square_wave.i)
- (modules/thermal_hydraulics/test/tests/problems/sedov_blast_wave/sedov_blast_wave.i)
- (modules/thermal_hydraulics/test/tests/problems/super_sonic_tube/test.i)
- (modules/thermal_hydraulics/test/tests/problems/double_rarefaction/1phase.i)
Formulation
In this formulation, no boundary data is supplied; thus, the boundary flux is evaluated entirely using the interior solution:
(modules/thermal_hydraulics/test/tests/problems/abrupt_area_change_liquid/without_junction.i)
# Version without junction
!include base_params.i
!include base.i
[Functions]
[A_fn]
type = PiecewiseConstant
axis = x
x = '0 ${xR}'
y = '${AL} ${AR}'
[]
[]
[Components]
[pipe]
type = FlowChannel1Phase
position = '0 0 0'
orientation = '1 0 0'
length = ${fparse lengthL + lengthR}
n_elems = ${fparse NL + NR}
A = A_fn
[]
[left_boundary]
type = FreeBoundary1Phase
input = 'pipe:in'
[]
[right_boundary]
type = FreeBoundary1Phase
input = 'pipe:out'
[]
[]
[VectorPostprocessors]
[vpp]
type = ADSampler1DReal
block = 'pipe'
property = 'p vel'
sort_by = x
execute_on = 'FINAL'
[]
[]
[Outputs]
file_base = 'without_junction'
[]
(modules/thermal_hydraulics/test/tests/problems/abrupt_area_change_liquid/with_junction.i)
# Version with junction
!include base_params.i
# For equivalent comparison to the unsplit case, we reduce the length of the
# right pipe by one element
dx = ${fparse lengthL / NL}
NR_minus_junction = ${fparse NR - 1}
lengthR_minus_junction = ${fparse lengthR - dx}
xR_minus_junction = ${fparse xR + dx}
xJ = ${fparse lengthL + 0.5 * dx}
AJ = ${fparse AL + AR}
RJ = ${fparse sqrt(AJ / (4 * pi))} # A = 4 pi R^2
VJ = ${fparse 4/3 * pi * RJ^3}
!include base.i
[Components]
[pipeL]
type = FlowChannel1Phase
position = '0 0 0'
orientation = '1 0 0'
length = ${lengthL}
n_elems = ${NL}
A = ${AL}
[]
[pipeR]
type = FlowChannel1Phase
position = '${xR_minus_junction} 0 0'
orientation = '1 0 0'
length = ${lengthR_minus_junction}
n_elems = ${NR_minus_junction}
A = ${AR}
[]
[junction]
type = VolumeJunction1Phase
connections = 'pipeL:out pipeR:in'
position = '${xJ} 0 0'
volume = ${VJ}
initial_vel_x = 0
initial_vel_y = 0
initial_vel_z = 0
scaling_factor_rhoEV = 1e-5
apply_velocity_scaling = true
[]
[left_boundary]
type = FreeBoundary1Phase
input = 'pipeL:in'
[]
[right_boundary]
type = FreeBoundary1Phase
input = 'pipeR:out'
[]
[]
[VectorPostprocessors]
[vpp]
type = ADSampler1DReal
block = 'pipeL pipeR'
property = 'p vel'
sort_by = x
execute_on = 'FINAL'
[]
[]
[Outputs]
file_base = 'with_junction'
[]
(modules/thermal_hydraulics/test/tests/components/free_boundary_1phase/phy.conservation_free_boundary_1phase.i)
# This test tests conservation of mass, momentum, and energy on a transient
# problem with an inlet and outlet (using free boundaries for each). This test
# takes 1 time step with Crank-Nicolson and some boundary flux integral
# post-processors needed for the full conservation statement. Lastly, the
# conservation quantities are shown on the console, which should ideally be zero
# for full conservation.
[GlobalParams]
gravity_vector = '0 0 0'
scaling_factor_1phase = '1 1 1e-6'
closures = simple_closures
[]
[Functions]
[T_fn]
type = ParsedFunction
expression = '300 + 10 * (cos(2*pi*x + pi))'
[]
[]
[FluidProperties]
[fp]
type = StiffenedGasFluidProperties
gamma = 2.35
cv = 1816.0
q = -1.167e6
p_inf = 1.0e9
q_prime = 0
[]
[]
[Closures]
[simple_closures]
type = Closures1PhaseSimple
[]
[]
[Components]
[inlet]
type = FreeBoundary1Phase
input = pipe:in
[]
[pipe]
type = FlowChannel1Phase
position = '0 0 0'
orientation = '1 0 0'
length = 1.0
n_elems = 10
A = 1.0
initial_T = T_fn
initial_p = 1e5
initial_vel = 1
f = 0
fp = fp
[]
[outlet]
type = FreeBoundary1Phase
input = pipe:out
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
scheme = crank-nicolson
start_time = 0.0
end_time = 0.01
dt = 0.01
abort_on_solve_fail = true
solve_type = 'PJFNK'
line_search = 'basic'
nl_rel_tol = 1e-6
nl_abs_tol = 1e-4
nl_max_its = 10
l_tol = 1e-2
l_max_its = 20
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
[]
[Postprocessors]
# MASS
[massflux_left]
type = MassFluxIntegral
boundary = inlet
arhouA = rhouA
[]
[massflux_right]
type = MassFluxIntegral
boundary = outlet
arhouA = rhouA
[]
[massflux_difference]
type = DifferencePostprocessor
value1 = massflux_right
value2 = massflux_left
[]
[massflux_integral]
type = TimeIntegratedPostprocessor
value = massflux_difference
[]
[mass]
type = ElementIntegralVariablePostprocessor
variable = rhoA
execute_on = 'initial timestep_end'
[]
[mass_change]
type = ChangeOverTimePostprocessor
postprocessor = mass
change_with_respect_to_initial = true
execute_on = 'initial timestep_end'
[]
[mass_conservation]
type = SumPostprocessor
values = 'mass_change massflux_integral'
[]
# MOMENTUM
[momentumflux_left]
type = MomentumFluxIntegral
boundary = inlet
arhouA = rhouA
vel = vel
p = p
A = A
[]
[momentumflux_right]
type = MomentumFluxIntegral
boundary = outlet
arhouA = rhouA
vel = vel
p = p
A = A
[]
[momentumflux_difference]
type = DifferencePostprocessor
value1 = momentumflux_right
value2 = momentumflux_left
[]
[momentumflux_integral]
type = TimeIntegratedPostprocessor
value = momentumflux_difference
[]
[momentum]
type = ElementIntegralVariablePostprocessor
variable = rhouA
execute_on = 'initial timestep_end'
[]
[momentum_change]
type = ChangeOverTimePostprocessor
postprocessor = momentum
change_with_respect_to_initial = true
execute_on = 'initial timestep_end'
[]
[momentum_conservation]
type = SumPostprocessor
values = 'momentum_change momentumflux_integral'
[]
# ENERGY
[energyflux_left]
type = EnergyFluxIntegral
boundary = inlet
arhouA = rhouA
H = H
[]
[energyflux_right]
type = EnergyFluxIntegral
boundary = outlet
arhouA = rhouA
H = H
[]
[energyflux_difference]
type = DifferencePostprocessor
value1 = energyflux_right
value2 = energyflux_left
[]
[energyflux_integral]
type = TimeIntegratedPostprocessor
value = energyflux_difference
[]
[energy]
type = ElementIntegralVariablePostprocessor
variable = rhoEA
execute_on = 'initial timestep_end'
[]
[energy_change]
type = ChangeOverTimePostprocessor
postprocessor = energy
change_with_respect_to_initial = true
execute_on = 'initial timestep_end'
[]
[energy_conservation]
type = SumPostprocessor
values = 'energy_change energyflux_integral'
[]
[]
[Outputs]
[console]
type = Console
show = 'mass_conservation momentum_conservation energy_conservation'
[]
velocity_as_vector = false
[]
(modules/thermal_hydraulics/test/tests/components/junction_one_to_one_1phase/junction_one_to_one_1phase.i)
# This input file simulates the Sod shock tube using a junction in the middle
# of the domain. The solution should be exactly equivalent to the problem with
# no junction. This test examines the solutions at the junction connections
# and compares them to gold values generated from a version of this input file
# that has no junction.
[GlobalParams]
gravity_vector = '0 0 0'
closures = simple_closures
[]
[Functions]
[p_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = '0.5 1.0'
y = '1.0 0.1'
[]
[T_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = '0.5 1.0'
y = '1.4 1.12'
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
gamma = 1.4
molar_mass = 11.64024372
[]
[]
[Closures]
[simple_closures]
type = Closures1PhaseSimple
[]
[]
[Components]
[left_boundary]
type = FreeBoundary1Phase
input = 'left_channel:in'
[]
[left_channel]
type = FlowChannel1Phase
fp = fp
position = '0 0 0'
orientation = '1 0 0'
length = 0.5
n_elems = 50
A = 1.0
initial_T = T_ic_fn
initial_p = p_ic_fn
initial_vel = 0
f = 0
[]
[junction]
type = JunctionOneToOne1Phase
connections = 'left_channel:out right_channel:in'
[]
[right_channel]
type = FlowChannel1Phase
fp = fp
position = '0.5 0 0'
orientation = '1 0 0'
length = 0.5
n_elems = 50
A = 1.0
initial_T = T_ic_fn
initial_p = p_ic_fn
initial_vel = 0
f = 0
[]
[right_boundary]
type = FreeBoundary1Phase
input = 'right_channel:out'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
scheme = bdf2
solve_type = NEWTON
nl_rel_tol = 1e-10
nl_abs_tol = 1e-8
nl_max_its = 60
l_tol = 1e-4
start_time = 0.0
dt = 1e-3
num_steps = 5
abort_on_solve_fail = true
[]
[Postprocessors]
[rhoA_left]
type = SideAverageValue
variable = rhoA
boundary = left_channel:out
execute_on = 'initial timestep_end'
[]
[rhouA_left]
type = SideAverageValue
variable = rhouA
boundary = left_channel:out
execute_on = 'initial timestep_end'
[]
[rhoEA_left]
type = SideAverageValue
variable = rhoEA
boundary = left_channel:out
execute_on = 'initial timestep_end'
[]
[rhoA_right]
type = SideAverageValue
variable = rhoA
boundary = right_channel:in
execute_on = 'initial timestep_end'
[]
# rhouA_right is added by tests file
[rhoEA_right]
type = SideAverageValue
variable = rhoEA
boundary = right_channel:in
execute_on = 'initial timestep_end'
[]
# This is present to test that junction sidesets work properly
[p_avg_junction]
type = SideAverageValue
boundary = 'junction'
variable = p
execute_on = 'initial timestep_end'
[]
[]
[Outputs]
csv = true
show = 'rhoA_left rhouA_left rhoEA_left rhoA_right rhouA_right rhoEA_right'
execute_on = 'initial timestep_end'
[]
(modules/thermal_hydraulics/test/tests/problems/lax_shock_tube/lax_shock_tube.i)
# This test problem is the Lax shock tube test problem,
# which is a Riemann problem with the following parameters:
# * domain = (0,1)
# * gravity = 0
# * EoS: Ideal gas EoS with gamma = 1.4, R = 0.71428571428571428571
# * interface: x = 0.5
# * typical end time: 0.15
# Left initial values:
# * rho = 0.445
# * vel = 0.692
# * p = 3.52874226
# Right initial values:
# * rho = 0.5
# * vel = 0
# * p = 0.571
[GlobalParams]
gravity_vector = '0 0 0'
rdg_slope_reconstruction = minmod
closures = simple_closures
[]
[Functions]
[p_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = '0.5 1.0'
y = '3.52874226 0.571'
[]
[T_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = '0.5 1.0'
y = '11.1016610426966 1.5988'
[]
[vel_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = '0.5 1.0'
y = '0.692 0.0'
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
gamma = 1.4
molar_mass = 11.64024372
[]
[]
[Closures]
[simple_closures]
type = Closures1PhaseSimple
[]
[]
[Components]
[pipe]
type = FlowChannel1Phase
fp = fp
# geometry
position = '0 0 0'
orientation = '1 0 0'
length = 1.0
n_elems = 100
A = 1.0
# IC
initial_T = T_ic_fn
initial_p = p_ic_fn
initial_vel = vel_ic_fn
f = 0
[]
[left_boundary]
type = FreeBoundary1Phase
input = 'pipe:in'
[]
[right_boundary]
type = FreeBoundary1Phase
input = 'pipe:out'
[]
[]
[Executioner]
type = Transient
[TimeIntegrator]
type = ExplicitSSPRungeKutta
# add order via 'cli_args' in 'tests'
[]
solve_type = LINEAR
l_tol = 1e-4
nl_rel_tol = 1e-8
nl_abs_tol = 1e-8
nl_max_its = 60
# run to t = 0.15
start_time = 0.0
dt = 1e-3
num_steps = 150
abort_on_solve_fail = true
[]
[Outputs]
file_base = 'lax_shock_tube'
velocity_as_vector = false
execute_on = 'initial timestep_end'
[out]
type = Exodus
show = 'rho p vel'
[]
[]
(modules/thermal_hydraulics/test/tests/components/junction_one_to_one_1phase/no_junction_1phase.i)
# This input file is used to generate gold values for the junction_one_to_one_1phase.i
# test. Unlike junction_one_to_one_1phase.i, this file has no junction in the
# middle of the domain. In junction_one_to_one_1phase.i, the post-processors are
# side post-processors, but in this input file, side post-processors cannot be
# used to obtain the solution at these positions since there are no sides there.
# Therefore, the solution is sampled at points just to the left and right of
# the middle to obtain the piecewise constant solution values to either side of
# the interface.
[GlobalParams]
gravity_vector = '0 0 0'
closures = simple_closures
[]
[Functions]
[p_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = '0.5 1.0'
y = '1.0 0.1'
[]
[T_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = '0.5 1.0'
y = '1.4 1.12'
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
gamma = 1.4
molar_mass = 11.64024372
[]
[]
[Closures]
[simple_closures]
type = Closures1PhaseSimple
[]
[]
[Components]
[left_boundary]
type = FreeBoundary1Phase
input = 'channel:in'
[]
[channel]
type = FlowChannel1Phase
fp = fp
position = '0 0 0'
orientation = '1 0 0'
length = 1.0
n_elems = 100
A = 1.0
initial_T = T_ic_fn
initial_p = p_ic_fn
initial_vel = 0
f = 0
[]
[right_boundary]
type = FreeBoundary1Phase
input = 'channel:out'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
scheme = bdf2
solve_type = NEWTON
nl_rel_tol = 1e-10
nl_abs_tol = 1e-8
nl_max_its = 60
l_tol = 1e-4
start_time = 0.0
dt = 1e-3
num_steps = 5
abort_on_solve_fail = true
[]
[Postprocessors]
[rhoA_left]
type = PointValue
variable = rhoA
point = '0.4999 0 0'
execute_on = 'initial timestep_end'
[]
[rhouA_left]
type = PointValue
variable = rhouA
point = '0.4999 0 0'
execute_on = 'initial timestep_end'
[]
[rhoEA_left]
type = PointValue
variable = rhoEA
point = '0.4999 0 0'
execute_on = 'initial timestep_end'
[]
[rhoA_right]
type = PointValue
variable = rhoA
point = '0.5001 0 0'
execute_on = 'initial timestep_end'
[]
[rhouA_right]
type = PointValue
variable = rhouA
point = '0.5001 0 0'
execute_on = 'initial timestep_end'
[]
[rhoEA_right]
type = PointValue
variable = rhoEA
point = '0.5001 0 0'
execute_on = 'initial timestep_end'
[]
[]
[Outputs]
csv = true
file_base = 'junction_one_to_one_1phase_out'
execute_on = 'initial timestep_end'
[]
(modules/thermal_hydraulics/test/tests/problems/sod_shock_tube/sod_shock_tube.i)
# This test problem is the classic Sod shock tube test problem,
# which is a Riemann problem with the following parameters:
# * domain = (0,1)
# * gravity = 0
# * EoS: Ideal gas EoS with gamma = 1.4, R = 0.71428571428571428571
# * interface: x = 0.5
# * typical end time: 0.2
# Left initial values:
# * rho = 1
# * vel = 0
# * p = 1
# Right initial values:
# * rho = 0.125
# * vel = 0
# * p = 0.1
#
# The output can be viewed by opening Paraview with the state file `plot.pvsm`:
# paraview --state=plot.pvsm
# This will plot the numerical solution against the analytical solution
[GlobalParams]
gravity_vector = '0 0 0'
rdg_slope_reconstruction = minmod
closures = simple_closures
[]
[Functions]
[p_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = '0.5 1.0'
y = '1.0 0.1'
[]
[T_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = '0.5 1.0'
y = '1.4 1.12'
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
gamma = 1.4
molar_mass = 11.64024372
[]
[]
[Closures]
[simple_closures]
type = Closures1PhaseSimple
[]
[]
[Components]
[pipe]
type = FlowChannel1Phase
fp = fp
# geometry
position = '0 0 0'
orientation = '1 0 0'
length = 1.0
n_elems = 100
A = 1.0
# IC
initial_T = T_ic_fn
initial_p = p_ic_fn
initial_vel = 0
f = 0
[]
[left_boundary]
type = FreeBoundary1Phase
input = 'pipe:in'
[]
[right_boundary]
type = FreeBoundary1Phase
input = 'pipe:out'
[]
[]
[Executioner]
type = Transient
[TimeIntegrator]
type = ExplicitSSPRungeKutta
[]
solve_type = LINEAR
l_tol = 1e-4
nl_rel_tol = 1e-20
nl_abs_tol = 1e-8
nl_max_its = 60
# run to t = 0.2
start_time = 0.0
dt = 1e-3
num_steps = 200
abort_on_solve_fail = true
[]
[Outputs]
file_base = 'sod_shock_tube'
velocity_as_vector = false
execute_on = 'initial timestep_end'
[out]
type = Exodus
show = 'rho p vel'
[]
[]
(modules/thermal_hydraulics/test/tests/problems/square_wave/square_wave.i)
# Square wave problem
[GlobalParams]
gravity_vector = '0 0 0'
rdg_slope_reconstruction = minmod
closures = simple_closures
[]
[Functions]
[T_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = '0.1 0.6 1.0'
y = '2.8 1.4 2.8'
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
gamma = 1.4
molar_mass = 11.64024372
[]
[]
[Closures]
[simple_closures]
type = Closures1PhaseSimple
[]
[]
[Components]
[pipe]
type = FlowChannel1Phase
fp = fp
# geometry
position = '0 0 0'
orientation = '1 0 0'
length = 1.0
n_elems = 400
A = 1.0
# IC
initial_T = T_ic_fn
initial_p = 1
initial_vel = 1
f = 0
[]
[left_boundary]
type = FreeBoundary1Phase
input = 'pipe:in'
[]
[right_boundary]
type = FreeBoundary1Phase
input = 'pipe:out'
[]
[]
[Executioner]
type = Transient
[TimeIntegrator]
type = ExplicitSSPRungeKutta
order = 2
[]
solve_type = LINEAR
l_tol = 1e-4
nl_rel_tol = 1e-20
nl_abs_tol = 1e-8
nl_max_its = 60
# run to t = 0.3
start_time = 0.0
dt = 2e-4
num_steps = 1500
abort_on_solve_fail = true
[]
[Outputs]
file_base = 'square_wave'
velocity_as_vector = false
execute_on = 'initial timestep_end'
[out]
type = Exodus
show = 'p T vel'
[]
[]
(modules/thermal_hydraulics/test/tests/problems/sedov_blast_wave/sedov_blast_wave.i)
# This test problem is the Sedov blast wave test problem,
# which is a Riemann problem with the following parameters:
# * domain = (0,1)
# * gravity = 0
# * EoS: Ideal gas EoS with gamma = 1.4, R = 0.71428571428571428571
# * interface: x = 0.5
# * typical end time: 0.15
# Left initial values:
# * rho = 0.445
# * vel = 0.692
# * p = 3.52874226
# Right initial values:
# * rho = 0.5
# * vel = 0
# * p = 0.571
[GlobalParams]
gravity_vector = '0 0 0'
closures = simple_closures
[]
[Functions]
[p_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = '0.0025 1'
y = '1.591549333333333e+06 6.666666666666668e-09'
[]
[T_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = '0.0025 1'
y = '2.228169066666667e+06 9.333333333333334e-09'
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
gamma = 1.66666666666666666667
molar_mass = 11.64024372
[]
[]
[Closures]
[simple_closures]
type = Closures1PhaseSimple
[]
[]
[Components]
[pipe]
type = FlowChannel1Phase
fp = fp
# geometry
position = '0 0 0'
orientation = '1 0 0'
length = 1.0
n_elems = 400
A = 1.0
# IC
initial_T = T_ic_fn
initial_p = p_ic_fn
initial_vel = 0
f = 0
[]
[left_boundary]
type = SolidWall1Phase
input = 'pipe:in'
[]
[right_boundary]
type = FreeBoundary1Phase
input = 'pipe:out'
[]
[]
[Executioner]
type = Transient
[TimeIntegrator]
type = ExplicitSSPRungeKutta
order = 2
[]
solve_type = LINEAR
l_tol = 1e-4
nl_rel_tol = 1e-20
nl_abs_tol = 1e-8
nl_max_its = 60
# run to t = 0.005
start_time = 0.0
dt = 1e-6
num_steps = 5000
abort_on_solve_fail = true
[]
[Outputs]
file_base = 'sedov_blast_wave'
velocity_as_vector = false
execute_on = 'initial timestep_end'
[out]
type = Exodus
show = 'p T vel'
[]
[]
(modules/thermal_hydraulics/test/tests/problems/super_sonic_tube/test.i)
[GlobalParams]
gravity_vector = '0 0 0'
scaling_factor_1phase = '1 1e-2 1e-4'
initial_p = 101325
initial_T = 300
initial_vel = 522.676
closures = simple_closures
spatial_discretization = cg
[]
[FluidProperties]
[ig]
type = IdealGasFluidProperties
gamma = 1.41
molar_mass = 0.028966206103678928
[]
[]
[Closures]
[simple_closures]
type = Closures1PhaseSimple
[]
[]
[Components]
[pipe]
type = FlowChannel1Phase
fp = ig
# geometry
position = '0 0 0'
orientation = '1 0 0'
A = 1.
D_h = 1.12837916709551
f = 0.0
length = 1
n_elems = 100
[]
[inlet]
type = SupersonicInlet
input = 'pipe:in'
p = 101325
T = 300.0
vel = 522.676
[]
[outlet]
type = FreeBoundary1Phase
input = 'pipe:out'
[]
[]
[Preconditioning]
[SMP_PJFNK]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
scheme = 'bdf2'
start_time = 0
dt = 1e-5
num_steps = 10
abort_on_solve_fail = true
solve_type = 'PJFNK'
nl_rel_tol = 1e-9
nl_abs_tol = 1e-8
nl_max_its = 30
l_tol = 1e-3
l_max_its = 100
[Quadrature]
type = TRAP
order = FIRST
[]
[]
[Outputs]
[out]
type = Exodus
[]
[]
(modules/thermal_hydraulics/test/tests/problems/double_rarefaction/1phase.i)
# Riemann problem that has a double-rarefaction solution
[GlobalParams]
gravity_vector = '0 0 0'
rdg_slope_reconstruction = minmod
closures = simple_closures
[]
[Functions]
[vel_ic_fn]
type = PiecewiseConstant
axis = x
direction = right
x = ' 0.0 0.1'
y = '-1.0 1.0'
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
gamma = 1.4
molar_mass = 11.64024372
[]
[]
[Closures]
[simple_closures]
type = Closures1PhaseSimple
[]
[]
[Components]
[pipe]
type = FlowChannel1Phase
fp = fp
# geometry
position = '-1 0 0'
orientation = '1 0 0'
length = 2.0
n_elems = 100
A = 1.0
# IC
initial_T = 0.04
initial_p = 0.2
initial_vel = vel_ic_fn
f = 0
[]
[left_boundary]
type = FreeBoundary1Phase
input = 'pipe:in'
[]
[right_boundary]
type = FreeBoundary1Phase
input = 'pipe:out'
[]
[]
[Executioner]
type = Transient
[TimeIntegrator]
type = ExplicitSSPRungeKutta
order = 2
[]
solve_type = LINEAR
l_tol = 1e-4
nl_rel_tol = 1e-20
nl_abs_tol = 1e-8
nl_max_its = 60
# run to t = 0.6
start_time = 0.0
dt = 1e-3
num_steps = 600
abort_on_solve_fail = true
[]
[Outputs]
file_base = '1phase'
velocity_as_vector = false
execute_on = 'initial timestep_end'
[out]
type = Exodus
show = 'p T vel'
[]
[]