- rhie_chow_user_objectThe rhie-chow user-object which is used to determine the face velocity.
C++ Type:UserObjectName
Controllable:No
Description:The rhie-chow user-object which is used to determine the face velocity.
- variableThe name of the variable whose linear system this object contributes to
C++ Type:LinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable whose linear system this object contributes to
LinearFVEnergyAdvection
This kernel adds the contributions of the energy advection term to the matrix and right hand side of the energy equation system for the finite volume SIMPLE segregated solver SIMPLE.
This kernel currently supports the advection of specific enthalpy h or temperature T. Important consideration: Temperature advection is only supported for constant specific heat, where h can be defined as h=cpT. For variable cp, the user should use the enthalpy formulation. Parameter "advected_quantity" lets the user select "enthalpy" or "temperature".
This term is described by ∇⋅(ρuh) for enthalpy or ∇⋅(ρucpT) for constant specific heat. This term is present in the energy equation conservation for an incompressible/weakly-compressible formulation.
For FV, the integral of the advection term over a cell can be expressed as:
VC∫∇⋅(ρuh)dV≈f∑(ρu⋅n)RChf∣Sf∣
where hf is a face enthalpy. The enthalpy acts as the advected quantity and an interpolation scheme (e.g. upwind) can be used to compute the face value. This kernel adds the face contribution for each face f to the right hand side and matrix.
The face mass flux (ρu⋅n)RC is provided by the RhieChowMassFlux object which uses pressure gradients and the discrete momentum equation to compute face velocities and mass fluxes. For more information on the expression that is used, see SIMPLE.
Input Parameters
- advected_interp_methodupwindThe interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'.
Default:upwind
C++ Type:MooseEnum
Options:average, upwind, sou, min_mod, vanLeer, quick, venkatakrishnan, skewness-corrected
Controllable:No
Description:The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'.
- advected_quantityenthalpyThe advected quantity
Default:enthalpy
C++ Type:MooseEnum
Options:enthalpy, temperature
Controllable:No
Description:The advected quantity
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- cpConstant specific heat value
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:Constant specific heat value
- force_boundary_executionFalseWhether to force execution of this object on all external boundaries.
Default:False
C++ Type:bool
Controllable:No
Description:Whether to force execution of this object on all external boundaries.
- 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
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Controllable:No
Description:The seed for the master random number generator
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
- ghost_layers1The number of layers of elements to ghost.
Default:1
C++ Type:unsigned short
Controllable:No
Description:The number of layers of elements to ghost.
- use_point_neighborsFalseWhether to use point neighbors, which introduces additional ghosting to that used for simple face neighbors.
Default:False
C++ Type:bool
Controllable:No
Description:Whether to use point neighbors, which introduces additional ghosting to that used for simple face neighbors.
Parallel Ghosting Parameters
Input Files
- (modules/navier_stokes/test/tests/finite_volume/wcns/enthalpy_equation/1d_test_h.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/cht/flow-around-square-linear.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-boussinesq-transient.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-boussinesq.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d-heated/fluid.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/steady-transient-compare/common-blocks.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/natural_convection/linear_segregated/2d/diff_heated_cavity_linear_segregated.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/cht/flow-around-square-linear-fluidonly.i)
- (modules/navier_stokes/test/tests/finite_volume/wcns/enthalpy_equation/1d_test_h_fp.i)
- (modules/navier_stokes/test/tests/finite_volume/wcns/enthalpy_equation/enthalpy_equation.i)
advected_quantity
Default:enthalpy
C++ Type:MooseEnum
Options:enthalpy, temperature
Controllable:No
Description:The advected quantity
(modules/navier_stokes/test/tests/finite_volume/wcns/enthalpy_equation/1d_test_h.i)
L = 30
nx = 600
bulk_u = 0.01
q_source = 50000.
A_cp = 976.78
B_cp = 1.0634
T_in = 860.
p_ref = 101325.0
rho = 2000.
advected_interp_method = 'upwind'
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 1
xmin = 0
xmax = ${L}
nx = ${nx}
[]
allow_renumbering = false
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = ${advected_interp_method}
u = vel_x
[]
[Problem]
linear_sys_names = 'u_system pressure_system energy_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = RhieChowMassFlux
u = vel_x
pressure = pressure
rho = 'rho'
p_diffusion_kernel = p_diffusion
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
solver_sys = u_system
initial_condition = ${bulk_u}
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = ${p_ref}
[]
[h]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = ${fparse 860.*1900.}
[]
[]
[AuxVariables]
[rho_var]
type = MooseLinearVariableFVReal
[]
[cp_var]
type = MooseLinearVariableFVReal
[]
[mu_var]
type = MooseLinearVariableFVReal
[]
[k_var]
type = MooseLinearVariableFVReal
[]
[alpha_var]
type = MooseLinearVariableFVReal
[]
[T]
type = MooseLinearVariableFVReal
initial_condition = 860.
[]
[h_aux]
type = MooseLinearVariableFVReal
[]
[]
[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'
[]
[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
[]
[temp_advection]
type = LinearFVEnergyAdvection
variable = h
[]
[source]
type = LinearFVSource
variable = h
source_density = source_func
[]
[]
[LinearFVBCs]
[inlet_u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_x
functor = ${bulk_u} #${bulk_u} #'fully_developed_velocity'
[]
[inlet_h]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = h
boundary = 'left'
functor = 'h_from_p_T'
[]
[inlet_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
boundary = 'left'
functor = ${T_in}
[]
[outlet_p]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'right'
variable = pressure
functor = ${p_ref}
[]
[outlet_h]
type = LinearFVAdvectionDiffusionOutflowBC
variable = h
use_two_term_expansion = false
boundary = 'right'
[]
[outlet_u]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_x
use_two_term_expansion = false
boundary = 'right'
[]
[]
[Functions]
[source_func]
type = ParsedFunction
expression = ${q_source}
[]
[T_analytical]
type = ParsedFunction
expression = ${fparse (-A_cp+sqrt(A_cp^2-2*B_cp*(-q_source/rho/bulk_u*L-A_cp*T_in-B_cp/2*T_in*T_in)))/B_cp}
[]
[]
[FunctorMaterials]
[enthalpy_material]
type = LinearFVEnthalpyFunctorMaterial
pressure = ${p_ref}
T_fluid = T
h = h
h_from_p_T_functor = h_from_p_T_functor
T_from_p_h_functor = T_from_p_h_functor
[]
[h_from_p_T_functor]
type = ParsedFunctorMaterial
property_name = 'h_from_p_T_functor'
functor_names = 'T'
expression = '${A_cp}*T+${B_cp}/2*(T^2)'
[]
[T_from_p_h_functor]
type = ParsedFunctorMaterial
property_name = 'T_from_p_h_functor'
functor_names = 'h'
expression = '(-${A_cp}+sqrt(${A_cp}^2+2*h*${B_cp}))/${B_cp}'
[]
[rho]
type = ADParsedFunctorMaterial
property_name = 'rho'
functor_names = 'T'
expression = ${rho}
[]
[cp]
type = ADParsedFunctorMaterial
property_name = 'cp'
functor_names = 'T'
expression = '${A_cp}+${B_cp}*T'
[]
[mu]
type = ADParsedFunctorMaterial
property_name = 'mu'
functor_names = 'T'
expression = '4.5e-3'
[]
[k]
type = ADParsedFunctorMaterial
property_name = 'k'
functor_names = 'T'
expression = 0.7
[]
[]
[AuxKernels]
[rho_out]
type = FunctorAux
functor = 'rho'
variable = 'rho_var'
execute_on = 'NONLINEAR'
[]
[cp_out]
type = FunctorAux
functor = 'cp'
variable = 'cp_var'
execute_on = 'NONLINEAR'
[]
[mu_out]
type = FunctorAux
functor = 'mu'
variable = 'mu_var'
execute_on = 'NONLINEAR'
[]
[k_out]
type = FunctorAux
functor = 'k'
variable = 'k_var'
execute_on = 'NONLINEAR'
[]
[T_from_h_functor_aux]
type = FunctorAux
functor = 'T_from_p_h'
variable = 'T'
execute_on = 'NONLINEAR'
[]
[h_from_T_functor_aux]
type = FunctorAux
functor = 'h_from_p_T'
variable = 'h_aux'
execute_on = 'NONLINEAR'
[]
[]
[Postprocessors]
[T_out_sim]
type = ElementalVariableValue
variable = T
elementid = ${fparse nx-1}
[]
[T_out_analytic]
type = FunctionValuePostprocessor
function = T_analytical
[]
[]
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-12
pressure_l_abs_tol = 1e-12
energy_l_abs_tol = 1e-12
momentum_l_tol = 0
pressure_l_tol = 0
energy_l_tol = 0
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system'
pressure_system = 'pressure_system'
energy_system = 'energy_system'
momentum_equation_relaxation = 0.7
pressure_variable_relaxation = 0.3
energy_equation_relaxation = 0.95
num_iterations = 100
pressure_absolute_tolerance = 1e-8
momentum_absolute_tolerance = 1e-8
energy_absolute_tolerance = 1e-6
print_fields = false
momentum_l_max_its = 200
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
[]
[Outputs]
[out]
type = CSV
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/cht/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/channel-flow/linear-segregated/2d/2d-boussinesq-transient.i)
mu = 2.6
rho = 1.0
advected_interp_method = 'average'
cp = 300
k = 10
alpha_b = 1e-4
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '1.'
dy = '0.2'
ix = '10'
iy = '5'
[]
[]
[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
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
initial_condition = 0.5
solver_sys = u_system
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
initial_condition = 0.0
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = 0.2
[]
[T]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = 300
[]
[]
[LinearFVKernels]
[u_time]
type = LinearFVTimeDerivative
variable = vel_x
factor = ${rho}
[]
[v_time]
type = LinearFVTimeDerivative
variable = vel_y
factor = ${rho}
[]
[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'
[]
[u_boussinesq]
type = LinearFVMomentumBoussinesq
variable = vel_x
rho = ${rho}
gravity = '0 -9.81 0'
alpha_name = ${alpha_b}
ref_temperature = 300.0
T_fluid = T
momentum_component = 'x'
[]
[v_boussinesq]
type = LinearFVMomentumBoussinesq
variable = vel_y
rho = ${rho}
gravity = '0 -9.81 0'
alpha_name = ${alpha_b}
ref_temperature = 300.0
T_fluid = T
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_time]
type = LinearFVTimeDerivative
variable = T
factor = ${fparse rho*cp}
[]
[h_advection]
type = LinearFVEnergyAdvection
variable = T
advected_quantity = temperature
cp = ${cp}
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
[]
[conduction]
type = LinearFVDiffusion
variable = T
diffusion_coeff = ${k}
use_nonorthogonal_correction = false
[]
[]
[FunctorMaterials]
[constant_functors]
type = GenericFunctorMaterial
prop_names = 'cp alpha_b'
prop_values = '${cp} ${alpha_b}'
[]
[]
[LinearFVBCs]
[inlet-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_x
functor = '1.1'
[]
[inlet-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_y
functor = '0.0'
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
variable = vel_x
functor = 0.0
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
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_top_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
functor = 300.0
boundary = 'left top'
[]
[bottom_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
functor = wall-temperature
boundary = bottom
[]
[outlet_T]
type = LinearFVAdvectionDiffusionOutflowBC
variable = T
use_two_term_expansion = false
boundary = right
[]
[]
[Functions]
[wall-temperature]
type = ParsedFunction
expression = '350 + 50 * sin(6.28*t)'
[]
[]
[Executioner]
type = PIMPLE
momentum_l_abs_tol = 1e-12
pressure_l_abs_tol = 1e-12
energy_l_abs_tol = 1e-12
momentum_l_tol = 1e-12
pressure_l_tol = 1e-12
energy_l_tol = 1e-12
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system v_system'
pressure_system = 'pressure_system'
energy_system = 'energy_system'
momentum_equation_relaxation = 0.8
pressure_variable_relaxation = 0.3
energy_equation_relaxation = 0.9
num_iterations = 100
pressure_absolute_tolerance = 1e-11
momentum_absolute_tolerance = 1e-11
energy_absolute_tolerance = 1e-11
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
dt = 0.01
num_steps = 6
num_piso_iterations = 0
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d/2d-boussinesq.i)
mu = 2.6
rho = 1.0
advected_interp_method = 'average'
cp = 300
k = 20
alpha_b = 1e-4
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '1.'
dy = '0.2'
ix = '10'
iy = '5'
[]
[]
[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
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
initial_condition = 0.5
solver_sys = u_system
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
initial_condition = 0.0
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = 0.2
[]
[T]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = 300
[]
[]
[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'
[]
[u_boussinesq]
type = LinearFVMomentumBoussinesq
variable = vel_x
rho = ${rho}
gravity = '0 -9.81 0'
alpha_name = ${alpha_b}
ref_temperature = 300.0
T_fluid = T
momentum_component = 'x'
[]
[v_boussinesq]
type = LinearFVMomentumBoussinesq
variable = vel_y
rho = ${rho}
gravity = '0 -9.81 0'
alpha_name = ${alpha_b}
ref_temperature = 300.0
T_fluid = T
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
advected_quantity = temperature
cp = ${cp}
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
[]
[conduction]
type = LinearFVDiffusion
variable = T
diffusion_coeff = ${k}
use_nonorthogonal_correction = false
[]
[]
[FunctorMaterials]
[constant_functors]
type = GenericFunctorMaterial
prop_names = 'cp alpha_b'
prop_values = '${cp} ${alpha_b}'
[]
[]
[LinearFVBCs]
[inlet-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_x
functor = '1.1'
[]
[inlet-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_y
functor = '0.0'
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
variable = vel_x
functor = 0.0
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
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_top_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
functor = ${fparse 300.0}
boundary = 'left top'
[]
[bottom_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
functor = ${fparse 350.0}
boundary = bottom
[]
[outlet_T]
type = LinearFVAdvectionDiffusionOutflowBC
variable = T
use_two_term_expansion = false
boundary = right
[]
[]
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-10
pressure_l_abs_tol = 1e-10
energy_l_abs_tol = 1e-10
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 = 0.9
pressure_variable_relaxation = 0.3
num_iterations = 200
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
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/linear-segregated/2d-heated/fluid.i)
mu = 2.6
rho = 1.0
advected_interp_method = 'upwind'
cp = 1000
k = 1
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '0.25 0.25'
dy = '0.2'
ix = '5 5'
iy = '5'
subdomain_id = '0 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
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
initial_condition = 0.5
solver_sys = u_system
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
initial_condition = 0.0
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = 0.2
[]
[T_fluid]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = 300
[]
[]
[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
[]
[heat_exchange]
type = LinearFVVolumetricHeatTransfer
variable = T_fluid
h_solid_fluid = 100
T_fluid = T_fluid
T_solid = T_solid
is_solid = false
block = 1
[]
[]
[FunctorMaterials]
[constant_functors]
type = GenericFunctorMaterial
prop_names = 'cp'
prop_values = '${cp}'
[]
[]
[LinearFVBCs]
[inlet-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_x
functor = '1.1'
[]
[inlet-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_y
functor = '0.0'
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
variable = vel_x
functor = 0.0
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
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_wall_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T_fluid
functor = 300
boundary = 'left top bottom'
[]
[outlet_T]
type = LinearFVAdvectionDiffusionOutflowBC
variable = T_fluid
use_two_term_expansion = false
boundary = right
[]
[]
[AuxVariables]
[T_solid]
type = MooseLinearVariableFVReal
initial_condition = 300
block = 1
[]
[]
[MultiApps]
inactive = 'solid'
[solid]
type = FullSolveMultiApp
input_files = solid.i
execute_on = timestep_begin
no_restore = true
[]
[]
[Transfers]
inactive = 'from_solid to_solid'
[from_solid]
type = MultiAppGeneralFieldShapeEvaluationTransfer
from_multi_app = solid
source_variable = 'T_solid'
variable = 'T_solid'
execute_on = timestep_begin
from_blocks = 1
[]
[to_solid]
type = MultiAppGeneralFieldShapeEvaluationTransfer
to_multi_app = solid
source_variable = 'T_fluid'
variable = 'T_fluid'
execute_on = timestep_begin
to_blocks = 1
[]
[]
[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 = 0.9
pressure_variable_relaxation = 0.3
num_iterations = 20
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/channel-flow/linear-segregated/steady-transient-compare/common-blocks.i)
mu = 2.6
rho = 1.0
advected_interp_method = 'average'
cp = 300
k = 20
alpha_b = 1e-4
[Mesh]
[mesh]
type = CartesianMeshGenerator
dim = 2
dx = '1.'
dy = '0.2'
ix = '10'
iy = '5'
[]
[]
[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
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
initial_condition = 0.5
solver_sys = u_system
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
initial_condition = 0.0
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = 0.2
[]
[T]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = 300
[]
[]
[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'
[]
[u_boussinesq]
type = LinearFVMomentumBoussinesq
variable = vel_x
rho = ${rho}
gravity = '0 -9.81 0'
alpha_name = ${alpha_b}
ref_temperature = 300.0
T_fluid = T
momentum_component = 'x'
[]
[v_boussinesq]
type = LinearFVMomentumBoussinesq
variable = vel_y
rho = ${rho}
gravity = '0 -9.81 0'
alpha_name = ${alpha_b}
ref_temperature = 300.0
T_fluid = T
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
advected_quantity = temperature
cp = ${cp}
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
[]
[conduction]
type = LinearFVDiffusion
variable = T
diffusion_coeff = ${k}
use_nonorthogonal_correction = false
[]
[]
[FunctorMaterials]
[constant_functors]
type = GenericFunctorMaterial
prop_names = 'cp alpha_b'
prop_values = '${cp} ${alpha_b}'
[]
[]
[LinearFVBCs]
[inlet-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_x
functor = '1.1'
[]
[inlet-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_y
functor = '0.0'
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
variable = vel_x
functor = 0.0
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'top bottom'
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_top_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
functor = ${fparse 300.0}
boundary = 'left top'
[]
[bottom_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
functor = ${fparse 350.0}
boundary = bottom
[]
[outlet_T]
type = LinearFVAdvectionDiffusionOutflowBC
variable = T
use_two_term_expansion = false
boundary = right
[]
[]
[Outputs]
exodus = true
execute_on = FINAL
[]
(modules/navier_stokes/test/tests/finite_volume/ins/natural_convection/linear_segregated/2d/diff_heated_cavity_linear_segregated.i)
################################################################################
# MATERIAL PROPERTIES
################################################################################
rho = 3279.
T_0 = 875.0
mu = 1.
k_cond = 38.0
cp = 640.
alpha = 3.26e-4
walls = 'right left top bottom'
[GlobalParams]
rhie_chow_user_object = 'ins_rhie_chow_interpolator'
advected_interp_method = 'upwind'
u = superficial_vel_x
v = superficial_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 = 0
xmax = 1
ymin = 0
ymax = 1
nx = 30
ny = 30
[]
[]
################################################################################
# EQUATIONS: VARIABLES, KERNELS & BCS
################################################################################
[UserObjects]
[ins_rhie_chow_interpolator]
type = RhieChowMassFlux
u = superficial_vel_x
v = superficial_vel_y
pressure = pressure
rho = ${rho}
p_diffusion_kernel = p_diffusion
[]
[]
[Variables]
[superficial_vel_x]
type = MooseLinearVariableFVReal
solver_sys = u_system
[]
[superficial_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 = superficial_vel_x
mu = ${mu}
momentum_component = 'x'
use_nonorthogonal_correction = true
[]
[u_pressure]
type = LinearFVMomentumPressure
variable = superficial_vel_x
pressure = pressure
momentum_component = 'x'
[]
[u_buoyancy]
type = LinearFVMomentumBoussinesq
variable = superficial_vel_x
T_fluid = T_fluid
gravity = '0 -9.81 0'
rho = ${rho}
ref_temperature = ${T_0}
alpha_name = ${alpha}
momentum_component = 'x'
[]
[v_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = superficial_vel_y
mu = ${mu}
momentum_component = 'y'
use_nonorthogonal_correction = true
[]
[v_pressure]
type = LinearFVMomentumPressure
variable = superficial_vel_y
pressure = pressure
momentum_component = 'y'
[]
[v_buoyancy]
type = LinearFVMomentumBoussinesq
variable = superficial_vel_y
T_fluid = T_fluid
gravity = '0 -9.81 0'
rho = ${rho}
ref_temperature = ${T_0}
alpha_name = ${alpha}
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
[]
####### FUEL ENERGY EQUATION #######
[heat_advection]
type = LinearFVEnergyAdvection
variable = T_fluid
advected_quantity = temperature
cp = ${cp}
[]
[conduction]
type = LinearFVDiffusion
variable = T_fluid
diffusion_coeff = ${fparse k_cond}
[]
[]
[LinearFVBCs]
[no-slip-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = superficial_vel_x
boundary = ${walls}
functor = 0
[]
[no-slip-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = superficial_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 = LinearFVAdvectionDiffusionExtrapolatedBC
variable = T_fluid
boundary = 'top bottom'
use_two_term_expansion = false
[]
[pressure-extrapolation]
type = LinearFVExtrapolatedPressureBC
boundary = ${walls}
variable = pressure
use_two_term_expansion = false
[]
[]
[FunctorMaterials]
[constant_functors]
type = GenericFunctorMaterial
prop_names = 'cp alpha_b'
prop_values = '${cp} ${alpha}'
[]
[]
################################################################################
# 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.3
pressure_variable_relaxation = 0.7
energy_equation_relaxation = 0.95
num_iterations = 3000
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 = '0.5 0.0 0.0'
# momentum_petsc_options = '-ksp_monitor'
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'
[]
################################################################################
# SIMULATION OUTPUTS
################################################################################
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/cht/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/wcns/enthalpy_equation/1d_test_h_fp.i)
L = 30
nx = 600
bulk_u = 0.01
p_ref = 101325.0
T_in = 860.
q_source = 20000.
advected_interp_method = 'upwind'
[Mesh]
[gmg]
type = GeneratedMeshGenerator
dim = 1
xmin = 0
xmax = ${L}
nx = ${nx}
[]
allow_renumbering = false
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = ${advected_interp_method}
u = vel_x
[]
[Problem]
linear_sys_names = 'u_system pressure_system energy_system'
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = RhieChowMassFlux
u = vel_x
pressure = pressure
rho = 'rho'
p_diffusion_kernel = p_diffusion
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
solver_sys = u_system
initial_condition = ${bulk_u}
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = ${p_ref}
[]
[h]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = ${fparse 860.*240.}
[]
[]
[AuxVariables]
[rho_var]
type = MooseLinearVariableFVReal
[]
[cp_var]
type = MooseLinearVariableFVReal
[]
[mu_var]
type = MooseLinearVariableFVReal
[]
[k_var]
type = MooseLinearVariableFVReal
[]
[alpha_var]
type = MooseLinearVariableFVReal
[]
[T]
type = MooseLinearVariableFVReal
initial_condition = ${T_in}
[]
[h_aux]
type = MooseLinearVariableFVReal
[]
[]
[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'
[]
[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
[]
[temp_advection]
type = LinearFVEnergyAdvection
variable = h
[]
[source]
type = LinearFVSource
variable = h
source_density = source_func
[]
[]
[LinearFVBCs]
[inlet_u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_x
functor = ${bulk_u}
[]
[inlet_h]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = h
boundary = 'left'
functor = 'h_from_p_T'
[]
[inlet_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
boundary = 'left'
functor = ${T_in}
[]
[outlet_p]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'right'
variable = pressure
functor = ${p_ref}
[]
[outlet_h]
type = LinearFVAdvectionDiffusionOutflowBC
variable = h
use_two_term_expansion = false
boundary = 'right'
[]
[outlet_u]
type = LinearFVAdvectionDiffusionOutflowBC
variable = vel_x
use_two_term_expansion = false
boundary = 'right'
[]
[]
[FluidProperties]
[lead]
type = LeadFluidProperties
[]
[]
[FunctorMaterials]
[enthalpy_material]
type = LinearFVEnthalpyFunctorMaterial
pressure = ${p_ref}
T_fluid = T
h = h
fp = lead
[]
[fluid_props_to_mat_props]
type = GeneralFunctorFluidProps
fp = lead
pressure = ${p_ref}
T_fluid = 'T'
speed = 1
porosity = 1
characteristic_length = 1
[]
[source_func]
type = ADParsedFunctorMaterial
property_name = source_func
functor_names = 'rho'
expression = ${q_source}
[]
[]
[AuxKernels]
[rho_out]
type = FunctorAux
functor = 'rho'
variable = 'rho_var'
execute_on = 'NONLINEAR'
[]
[cp_out]
type = FunctorAux
functor = 'cp'
variable = 'cp_var'
execute_on = 'NONLINEAR'
[]
[mu_out]
type = FunctorAux
functor = 'mu'
variable = 'mu_var'
execute_on = 'NONLINEAR'
[]
[k_out]
type = FunctorAux
functor = 'k'
variable = 'k_var'
execute_on = 'NONLINEAR'
[]
[T_from_h_functor_aux]
type = FunctorAux
functor = 'T_from_p_h'
variable = 'T'
execute_on = 'NONLINEAR'
[]
[h_from_T_functor_aux]
type = FunctorAux
functor = 'h_from_p_T'
variable = 'h_aux'
execute_on = 'NONLINEAR'
[]
[]
[Postprocessors]
[T_out_sim]
type = ElementalVariableValue
variable = T
elementid = ${fparse nx-1}
[]
[]
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-12
pressure_l_abs_tol = 1e-12
energy_l_abs_tol = 1e-12
momentum_l_tol = 0
pressure_l_tol = 0
energy_l_tol = 0
rhie_chow_user_object = 'rc'
momentum_systems = 'u_system'
pressure_system = 'pressure_system'
energy_system = 'energy_system'
momentum_equation_relaxation = 0.7
pressure_variable_relaxation = 0.3
energy_equation_relaxation = 0.95
num_iterations = 100
pressure_absolute_tolerance = 1e-8
momentum_absolute_tolerance = 1e-8
energy_absolute_tolerance = 1e-6
print_fields = false
momentum_l_max_its = 200
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
[]
[Outputs]
[out]
type = CSV
[]
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/enthalpy_equation/enthalpy_equation.i)
H = 0.015 #halfwidth of the channel, 10 cm of channel height
L = 1
bulk_u = 0.01
p_ref = 101325.0
advected_interp_method = 'upwind'
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${L}
ymin = -${H}
ymax = ${H}
nx = 30
ny = 15
[]
[]
[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
[]
[]
[Variables]
[vel_x]
type = MooseLinearVariableFVReal
solver_sys = u_system
initial_condition = ${bulk_u}
[]
[vel_y]
type = MooseLinearVariableFVReal
solver_sys = v_system
initial_condition = 0
[]
[pressure]
type = MooseLinearVariableFVReal
solver_sys = pressure_system
initial_condition = ${p_ref}
[]
[h]
type = MooseLinearVariableFVReal
solver_sys = energy_system
initial_condition = 44000 # 1900 is an approx of cp(T)
[]
[]
[AuxVariables]
[rho_var]
type = MooseLinearVariableFVReal
[]
[cp_var]
type = MooseLinearVariableFVReal
[]
[mu_var]
type = MooseLinearVariableFVReal
[]
[k_var]
type = MooseLinearVariableFVReal
[]
[T]
type = MooseLinearVariableFVReal
initial_condition = 777.
[]
[]
[LinearFVKernels]
[u_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_x
mu = 'mu'
momentum_component = 'x'
use_nonorthogonal_correction = false
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
u = vel_x
v = vel_y
[]
[u_pressure]
type = LinearFVMomentumPressure
variable = vel_x
pressure = pressure
momentum_component = 'x'
[]
[v_advection_stress]
type = LinearWCNSFVMomentumFlux
variable = vel_y
mu = 'mu'
momentum_component = 'y'
use_nonorthogonal_correction = false
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
u = vel_x
v = vel_y
[]
[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
[]
[temp_conduction]
type = LinearFVDiffusion
diffusion_coeff = 'alpha'
variable = h
[]
[temp_advection]
type = LinearFVEnergyAdvection
variable = h
advected_interp_method = ${advected_interp_method}
rhie_chow_user_object = 'rc'
[]
[]
[LinearFVBCs]
[inlet_u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_x
functor = ${bulk_u} #${bulk_u} #'fully_developed_velocity'
[]
[inlet-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'left'
variable = vel_y
functor = 0
[]
[inlet_h]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = h
boundary = 'left'
functor = h_from_p_T # ${fparse 1900.*860.}
[]
[inlet_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
boundary = 'left'
functor = 860.
[]
[walls-u]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = vel_x
boundary = 'top bottom'
functor = 0.
[]
[walls-v]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = vel_y
boundary = 'top bottom'
functor = 0.
[]
[walls_h]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = h
boundary = 'top bottom'
functor = h_from_p_T # ${fparse 1900. * 950}
[]
[walls_T]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
variable = T
boundary = 'top bottom'
functor = 950.
[]
[walls_p]
type = LinearFVExtrapolatedPressureBC
boundary = 'top bottom'
variable = pressure
use_two_term_expansion = false
[]
[outlet_p]
type = LinearFVAdvectionDiffusionFunctorDirichletBC
boundary = 'right'
variable = pressure
functor = ${p_ref}
[]
[outlet_h]
type = LinearFVAdvectionDiffusionOutflowBC
variable = h
use_two_term_expansion = false
boundary = 'right'
[]
[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
[]
[]
[FluidProperties]
[lead]
type = LeadFluidProperties
[]
[]
[FunctorMaterials]
[fluid_props_to_mat_props]
type = GeneralFunctorFluidProps
fp = lead
pressure = ${p_ref}
T_fluid = 'T'
speed = 1
porosity = 1
characteristic_length = 1
[]
[alpha]
type = ADParsedFunctorMaterial
property_name = 'alpha'
functor_names = 'k cp'
expression = 'k/cp'
[]
[enthalpy_material]
type = LinearFVEnthalpyFunctorMaterial
pressure = ${p_ref}
T_fluid = T
h = h
fp = lead
[]
[]
[AuxKernels]
[rho_out]
type = FunctorAux
functor = 'rho'
variable = 'rho_var'
execute_on = 'NONLINEAR'
[]
[cp_out]
type = FunctorAux
functor = 'cp'
variable = 'cp_var'
execute_on = 'NONLINEAR'
[]
[mu_out]
type = FunctorAux
functor = 'mu'
variable = 'mu_var'
execute_on = 'NONLINEAR'
[]
[k_out]
type = FunctorAux
functor = 'k'
variable = 'k_var'
execute_on = 'NONLINEAR'
[]
[T_from_h_functor]
type = FunctorAux
functor = 'T_from_p_h'
variable = 'T'
execute_on = 'NONLINEAR'
[]
[]
[Executioner]
type = SIMPLE
momentum_l_abs_tol = 1e-6
pressure_l_abs_tol = 1e-6
energy_l_abs_tol = 1e-8
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.7
pressure_variable_relaxation = 0.3
energy_equation_relaxation = 0.9
num_iterations = 200
pressure_absolute_tolerance = 1e-6
momentum_absolute_tolerance = 1e-6
energy_absolute_tolerance = 1e-6
print_fields = false
momentum_l_max_its = 1000
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
[]
[Outputs]
exodus = true
execute_on = 'TIMESTEP_BEGIN FINAL'
[]