- 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
- functorThe diffusive flux functor for this boundary condition. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:The diffusive flux functor for this boundary condition. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
- variableThe name of the variable that this boundary condition applies to
C++ Type:LinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable that this boundary condition applies to
LinearFVAdvectionDiffusionFunctorNeumannBC
Description
LinearFVAdvectionDiffusionFunctorNeumannBC specifies the diffusive flux at the boundary. The value will be determined by a Functor (through the "functor" parameter).
Specifically, the specified functor defines the following boundary flux for variable
where the is the boundary normal unit vector and is the diffusion coefficient at the boundary.
This boundary condition should only be used for problems which involve advection and/or diffusion problems.
Input Parameters
- diffusion_coeff1The diffusion coefficient. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
Default:1
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:The diffusion 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_tagsrhsThe tag for the vectors this Kernel should fill
Default:rhs
C++ Type:MultiMooseEnum
Options:rhs, 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.
- 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
- search_methodnearest_node_connected_sidesChoice of search algorithm. All options begin by finding the nearest node in the primary boundary to a query point in the secondary boundary. In the default nearest_node_connected_sides algorithm, primary boundary elements are searched iff that nearest node is one of their nodes. This is fast to determine via a pregenerated node-to-elem map and is robust on conforming meshes. In the optional all_proximate_sides algorithm, primary boundary elements are searched iff they touch that nearest node, even if they are not topologically connected to it. This is more CPU-intensive but is necessary for robustness on any boundary surfaces which has disconnections (such as Flex IGA meshes) or non-conformity (such as hanging nodes in adaptively h-refined meshes).
Default:nearest_node_connected_sides
C++ Type:MooseEnum
Options:nearest_node_connected_sides, all_proximate_sides
Controllable:No
Description:Choice of search algorithm. All options begin by finding the nearest node in the primary boundary to a query point in the secondary boundary. In the default nearest_node_connected_sides algorithm, primary boundary elements are searched iff that nearest node is one of their nodes. This is fast to determine via a pregenerated node-to-elem map and is robust on conforming meshes. In the optional all_proximate_sides algorithm, primary boundary elements are searched iff they touch that nearest node, even if they are not topologically connected to it. This is more CPU-intensive but is necessary for robustness on any boundary surfaces which has disconnections (such as Flex IGA meshes) or non-conformity (such as hanging nodes in adaptively h-refined meshes).
Advanced Parameters
Input Files
- (modules/navier_stokes/test/tests/finite_volume/ins/natural_convection/linear_segregated/2d/diff_heated_cavity_linear_buoyancy.i)
- (test/tests/linearfvkernels/diffusion-reaction-advection/advection-diffusion-reaction-1d.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/cht/conjugate_heat_transfer/cht_neu-dir.i)
- (test/tests/linearfvkernels/diffusion-reaction-advection/advection-diffusion-reaction-2d.i)
- (test/tests/linearfvkernels/diffusion/diffusion-1d_neumann.i)
- (test/tests/linearfvkernels/diffusion/diffusion-2d_neumann.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/cht/bulk_heat_transfer/flow-around-square-linear-fluidonly.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/cht/bulk_heat_transfer/flow-around-square-linear.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/cht/conjugate_heat_transfer/cht_rob-rob.i)
functor
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:The diffusive flux functor for this boundary condition. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number.
(modules/navier_stokes/test/tests/finite_volume/ins/natural_convection/linear_segregated/2d/diff_heated_cavity_linear_buoyancy.i)
################################################################################
# MATERIAL PROPERTIES
################################################################################
rho_0 = 3279.
mu = 1.0
k_cond = 38.0
cp = ${fparse 640}
alpha_b = 3.26e-5
T_0 = 875.0
walls = 'right left top bottom'
[GlobalParams]
rhie_chow_user_object = 'ins_rhie_chow_interpolator'
advected_interp_method = 'upwind'
u = vel_x
v = vel_y
[]
[Problem]
linear_sys_names = 'u_system v_system pressure_system energy_system'
previous_nl_solution_required = true
[]
################################################################################
# GEOMETRY
################################################################################
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 1
xmax = 2.0
ymin = 1.0
ymax = 2.0
nx = 15
ny = 15
[]
[]
################################################################################
# EQUATIONS: VARIABLES, KERNELS & BCS
################################################################################
[UserObjects]
[ins_rhie_chow_interpolator]
type = RhieChowMassFlux
u = vel_x
v = vel_y
pressure = pressure
rho = 'rho'
p_diffusion_kernel = p_diffusion
body_force_kernel_names = 'u_buoyancy; v_buoyancy'
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
solver_sys = u_system
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
[]
[pressure]
type = MooseLinearVariableFVReal
initial_condition = 0
solver_sys = pressure_system
[]
[T_fluid]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = 875
[]
[]
[LinearFVKernels]
[u_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_x
mu = ${mu}
momentum_component = 'x'
use_nonorthogonal_correction = false
[]
[u_pressure]
type = LinearFVMomentumPressure
variable = vel_x
pressure = pressure
momentum_component = 'x'
[]
[u_buoyancy]
type = LinearFVMomentumBuoyancy
variable = vel_x
rho = 'rho'
reference_rho = ${rho_0}
gravity = '0 -9.81 0'
momentum_component = 'x'
[]
[v_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_y
mu = ${mu}
momentum_component = 'y'
use_nonorthogonal_correction = false
[]
[v_pressure]
type = LinearFVMomentumPressure
variable = vel_y
pressure = pressure
momentum_component = 'y'
[]
[v_buoyancy]
type = LinearFVMomentumBuoyancy
variable = vel_y
rho = 'rho'
reference_rho = ${rho_0}
gravity = '0 -9.81 0'
momentum_component = 'y'
[]
[p_diffusion]
type = LinearFVAnisotropicDiffusion
variable = pressure
diffusion_tensor = Ainv
use_nonorthogonal_correction = false
[]
[HbyA_divergence]
type = LinearFVDivergence
variable = pressure
face_flux = HbyA
force_boundary_execution = false
[]
####### FUEL ENERGY EQUATION #######
[heat_advection]
type = LinearFVEnergyAdvection
variable = T_fluid
advected_quantity = temperature
cp = ${cp}
[]
[conduction]
type = LinearFVDiffusion
variable = T_fluid
diffusion_coeff = ${k_cond}
[]
[]
[LinearFVBCs]
[no-slip-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = vel_x
boundary = ${walls}
functor = 0
[]
[no-slip-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = vel_y
boundary = ${walls}
functor = 0
[]
[T_cold]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T_fluid
boundary = 'right'
functor = 870.0
[]
[T_hot]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T_fluid
boundary = 'left'
functor = 880.0
[]
[T_all]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
variable = T_fluid
boundary = 'top bottom'
functor = 0.0
[]
[pressure]
type = LinearFVPressureFluxBC
boundary = 'top bottom left right'
variable = pressure
HbyA_flux = HbyA
Ainv = Ainv
[]
[]
[FunctorMaterials]
[rho_function]
type = ParsedFunctorMaterial
property_name = 'rho'
functor_names = 'T_fluid'
expression = '${rho_0}*(1-${alpha_b}*(T_fluid-${T_0})) '
[]
[]
################################################################################
# EXECUTION / SOLVE
################################################################################
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-11
pressure_l_abs_tol = 1e-11
energy_l_abs_tol = 1e-11
momentum_l_tol = 0
pressure_l_tol = 0
energy_l_tol = 0
rhie_chow_user_object = 'ins_rhie_chow_interpolator'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
energy_system = 'energy_system'
momentum_equation_relaxation = 0.7
pressure_variable_relaxation = 0.3
energy_equation_relaxation = 0.9
num_iterations = 1500
pressure_absolute_tolerance = 1e-8
momentum_absolute_tolerance = 1e-8
energy_absolute_tolerance = 1e-8
print_fields = false
momentum_l_max_its = 300
pin_pressure = true
pressure_pin_value = 0.0
pressure_pin_point = '1.5 1.5 0.0'
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
energy_petsc_options_iname = '-pc_type -pc_hypre_type'
energy_petsc_options_value = 'hypre boomeramg'
continue_on_max_its = true
[]
################################################################################
# SIMULATION OUTPUTS
################################################################################
[Outputs]
exodus = true
[]
(test/tests/linearfvkernels/diffusion-reaction-advection/advection-diffusion-reaction-1d.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 1
nx = 2
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = diff_coeff_func
use_nonorthogonal_correction = false
[]
[advection]
type = LinearFVAdvection
variable = u
velocity = "0.5 0 0"
advected_interp_method = average
[]
[reaction]
type = LinearFVReaction
variable = u
coeff = coeff_func
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
inactive = "outflow neumann"
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right"
functor = analytic_solution
[]
[outflow]
type = LinearFVAdvectionDiffusionOutflowBC
variable = u
boundary = "right"
use_two_term_expansion = true
[]
[neumann]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
variable = u
boundary = "left"
functor = analytic_solution_neumann_left
diffusion_coeff = diff_coeff_func
[]
[]
[Functions]
[diff_coeff_func]
type = ParsedFunction
expression = '1+0.5*x'
[]
[coeff_func]
type = ParsedFunction
expression = '1+1/(1+x)'
[]
[source_func]
type = ParsedFunction
expression = '(1+1/(x+1))*(sin(pi/2*x)+1.5)+0.25*pi*pi*(0.5*x+1)*sin(pi/2*x)'
[]
[analytic_solution]
type = ParsedFunction
expression = 'sin(pi/2*x)+1.5'
[]
[analytic_solution_neumann_left]
type = ParsedFunction
expression = '-(1+0.5*x)*cos(pi/2*x)*pi/2'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[]
[Convergence]
[linear]
type = IterationCountConvergence
max_iterations = 1
converge_at_max_iterations = true
[]
[]
[Executioner]
type = Steady
system_names = u_sys
multi_system_fixed_point=true
multi_system_fixed_point_convergence=linear
petsc_options_iname = '-pc_type -pc_hypre_type -ksp_rtol'
petsc_options_value = 'hypre boomeramg 1e-10'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/cht/conjugate_heat_transfer/cht_neu-dir.i)
### benchmark sources:
### https://doi.org/10.1016/j.compfluid.2018.06.016
### https://doi.org/10.1016/0017-9310(74)90087-8
b = 0.01 # plate thickness
l = 0.2 # plate length
nxi = 24 # nx in the inlet/entrance region
nyf = 18 # ny in fluid
nxf = 24 # nx in the main fluid region
nys = 8 # ny in the solid domain
fx1_bias = 1.00 # bdry layer bias - fluid
fx2_bias = '${fparse 1.0/1.00}' # bdry layer bias - solid
fy_bias = 1.20 # bdry layer bias - fluid
sy_bias = '${fparse 1.0/1.05}' # bdry layer bias - solid
k_s = 0.2876
rho = 0.3525
mu = 3.95e-5
k = 0.06808
cp = 1142.6
vin = 12.0
Tin = 1000.0
T_s_bottom = 600.0
P_out = 1.03e5
h_s = 0.0
advected_interp_method = 'upwind'
[Mesh]
[fluid_channel]
type = GeneratedMeshGenerator
dim = 2
nx = ${nxf}
ny = ${nyf}
xmin = 0
xmax = ${l}
ymin = 0
ymax = '${fparse 10.0*b}'
subdomain_ids = '1'
subdomain_name = 'fluid'
bias_x = '${fx1_bias}'
bias_y = '${fparse fy_bias}'
boundary_name_prefix = 'fluid'
[]
[solid_base]
type = GeneratedMeshGenerator
dim = 2
nx = ${nxf}
ny = ${nys}
xmin = 0
xmax = ${l}
ymin = '${fparse -b}'
ymax = 0
subdomain_ids = '2'
subdomain_name = 'solid'
bias_x = ${fx1_bias}
bias_y = '${fparse sy_bias}'
boundary_id_offset = 10
boundary_name_prefix = 'solid'
[]
[entrance]
type = GeneratedMeshGenerator
dim = 2
nx = '${fparse 2.0*nxi}'
ny = ${nyf}
xmin = '${fparse -2.0*l}'
xmax = 0
ymin = 0
ymax = '${fparse 10.0*b}'
subdomain_ids = '0'
subdomain_name = 'entrance'
bias_x = ${fx2_bias}
bias_y = '${fparse fy_bias}'
boundary_id_offset = 20
boundary_name_prefix = 'ent'
[]
[smg]
type = StitchedMeshGenerator
inputs = 'entrance fluid_channel solid_base'
stitch_boundaries_pairs = 'ent_right fluid_left;
fluid_bottom solid_top'
prevent_boundary_ids_overlap = false
[]
[interface]
type = SideSetsBetweenSubdomainsGenerator
input = 'smg'
primary_block = 'fluid'
paired_block = 'solid'
new_boundary = interface
[]
[symmetry_transform]
type = SymmetryTransformGenerator
input = interface
mirror_point = '0 0 0'
mirror_normal_vector = '0 1 0'
[]
inactive = 'symmetry_transform'
[]
[Problem]
linear_sys_names = 'u_system v_system pressure_system energy_system solid_energy_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = RhieChowMassFlux
u = vel_x
v = vel_y
pressure = pressure
rho = ${rho}
p_diffusion_kernel = p_diffusion
block = '0 1'
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
initial_condition = ${vin}
solver_sys = u_system
block = '0 1'
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
initial_condition = 0.0
block = '0 1'
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = ${P_out}
block = '0 1'
[]
[T_fluid]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = ${Tin}
block = '0 1'
[]
[T_solid]
type = MooseLinearVariableFVReal
solver_sys = solid_energy_system
initial_condition = ${T_s_bottom}
block = 2
[]
[]
[LinearFVKernels]
[u_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_x
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'x'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = false
[]
[v_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_y
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'y'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = false
[]
[u_pressure]
type = LinearFVMomentumPressure
variable = vel_x
pressure = pressure
momentum_component = 'x'
[]
[v_pressure]
type = LinearFVMomentumPressure
variable = vel_y
pressure = pressure
momentum_component = 'y'
[]
[p_diffusion]
type = LinearFVAnisotropicDiffusion
variable = pressure
diffusion_tensor = Ainv
use_nonorthogonal_correction = false
[]
[HbyA_divergence]
type = LinearFVDivergence
variable = pressure
face_flux = HbyA
force_boundary_execution = true
[]
[h_advection]
type = LinearFVEnergyAdvection
variable = T_fluid
advected_quantity = temperature
cp = ${cp}
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
[]
[conduction]
type = LinearFVDiffusion
variable = T_fluid
diffusion_coeff = ${k}
use_nonorthogonal_correction = false
[]
[solid-conduction]
type = LinearFVDiffusion
variable = T_solid
diffusion_coeff = ${k_s}
use_nonorthogonal_correction = false
[]
[]
[LinearFVBCs]
# velocity BCs
[inlet-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'ent_left'
variable = vel_x
functor = ${vin}
[]
[inlet-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'ent_left'
variable = vel_y
functor = '0.000'
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'ent_bottom interface'
variable = vel_x
functor = 0.0
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'ent_bottom interface'
variable = vel_y
functor = 0.0
[]
[outlet_p]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'fluid_right'
variable = pressure
functor = ${P_out}
[]
[outlet_u]
type = LinearFVAdvectionDiffusionOutflowBC
boundary = 'fluid_right'
variable = vel_x
use_two_term_expansion = false
[]
[outlet_v]
type = LinearFVAdvectionDiffusionOutflowBC
boundary = 'fluid_right'
variable = vel_y
use_two_term_expansion = false
[]
# freestream BCs for top of fluid domain
[freestream_u]
type = LinearFVAdvectionDiffusionOutflowBC
boundary = 'fluid_top ent_top'
variable = vel_x
use_two_term_expansion = false
[]
[freestream_v]
type = LinearFVAdvectionDiffusionOutflowBC
boundary = 'fluid_top ent_top'
variable = vel_y
use_two_term_expansion = false
[]
[freestream_p]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
boundary = 'fluid_top ent_top'
variable = pressure
functor = 0
[]
# temperature BCs
[inlet_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T_fluid
functor = ${Tin}
boundary = 'ent_left'
[]
[heated_wall_solid]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T_solid
functor = ${T_s_bottom}
boundary = 'solid_bottom'
[]
[insulated_fluid]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
variable = T_fluid
functor = 0
boundary = 'ent_bottom ent_top fluid_top'
[]
[insulated_solid]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
variable = T_solid
functor = 0
boundary = 'solid_left solid_right'
[]
[outlet_T]
type = LinearFVAdvectionDiffusionOutflowBC
variable = T_fluid
use_two_term_expansion = false
boundary = 'fluid_right'
[]
[fluid_solid]
type = LinearFVDirichletCHTBC
variable = T_fluid
boundary = interface
functor = interface_temperature_solid_interface
[]
[solid_fluid]
type = LinearFVRobinCHTBC
variable = T_solid
boundary = interface
h = ${h_s}
thermal_conductivity = ${k_s}
incoming_flux = heat_flux_to_solid_interface
surface_temperature = interface_temperature_fluid_interface
[]
[]
[FunctorMaterials]
[rhocpT]
property_name = 'rhocpT'
type = ParsedFunctorMaterial
functor_names = 'T_fluid'
expression = '${rho}*${cp}*T_fluid'
[]
[]
[Postprocessors]
[h_in]
type = VolumetricFlowRate
boundary = 'ent_left'
vel_x = vel_x
vel_y = vel_y
rhie_chow_user_object = rc
advected_quantity = 'rhocpT'
subtract_mesh_velocity = false
[]
[h_out]
type = VolumetricFlowRate
boundary = 'fluid_right fluid_top ent_top interface'
vel_x = vel_x
vel_y = vel_y
rhie_chow_user_object = rc
advected_quantity = 'rhocpT'
advected_interp_method = upwind
subtract_mesh_velocity = false
[]
[]
[VectorPostprocessors]
[y_vs_ts]
type = LineValueSampler
variable = 'T_solid'
start_point = '0.05 -1e-9 0' # making sure we are always in the domain
end_point = '0.05 ${fparse -b+1e-9} 0'
num_points = 8
sort_by = id
warn_discontinuous_face_values = false
[]
[y_vs_tf]
type = LineValueSampler
variable = 'T_fluid'
start_point = '0.05 1e-9 0' # making sure we are always in the domain
end_point = '0.05 ${fparse b-1e-9} 0'
num_points = 12
sort_by = id
warn_discontinuous_face_values = false
[]
[]
[Executioner]
type = SIMPLE
num_iterations = 1000
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
rhie_chow_user_object = 'rc'
momentum_l_abs_tol = 1e-10
pressure_l_abs_tol = 1e-10
momentum_l_tol = 0
pressure_l_tol = 0
momentum_equation_relaxation = 0.9
pressure_variable_relaxation = 0.3
momentum_absolute_tolerance = 1e-7
pressure_absolute_tolerance = 1e-7
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
energy_system = 'energy_system'
solid_energy_system = 'solid_energy_system'
energy_l_abs_tol = 1e-10
solid_energy_l_abs_tol = 1e-10
energy_l_tol = 0
solid_energy_l_tol = 0
energy_equation_relaxation = 1.0
energy_absolute_tolerance = 1e-7
solid_energy_absolute_tolerance = 1e-7
energy_petsc_options_iname = '-pc_type -pc_hypre_type'
energy_petsc_options_value = 'hypre boomeramg'
solid_energy_petsc_options_iname = '-pc_type -pc_hypre_type'
solid_energy_petsc_options_value = 'hypre boomeramg'
cht_interfaces = 'interface'
cht_solid_flux_relaxation = 0.4
cht_fluid_flux_relaxation = 0.4
cht_solid_temperature_relaxation = 0.4
cht_fluid_temperature_relaxation = 0.4
max_cht_fpi = 3
print_fields = false
[]
[Outputs]
csv = true
execute_on = timestep_end
[]
(test/tests/linearfvkernels/diffusion-reaction-advection/advection-diffusion-reaction-2d.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 2
ny = 1
ymax = 0.5
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = diff_coeff_func
use_nonorthogonal_correction = false
[]
[advection]
type = LinearFVAdvection
variable = u
velocity = "0.5 0 0"
advected_interp_method = average
[]
[reaction]
type = LinearFVReaction
variable = u
coeff = coeff_func
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
inactive = "outflow neumann"
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right top bottom"
functor = analytic_solution
[]
[outflow]
type = LinearFVAdvectionDiffusionOutflowBC
variable = u
boundary = "right"
use_two_term_expansion = true
[]
[neumann]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
variable = u
boundary = "top"
functor = analytic_solution_neumann_top
diffusion_coeff = diff_coeff_func
[]
[]
[Functions]
[diff_coeff_func]
type = ParsedFunction
expression = '1.0+0.5*x*y'
[]
[coeff_func]
type = ParsedFunction
expression = '1.0+1.0/(1+x*y)'
[]
[source_func]
type = ParsedFunction
expression = '-1.0*x*pi*sin((1/2)*x*pi)*cos(2*y*pi) - 0.25*y*pi*sin(2*y*pi)*cos((1/2)*x*pi) + (1.0 + 1.0/(x*y + 1))*(sin((1/2)*x*pi)*sin(2*y*pi) + 1.5) + (17/4)*pi^2*(0.5*x*y + 1.0)*sin((1/2)*x*pi)*sin(2*y*pi) + 0.25*pi*sin(2*y*pi)*cos((1/2)*x*pi)'
[]
[analytic_solution]
type = ParsedFunction
expression = 'sin((1/2)*x*pi)*sin(2*y*pi) + 1.5'
[]
[analytic_solution_neumann_top]
type = ParsedFunction
expression = '(1.0+0.5*x*y)*sin((1/2)*x*pi)*cos(2*y*pi)*2*pi'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[]
[Convergence]
[linear]
type = IterationCountConvergence
max_iterations = 1
converge_at_max_iterations = true
[]
[]
[Executioner]
type = Steady
system_names = u_sys
l_tol = 1e-10
multi_system_fixed_point=true
multi_system_fixed_point_convergence=linear
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(test/tests/linearfvkernels/diffusion/diffusion-1d_neumann.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 1
nx = 2
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = coeff_func
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "right"
functor = analytic_solution
[]
[neu]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
variable = u
boundary = "left"
functor = analytic_solution_neumann
[]
[]
[Functions]
[coeff_func]
type = ParsedFunction
expression = '0.5*x'
[]
[source_func]
type = ParsedFunction
expression = '2*x'
[]
[analytic_solution]
type = ParsedFunction
expression = '1-x*x'
[]
[analytic_solution_neumann]
type = ParsedFunction
expression = '-(0.5*x)*(-2*x)'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[]
[Executioner]
type = Steady
system_names = u_sys
l_tol = 1e-10
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(test/tests/linearfvkernels/diffusion/diffusion-2d_neumann.i)
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 2
nx = 2
ny = 1
ymax = 0.5
[]
[]
[Problem]
linear_sys_names = 'u_sys'
[]
[Variables]
[u]
type = MooseLinearVariableFVReal
solver_sys = 'u_sys'
initial_condition = 1.0
[]
[]
[LinearFVKernels]
[diffusion]
type = LinearFVDiffusion
variable = u
diffusion_coeff = coeff_func
use_nonorthogonal_correction = false
[]
[source]
type = LinearFVSource
variable = u
source_density = source_func
[]
[]
[LinearFVBCs]
[dir]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = u
boundary = "left right"
functor = analytic_solution
[]
[neu_bottom]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
variable = u
boundary = "bottom"
functor = analytic_solution_neumann_bottom
diffusion_coeff = coeff_func
[]
[neu_top]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
variable = u
boundary = "top"
functor = analytic_solution_neumann_top
diffusion_coeff = coeff_func
[]
[]
[Functions]
[coeff_func]
type = ParsedFunction
expression = '1+0.5*x*y'
[]
[source_func]
type = ParsedFunction
expression = '2*(1.5-y*y)+2*x*y*(1.5-y*y)+2*(1.5-x*x)+2*x*y*(1.5-x*x)'
[]
[analytic_solution]
type = ParsedFunction
expression = '(1.5-x*x)*(1.5-y*y)'
[]
[analytic_solution_neumann_bottom]
type = ParsedFunction
expression = '-(1+0.5*x*y)*(1.5-x*x)*(-2*y)'
[]
[analytic_solution_neumann_top]
type = ParsedFunction
expression = '(1+0.5*x*y)*(1.5-x*x)*(-2*y)'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = FINAL
[]
[error]
type = ElementL2FunctorError
approximate = u
exact = analytic_solution
execute_on = FINAL
[]
[]
[Convergence]
[linear]
type = IterationCountConvergence
max_iterations = 1
converge_at_max_iterations = true
[]
[]
[Executioner]
type = Steady
system_names = u_sys
l_tol = 1e-10
multi_system_fixed_point=true
multi_system_fixed_point_convergence=linear
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
[csv]
type = CSV
execute_on = FINAL
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/cht/bulk_heat_transfer/flow-around-square-linear-fluidonly.i)
mu = 0.01
rho = 1.1
k = 0.0005
cp = 10
h_conv = 5
advected_interp_method = 'upwind'
[Mesh]
[generated_mesh]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 10
xmin = 0
ymin = 0
ymax = 0.1
xmax = 0.1
[]
[subdomain1]
type = SubdomainBoundingBoxGenerator
input = generated_mesh
block_name = subdomain1
bottom_left = '0.04 0.04 0'
block_id = 1
top_right = '0.06 0.06 0'
[]
[interface]
type = SideSetsBetweenSubdomainsGenerator
input = subdomain1
primary_block = 0
paired_block = 1
new_boundary = interface
[]
[delete]
type = BlockDeletionGenerator
input = interface
block = 1
[]
[]
[Problem]
linear_sys_names = 'u_system v_system pressure_system energy_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = RhieChowMassFlux
u = vel_x
v = vel_y
pressure = pressure
rho = ${rho}
p_diffusion_kernel = p_diffusion
block = 0
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
initial_condition = 0.1
solver_sys = u_system
block = 0
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
initial_condition = 0.01
block = 0
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = 0.2
block = 0
[]
[T_fluid]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = 300
block = 0
[]
[]
[LinearFVKernels]
[u_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_x
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'x'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = true
[]
[v_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_y
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'y'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = true
[]
[u_pressure]
type = LinearFVMomentumPressure
variable = vel_x
pressure = pressure
momentum_component = 'x'
[]
[v_pressure]
type = LinearFVMomentumPressure
variable = vel_y
pressure = pressure
momentum_component = 'y'
[]
[p_diffusion]
type = LinearFVAnisotropicDiffusion
variable = pressure
diffusion_tensor = Ainv
use_nonorthogonal_correction = true
[]
[HbyA_divergence]
type = LinearFVDivergence
variable = pressure
face_flux = HbyA
force_boundary_execution = true
[]
[h_advection]
type = LinearFVEnergyAdvection
variable = T_fluid
advected_quantity = temperature
cp = ${cp}
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
[]
[conduction]
type = LinearFVDiffusion
variable = T_fluid
diffusion_coeff = ${k}
use_nonorthogonal_correction = true
[]
[]
[LinearFVBCs]
[inlet-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_x
functor = '0.1'
[]
[inlet-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_y
functor = '0.0'
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom interface'
variable = vel_x
functor = 0.0
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom interface'
variable = vel_y
functor = 0.0
[]
[outlet_p]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'right'
variable = pressure
functor = 1.4
[]
[outlet_u]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_x
use_two_term_expansion = false
boundary = right
[]
[outlet_v]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_y
use_two_term_expansion = false
boundary = right
[]
[inlet_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T_fluid
functor = 300
boundary = 'left'
[]
[walls_T]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
variable = T_fluid
functor = 0.0
boundary = 'top bottom'
[]
[outlet_T]
type = LinearFVAdvectionDiffusionOutflowBC
variable = T_fluid
use_two_term_expansion = false
boundary = right
[]
[fluid_solid]
type = LinearFVConvectiveHeatTransferBC
variable = T_fluid
T_solid = boundary_value
T_fluid = T_fluid
boundary = interface
h = ${h_conv}
[]
[]
[FunctorMaterials]
[rhocpT]
property_name = 'rhocpT'
type = ParsedFunctorMaterial
functor_names = 'T_fluid'
expression = '${rho}*${cp}*T_fluid'
[]
[]
[Functions]
[boundary_value]
type = ConstantFunction
value = 350
[]
[]
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-13
pressure_l_abs_tol = 1e-13
energy_l_abs_tol = 1e-13
momentum_l_tol = 0
pressure_l_tol = 0
energy_l_tol = 0
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
energy_system = 'energy_system'
momentum_equation_relaxation = 0.8
energy_equation_relaxation = 1.0
pressure_variable_relaxation = 0.3
num_iterations = 1000
pressure_absolute_tolerance = 1e-10
momentum_absolute_tolerance = 1e-10
energy_absolute_tolerance = 1e-10
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
energy_petsc_options_iname = '-pc_type -pc_hypre_type'
energy_petsc_options_value = 'hypre boomeramg'
print_fields = false
continue_on_max_its = true
[]
[Outputs]
exodus = true
execute_on = timestep_end
[]
(modules/navier_stokes/test/tests/finite_volume/ins/cht/bulk_heat_transfer/flow-around-square-linear.i)
mu = 0.01
rho = 1.1
k = 0.0005
cp = 10
k_s = 3.0
h_conv = 5
power_density = 10000
advected_interp_method = 'upwind'
[Mesh]
[generated_mesh]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 10
xmin = 0
ymin = 0
ymax = 0.1
xmax = 0.1
[]
[subdomain1]
type = SubdomainBoundingBoxGenerator
input = generated_mesh
block_name = subdomain1
bottom_left = '0.04 0.04 0'
block_id = 1
top_right = '0.06 0.06 0'
[]
[interface]
type = SideSetsBetweenSubdomainsGenerator
input = subdomain1
primary_block = 0
paired_block = 1
new_boundary = interface
[]
[]
[Problem]
linear_sys_names = 'u_system v_system pressure_system energy_system solid_energy_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = RhieChowMassFlux
u = vel_x
v = vel_y
pressure = pressure
rho = ${rho}
p_diffusion_kernel = p_diffusion
block = 0
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
initial_condition = 0.1
solver_sys = u_system
block = 0
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
initial_condition = 0.01
block = 0
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = 0.2
block = 0
[]
[T_fluid]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = 300
block = 0
[]
[T_solid]
type = MooseLinearVariableFVReal
solver_sys = solid_energy_system
initial_condition = 500
block = 1
[]
[]
[LinearFVKernels]
[u_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_x
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'x'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = true
[]
[v_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_y
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'y'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = true
[]
[u_pressure]
type = LinearFVMomentumPressure
variable = vel_x
pressure = pressure
momentum_component = 'x'
[]
[v_pressure]
type = LinearFVMomentumPressure
variable = vel_y
pressure = pressure
momentum_component = 'y'
[]
[p_diffusion]
type = LinearFVAnisotropicDiffusion
variable = pressure
diffusion_tensor = Ainv
use_nonorthogonal_correction = true
[]
[HbyA_divergence]
type = LinearFVDivergence
variable = pressure
face_flux = HbyA
force_boundary_execution = true
[]
[h_advection]
type = LinearFVEnergyAdvection
variable = T_fluid
advected_quantity = temperature
cp = ${cp}
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
[]
[conduction]
type = LinearFVDiffusion
variable = T_fluid
diffusion_coeff = ${k}
use_nonorthogonal_correction = true
[]
[solid-conduction]
type = LinearFVDiffusion
variable = T_solid
diffusion_coeff = ${k_s}
use_nonorthogonal_correction = true
[]
[solid-source]
type = LinearFVSource
variable = T_solid
source_density = ${power_density}
[]
[]
[LinearFVBCs]
[inlet-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_x
functor = '0.1'
[]
[inlet-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_y
functor = '0.0'
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom interface'
variable = vel_x
functor = 0.0
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom interface'
variable = vel_y
functor = 0.0
[]
[outlet_p]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'right'
variable = pressure
functor = 1.4
[]
[outlet_u]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_x
use_two_term_expansion = false
boundary = right
[]
[outlet_v]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_y
use_two_term_expansion = false
boundary = right
[]
[inlet_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T_fluid
functor = 300
boundary = 'left'
[]
[walls_T]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
variable = T_fluid
functor = 0.0
boundary = 'top bottom'
[]
[outlet_T]
type = LinearFVAdvectionDiffusionOutflowBC
variable = T_fluid
use_two_term_expansion = false
boundary = right
[]
[fluid_solid]
type = LinearFVConvectiveHeatTransferBC
variable = T_fluid
T_solid = T_solid
T_fluid = T_fluid
boundary = interface
h = ${h_conv}
[]
[solid_fluid]
type = LinearFVConvectiveHeatTransferBC
variable = T_solid
T_solid = T_solid
T_fluid = T_fluid
boundary = interface
h = ${h_conv}
[]
[]
[FunctorMaterials]
[rhocpT]
property_name = 'rhocpT'
type = ParsedFunctorMaterial
functor_names = 'T_fluid'
expression = '${rho}*${cp}*T_fluid'
[]
[]
[Postprocessors]
[h_in]
type = VolumetricFlowRate
boundary = left
vel_x = vel_x
vel_y = vel_y
rhie_chow_user_object = rc
advected_quantity = 'rhocpT'
subtract_mesh_velocity = false
[]
[h_out]
type = VolumetricFlowRate
boundary = right
vel_x = vel_x
vel_y = vel_y
rhie_chow_user_object = rc
advected_quantity = 'rhocpT'
advected_interp_method = upwind
subtract_mesh_velocity = false
[]
[power]
type = ElementIntegralFunctorPostprocessor
functor = ${power_density}
block = 1
[]
[]
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-13
pressure_l_abs_tol = 1e-13
energy_l_abs_tol = 1e-13
solid_energy_l_abs_tol = 1e-13
momentum_l_tol = 0
pressure_l_tol = 0
energy_l_tol = 0
solid_energy_l_tol = 0
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
energy_system = 'energy_system'
solid_energy_system = 'solid_energy_system'
momentum_equation_relaxation = 0.8
energy_equation_relaxation = 1.0
pressure_variable_relaxation = 0.3
num_iterations = 1000
pressure_absolute_tolerance = 1e-10
momentum_absolute_tolerance = 1e-10
energy_absolute_tolerance = 1e-10
solid_energy_absolute_tolerance = 1e-10
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
energy_petsc_options_iname = '-pc_type -pc_hypre_type'
energy_petsc_options_value = 'hypre boomeramg'
solid_energy_petsc_options_iname = '-pc_type -pc_hypre_type'
solid_energy_petsc_options_value = 'hypre boomeramg'
print_fields = false
continue_on_max_its = true
[]
[Outputs]
exodus = true
execute_on = timestep_end
[]
(modules/navier_stokes/test/tests/finite_volume/ins/cht/conjugate_heat_transfer/cht_rob-rob.i)
### benchmark sources:
### https://doi.org/10.1016/j.compfluid.2018.06.016
### https://doi.org/10.1016/0017-9310(74)90087-8
b = 0.01 # plate thickness
l = 0.2 # plate length
nxi = 24 # nx in the inlet/entrance region
nyf = 18 # ny in fluid
nxf = 24 # nx in the main fluid region
nys = 8 # ny in the solid domain
fx1_bias = 1.00 # bdry layer bias - fluid
fx2_bias = '${fparse 1.0/1.00}' # bdry layer bias - solid
fy_bias = 1.20 # bdry layer bias - fluid
sy_bias = '${fparse 1.0/1.05}' # bdry layer bias - solid
k_s = 0.2876
rho = 0.3525
mu = 3.95e-5
k = 0.06808
cp = 1142.6
vin = 12.0
Tin = 1000.0
T_s_bottom = 600.0
P_out = 1.03e5
h_s = 1.0
h_f = 1.0
advected_interp_method = 'upwind'
[Mesh]
[fluid_channel]
type = GeneratedMeshGenerator
dim = 2
nx = ${nxf}
ny = ${nyf}
xmin = 0
xmax = ${l}
ymin = 0
ymax = '${fparse 10.0*b}'
subdomain_ids = '1'
subdomain_name = 'fluid'
bias_x = '${fx1_bias}'
bias_y = '${fparse fy_bias}'
boundary_name_prefix = 'fluid'
[]
[solid_base]
type = GeneratedMeshGenerator
dim = 2
nx = ${nxf}
ny = ${nys}
xmin = 0
xmax = ${l}
ymin = '${fparse -b}'
ymax = 0
subdomain_ids = '2'
subdomain_name = 'solid'
bias_x = ${fx1_bias}
bias_y = '${fparse sy_bias}'
boundary_id_offset = 10
boundary_name_prefix = 'solid'
[]
[entrance]
type = GeneratedMeshGenerator
dim = 2
nx = '${fparse 2.0*nxi}'
ny = ${nyf}
xmin = '${fparse -2.0*l}'
xmax = 0
ymin = 0
ymax = '${fparse 10.0*b}'
subdomain_ids = '0'
subdomain_name = 'entrance'
bias_x = ${fx2_bias}
bias_y = '${fparse fy_bias}'
boundary_id_offset = 20
boundary_name_prefix = 'ent'
[]
[smg]
type = StitchedMeshGenerator
inputs = 'entrance fluid_channel solid_base'
stitch_boundaries_pairs = 'ent_right fluid_left;
fluid_bottom solid_top'
prevent_boundary_ids_overlap = false
[]
[interface]
type = SideSetsBetweenSubdomainsGenerator
input = 'smg'
primary_block = 'fluid'
paired_block = 'solid'
new_boundary = interface
[]
[symmetry_transform]
type = SymmetryTransformGenerator
input = interface
mirror_point = '0 0 0'
mirror_normal_vector = '0 1 0'
[]
inactive = 'symmetry_transform'
[]
[Problem]
linear_sys_names = 'u_system v_system pressure_system energy_system solid_energy_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = RhieChowMassFlux
u = vel_x
v = vel_y
pressure = pressure
rho = ${rho}
p_diffusion_kernel = p_diffusion
block = '0 1'
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
initial_condition = ${vin}
solver_sys = u_system
block = '0 1'
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
initial_condition = 0.0
block = '0 1'
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = ${P_out}
block = '0 1'
[]
[T_fluid]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = ${Tin}
block = '0 1'
[]
[T_solid]
type = MooseLinearVariableFVReal
solver_sys = solid_energy_system
initial_condition = ${T_s_bottom}
block = 2
[]
[]
[LinearFVKernels]
[u_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_x
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'x'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = false
[]
[v_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_y
advected_interp_method = ${advected_interp_method}
mu = ${mu}
u = vel_x
v = vel_y
momentum_component = 'y'
rhie_chow_user_object = 'rc'
use_nonorthogonal_correction = false
[]
[u_pressure]
type = LinearFVMomentumPressure
variable = vel_x
pressure = pressure
momentum_component = 'x'
[]
[v_pressure]
type = LinearFVMomentumPressure
variable = vel_y
pressure = pressure
momentum_component = 'y'
[]
[p_diffusion]
type = LinearFVAnisotropicDiffusion
variable = pressure
diffusion_tensor = Ainv
use_nonorthogonal_correction = false
[]
[HbyA_divergence]
type = LinearFVDivergence
variable = pressure
face_flux = HbyA
force_boundary_execution = true
[]
[h_advection]
type = LinearFVEnergyAdvection
variable = T_fluid
advected_quantity = temperature
cp = ${cp}
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
[]
[conduction]
type = LinearFVDiffusion
variable = T_fluid
diffusion_coeff = ${k}
use_nonorthogonal_correction = false
[]
[solid-conduction]
type = LinearFVDiffusion
variable = T_solid
diffusion_coeff = ${k_s}
use_nonorthogonal_correction = false
[]
[]
[LinearFVBCs]
# velocity BCs
[inlet-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'ent_left'
variable = vel_x
functor = ${vin}
[]
[inlet-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'ent_left'
variable = vel_y
functor = '0.000'
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'ent_bottom interface'
variable = vel_x
functor = 0.0
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'ent_bottom interface'
variable = vel_y
functor = 0.0
[]
[outlet_p]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'fluid_right'
variable = pressure
functor = ${P_out}
[]
[outlet_u]
type = LinearFVAdvectionDiffusionOutflowBC
boundary = 'fluid_right'
variable = vel_x
use_two_term_expansion = false
[]
[outlet_v]
type = LinearFVAdvectionDiffusionOutflowBC
boundary = 'fluid_right'
variable = vel_y
use_two_term_expansion = false
[]
# freestream BCs for top of fluid domain
[freestream_u]
type = LinearFVAdvectionDiffusionOutflowBC
boundary = 'fluid_top ent_top'
variable = vel_x
use_two_term_expansion = false
[]
[freestream_v]
type = LinearFVAdvectionDiffusionOutflowBC
boundary = 'fluid_top ent_top'
variable = vel_y
use_two_term_expansion = false
[]
[freestream_p]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
boundary = 'fluid_top ent_top'
variable = pressure
functor = 0
[]
# temperature BCs
[inlet_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T_fluid
functor = ${Tin}
boundary = 'ent_left'
[]
[heated_wall_solid]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T_solid
functor = ${T_s_bottom}
boundary = 'solid_bottom'
[]
[insulated_fluid]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
variable = T_fluid
functor = 0
boundary = 'ent_top ent_bottom fluid_top'
[]
[insulated_solid]
type = LinearFVAdvectionDiffusionFunctorNeumannBC
variable = T_solid
functor = 0
boundary = 'solid_left solid_right'
[]
[outlet_T]
type = LinearFVAdvectionDiffusionOutflowBC
variable = T_fluid
use_two_term_expansion = false
boundary = 'fluid_right'
[]
[fluid_solid]
type = LinearFVRobinCHTBC
variable = T_fluid
boundary = interface
h = ${h_f}
incoming_flux = heat_flux_to_fluid_interface
surface_temperature = interface_temperature_solid_interface
thermal_conductivity = ${k}
[]
[solid_fluid]
type = LinearFVRobinCHTBC
variable = T_solid
boundary = interface
h = ${h_s}
incoming_flux = heat_flux_to_solid_interface
surface_temperature = interface_temperature_fluid_interface
thermal_conductivity = ${k_s}
[]
[]
[FunctorMaterials]
[rhocpT]
property_name = 'rhocpT'
type = ParsedFunctorMaterial
functor_names = 'T_fluid'
expression = '${rho}*${cp}*T_fluid'
[]
[]
[Postprocessors]
[h_in]
type = VolumetricFlowRate
boundary = 'ent_left'
vel_x = vel_x
vel_y = vel_y
rhie_chow_user_object = rc
advected_quantity = 'rhocpT'
subtract_mesh_velocity = false
[]
[h_out]
type = VolumetricFlowRate
boundary = 'fluid_right fluid_top ent_top interface'
vel_x = vel_x
vel_y = vel_y
rhie_chow_user_object = rc
advected_quantity = 'rhocpT'
advected_interp_method = upwind
subtract_mesh_velocity = false
[]
[]
[VectorPostprocessors]
[y_vs_ts]
type = LineValueSampler
variable = 'T_solid'
start_point = '0.05 -1e-9 0' # making sure we are always in the domain
end_point = '0.05 ${fparse -b+1e-9} 0'
num_points = 8
sort_by = id
warn_discontinuous_face_values = false
[]
[y_vs_tf]
type = LineValueSampler
variable = 'T_fluid'
start_point = '0.05 1e-9 0' # making sure we are always in the domain
end_point = '0.05 ${fparse b-1e-9} 0'
num_points = 12
sort_by = id
warn_discontinuous_face_values = false
[]
[]
[Executioner]
type = SIMPLE
num_iterations = 1000
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
rhie_chow_user_object = 'rc'
momentum_l_abs_tol = 1e-10
pressure_l_abs_tol = 1e-10
momentum_l_tol = 0
pressure_l_tol = 0
momentum_equation_relaxation = 0.9
pressure_variable_relaxation = 0.3
momentum_absolute_tolerance = 1e-7
pressure_absolute_tolerance = 1e-7
momentum_petsc_options_iname = '-pc_type -pc_hypre_type'
momentum_petsc_options_value = 'hypre boomeramg'
pressure_petsc_options_iname = '-pc_type -pc_hypre_type'
pressure_petsc_options_value = 'hypre boomeramg'
energy_system = 'energy_system'
solid_energy_system = 'solid_energy_system'
energy_l_abs_tol = 1e-10
solid_energy_l_abs_tol = 1e-10
energy_l_tol = 0
solid_energy_l_tol = 0
energy_equation_relaxation = 1.0
energy_absolute_tolerance = 1e-7
solid_energy_absolute_tolerance = 1e-7
energy_petsc_options_iname = '-pc_type -pc_hypre_type'
energy_petsc_options_value = 'hypre boomeramg'
solid_energy_petsc_options_iname = '-pc_type -pc_hypre_type'
solid_energy_petsc_options_value = 'hypre boomeramg'
cht_interfaces = 'interface'
cht_solid_flux_relaxation = 1.0
cht_fluid_flux_relaxation = 1.0
cht_solid_temperature_relaxation = 1.0
cht_fluid_temperature_relaxation = 1.0
max_cht_fpi = 2
print_fields = false
[]
[Outputs]
exodus = true
csv = true
execute_on = timestep_end
[]