- rhoDensity. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:Density. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable that this residual object operates on
INSFVEnergyTimeDerivative
Description
The INSFVEnergyTimeDerivative
kernel implements a time derivative for the domain given by
where is the material density, is the specific heat, is the fluid temperature and the second term on the left hand side corresponds to the strong forms of other kernels.
The Jacobian is computed with automatic differentiation.
Implementation
The derivative is obtained from the definition of the fluid energy. The isobaric and isochoric heat capacities being equal for incompressible fluids,
we take the partial derivative with regards to time:
The variation of the kinetic energy is not considered in this kernel.
The specific energy, , is currently approximated as .
Input Parameters
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Unit:(no unit assumed)
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- dh_dtdh_dtThe time derivative of the specific enthalpy. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
Default:dh_dt
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:The time derivative of the specific enthalpy. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- functorThe functor this kernel queries for the time derivative. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
C++ Type:MooseFunctorName
Unit:(no unit assumed)
Controllable:No
Description:The functor this kernel queries for the time derivative. A functor is any of the following: a variable, a functor material property, a function, a post-processor, or a number.
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
- use_interpolated_stateFalseFor the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:For the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
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>
Unit:(no unit assumed)
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>
Unit:(no unit assumed)
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>
Unit:(no unit assumed)
Controllable:No
Description:The extra tags for the vectors this Kernel should fill
- matrix_tagssystem timeThe tag for the matrices this Kernel should fill
Default:system time
C++ Type:MultiMooseEnum
Unit:(no unit assumed)
Options:nontime, system, time
Controllable:No
Description:The tag for the matrices this Kernel should fill
- vector_tagstimeThe tag for the vectors this Kernel should fill
Default:time
C++ Type:MultiMooseEnum
Unit:(no unit assumed)
Options:nontime, time
Controllable:No
Description:The tag for the vectors this Kernel should fill
Tagging Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable: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
Unit:(no unit assumed)
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
Unit:(no unit assumed)
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
Unit:(no unit assumed)
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
Unit:(no unit assumed)
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
Unit:(no unit assumed)
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/examples/solidification/gallium_melting.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/natural_convection/fuel_cavity.i)
- (modules/navier_stokes/examples/laser-welding/2d-fv.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-rc-transient.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/iks/flow-around-square/flow-around-square.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/solidification/solidification_no_advection.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/solidification/pipe_solidification.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/lid-driven/transient-lid-driven-with-energy.i)
Child Objects
(modules/navier_stokes/examples/solidification/gallium_melting.i)
##########################################################
# Simulation of Gallium Melting Experiment
# Ref: Gau, C., & Viskanta, R. (1986). Melting and solidification of a pure metal on a vertical wall.
# Key physics: melting/solidification, convective heat transfer, natural convection
##########################################################
mu = 1.81e-3
rho_solid = 6093
rho_liquid = 6093
k_solid = 32
k_liquid = 32
cp_solid = 381.5
cp_liquid = 381.5
L = 80160
alpha_b = 1.2e-4
T_solidus = 302.93
T_liquidus = '${fparse T_solidus + 0.1}'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
T_cold = 301.15
T_hot = 311.15
Nx = 100
Ny = 50
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 88.9e-3
ymin = 0
ymax = 63.5e-3
nx = ${Nx}
ny = ${Ny}
[]
[]
[AuxVariables]
[U]
type = MooseVariableFVReal
[]
[fl]
type = MooseVariableFVReal
initial_condition = 0.0
[]
[density]
type = MooseVariableFVReal
[]
[th_cond]
type = MooseVariableFVReal
[]
[cp_var]
type = MooseVariableFVReal
[]
[darcy_coef]
type = MooseVariableFVReal
[]
[fch_coef]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[mag]
type = VectorMagnitudeAux
variable = U
x = vel_x
y = vel_y
[]
[compute_fl]
type = NSLiquidFractionAux
variable = fl
temperature = T
T_liquidus = '${T_liquidus}'
T_solidus = '${T_solidus}'
execute_on = 'TIMESTEP_END'
[]
[rho_out]
type = FunctorAux
functor = 'rho_mixture'
variable = 'density'
[]
[th_cond_out]
type = FunctorAux
functor = 'k_mixture'
variable = 'th_cond'
[]
[cp_out]
type = FunctorAux
functor = 'cp_mixture'
variable = 'cp_var'
[]
[darcy_out]
type = FunctorAux
functor = 'Darcy_coefficient'
variable = 'darcy_coef'
[]
[fch_out]
type = FunctorAux
functor = 'Forchheimer_coefficient'
variable = 'fch_coef'
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0.0
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0.0
[]
[pressure]
type = INSFVPressureVariable
[]
[lambda]
family = SCALAR
order = FIRST
[]
[T]
type = INSFVEnergyVariable
initial_condition = '${T_cold}'
scaling = 1e-4
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = rho_mixture
[]
[mean_zero_pressure]
type = FVIntegralValueConstraint
variable = pressure
lambda = lambda
phi0 = 0.0
[]
[u_time]
type = INSFVMomentumTimeDerivative
variable = vel_x
rho = rho_mixture
momentum_component = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = rho_mixture
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[u_friction]
type = PINSFVMomentumFriction
variable = vel_x
momentum_component = 'x'
u = vel_x
v = vel_y
Darcy_name = 'Darcy_coeff'
Forchheimer_name = 'Forchheimer_coeff'
rho = ${rho_liquid}
mu = ${mu}
standard_friction_formulation = false
[]
[u_buoyancy]
type = INSFVMomentumBoussinesq
variable = vel_x
T_fluid = T
gravity = '0 -9.81 0'
rho = '${rho_liquid}'
ref_temperature = ${T_cold}
momentum_component = 'x'
[]
[u_gravity]
type = INSFVMomentumGravity
variable = vel_x
gravity = '0 -9.81 0'
rho = '${rho_liquid}'
momentum_component = 'x'
[]
[v_time]
type = INSFVMomentumTimeDerivative
variable = vel_y
rho = rho_mixture
momentum_component = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = rho_mixture
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[v_friction]
type = PINSFVMomentumFriction
variable = vel_y
momentum_component = 'y'
u = vel_x
v = vel_y
Darcy_name = 'Darcy_coeff'
Forchheimer_name = 'Forchheimer_coeff'
rho = ${rho_liquid}
mu = ${mu}
standard_friction_formulation = false
[]
[v_buoyancy]
type = INSFVMomentumBoussinesq
variable = vel_y
T_fluid = T
gravity = '0 -9.81 0'
rho = '${rho_liquid}'
ref_temperature = ${T_cold}
momentum_component = 'y'
[]
[v_gravity]
type = INSFVMomentumGravity
variable = vel_y
gravity = '0 -9.81 0'
rho = '${rho_liquid}'
momentum_component = 'y'
[]
[T_time]
type = INSFVEnergyTimeDerivative
variable = T
rho = rho_mixture
dh_dt = dh_dt
[]
[energy_advection]
type = INSFVEnergyAdvection
variable = T
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[energy_diffusion]
type = FVDiffusion
coeff = k_mixture
variable = T
[]
[energy_source]
type = NSFVPhaseChangeSource
variable = T
L = ${L}
liquid_fraction = fl
T_liquidus = ${T_liquidus}
T_solidus = ${T_solidus}
rho = 'rho_mixture'
[]
[]
[FVBCs]
[walls-u]
type = INSFVNoSlipWallBC
boundary = 'left right top bottom'
variable = vel_x
function = 0
[]
[walls-v]
type = INSFVNoSlipWallBC
boundary = 'left right top bottom'
variable = vel_y
function = 0
[]
[hot_wall]
type = FVDirichletBC
variable = T
value = '${T_hot}'
boundary = 'left'
[]
[cold_wall]
type = FVDirichletBC
variable = T
value = '${T_cold}'
boundary = 'right'
[]
[]
[FunctorMaterials]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
rho = rho_mixture
cp = cp_mixture
temperature = 'T'
[]
[eff_cp]
type = NSFVMixtureFunctorMaterial
phase_2_names = '${cp_solid} ${k_solid} ${rho_solid}'
phase_1_names = '${cp_liquid} ${k_liquid} ${rho_liquid}'
prop_names = 'cp_mixture k_mixture rho_mixture'
phase_1_fraction = fl
[]
[mushy_zone_resistance]
type = INSFVMushyPorousFrictionFunctorMaterial
liquid_fraction = 'fl'
mu = '${mu}'
rho_l = '${rho_liquid}'
dendrite_spacing_scaling = 1e-1
[]
[friction]
type = ADGenericVectorFunctorMaterial
prop_names = 'Darcy_coeff Forchheimer_coeff'
prop_values = 'darcy_coef darcy_coef darcy_coef fch_coef fch_coef fch_coef'
[]
[const_functor]
type = ADGenericFunctorMaterial
prop_names = 'alpha_b'
prop_values = '${alpha_b}'
[]
[]
[Executioner]
type = Transient
# Time-stepping parameters
start_time = 0.0
end_time = 200.0
num_steps = 2
[TimeStepper]
type = IterationAdaptiveDT
# Raise time step often but not by as much
# There's a rough spot for convergence near 10% fluid fraction
optimal_iterations = 15
growth_factor = 1.5
dt = 0.1
[]
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_rel_tol = 1e-6
nl_max_its = 30
line_search = 'none'
[]
[Postprocessors]
[ave_p]
type = ElementAverageValue
variable = 'pressure'
execute_on = 'INITIAL TIMESTEP_END'
[]
[ave_fl]
type = ElementAverageValue
variable = 'fl'
execute_on = 'INITIAL TIMESTEP_END'
[]
[ave_T]
type = ElementAverageValue
variable = 'T'
execute_on = 'INITIAL TIMESTEP_END'
[]
[]
[VectorPostprocessors]
[vel_x]
type = ElementValueSampler
variable = 'vel_x fl'
sort_by = 'x'
[]
[]
[Outputs]
exodus = true
csv = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/natural_convection/fuel_cavity.i)
# ========================================================================
# The purpose of this MOOSE scripts is to solve a 2-D axisymmetric
# problem with the following details:
# ------------------------------------------------------------------
# Physics: natural convection through a fluid and heat conduction
# in a solid and there is convective heat transfer from the
# solid to the liquid.
# ------------------------------------------------------------------
# Materials: the fluid is water and the solid is not specified.
# ------------------------------------------------------------------
# BCS: Inlet and outlet pressure with value of 0
# noslip conditions on the walls.
# Heat flux on the left wall with value of 40000 W/m^2
# ========================================================================
# ========================================================================
# Dimensions & Physical properties
# ========================================================================
Domain_length = 121.92e-2 # m
Solid_width = 0.7112e-3 # m
Liquid_width = 0.56261e-2 # m
mu = 0.00053157
rho = 987.27
k = 0.64247
k_solid = 15.0
cp = 4181.8
alpha_b = 210e-6
T_init = 300.0
input_heat_flux = 40000.0
# ========================================================================
# The main body of the script
# ========================================================================
[Mesh]
[cmg]
type = CartesianMeshGenerator
dim = 2
#dx = '0.7032625e-4 0.7112e-5'
dx = '${Liquid_width} ${Solid_width}'
ix = '10 3'
dy = '${fparse 1./5.*Domain_length} ${fparse 4./5.*Domain_length}'
iy = '30 10'
subdomain_id = '0 1
0 1'
[]
[interface]
type = SideSetsBetweenSubdomainsGenerator
input = 'cmg'
primary_block = 0
paired_block = 1
new_boundary = 'interface'
[]
[fluid_side]
type = BreakBoundaryOnSubdomainGenerator
input = 'interface'
boundaries = 'top bottom'
[]
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
block = 0
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
block = 0
initial_condition = 1e-6
[]
[vel_y]
type = INSFVVelocityVariable
block = 0
initial_condition = 1e-6
[]
[pressure]
type = INSFVPressureVariable
block = 0
[]
[T]
type = INSFVEnergyVariable
block = 0
initial_condition = ${T_init}
scaling = 1e-5
[]
[Ts]
type = INSFVEnergyVariable
block = 1
initial_condition = ${T_init}
scaling = 1e-3
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
rho = ${rho}
[]
[u_time]
type = INSFVMomentumTimeDerivative
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[u_buoyancy]
type = INSFVMomentumBoussinesq
variable = vel_x
T_fluid = T
gravity = '0 -9.81 0'
rho = ${rho}
ref_temperature = ${T_init}
momentum_component = 'x'
#alpha_name = ${alpha_b}
[]
[v_time]
type = INSFVMomentumTimeDerivative
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
#alpha_name = ${alpha_b}
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[v_buoyancy]
type = INSFVMomentumBoussinesq
variable = vel_y
T_fluid = T
gravity = '0 -9.81 0'
rho = ${rho}
ref_temperature = ${T_init}
momentum_component = 'y'
[]
[temp_time]
type = INSFVEnergyTimeDerivative
variable = T
rho = '${rho}'
dh_dt = dh_dt
[]
[temp_conduction]
type = FVDiffusion
coeff = 'k'
variable = T
[]
[temp_advection]
type = INSFVEnergyAdvection
variable = T
[]
[Ts_time]
type = INSFVEnergyTimeDerivative
variable = Ts
rho = '${rho}'
dh_dt = dh_solid_dt
[]
[solid_temp_conduction]
type = FVDiffusion
coeff = 'k_solid'
variable = Ts
[]
[]
[FVInterfaceKernels]
[convection]
type = FVConvectionCorrelationInterface
variable1 = T
variable2 = Ts
boundary = 'interface'
h = htc
T_solid = Ts
T_fluid = T
subdomain1 = 0
subdomain2 = 1
wall_cell_is_bulk = true
[]
[]
[FVBCs]
[walls_u]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'interface left bottom_to_0'
function = 0
[]
[walls_v]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'interface left bottom_to_0'
function = 0
[]
[outlet]
type = INSFVOutletPressureBC
variable = pressure
boundary = 'top_to_0'
function = 0.0
[]
[outlet_T]
type = NSFVOutflowTemperatureBC
variable = T
boundary = 'top_to_0'
u = vel_x
v = vel_y
rho = ${rho}
cp = '${cp}'
backflow_T = ${T_init}
[]
[Insulator]
type = FVNeumannBC
variable = 'T'
boundary = 'left'
value = 0.0
[]
[heater]
type = FVNeumannBC
variable = 'Ts'
boundary = 'right'
value = '${fparse input_heat_flux}'
[]
[Insulator_solid]
type = FVNeumannBC
variable = 'Ts'
boundary = 'top_to_1'
value = 0.0
[]
[inlet_T_1]
type = FVDirichletBC
variable = Ts
boundary = 'bottom_to_1'
value = ${T_init}
[]
[]
[AuxVariables]
[Ra]
type = INSFVScalarFieldVariable
initial_condition = 1000.0
[]
[htc]
type = INSFVScalarFieldVariable
initial_condition = 0.0
[]
[]
[AuxKernels]
[compute_Ra]
type = ParsedAux
variable = Ra
coupled_variables = 'T'
constant_names = 'g beta T_init width nu alpha'
constant_expressions = '9.81 ${alpha_b} ${T_init} ${Liquid_width} ${fparse mu/rho} ${fparse k/(rho*cp)}'
expression = 'g * beta * (T - T_init) * pow(width, 3) / (nu*alpha) + 1.0'
block = 0
[]
[htc]
type = ParsedAux
variable = htc
coupled_variables = 'Ra'
constant_names = 'Pr'
constant_expressions = '${fparse cp*mu/k}'
expression = '${k}* (0.68 + 0.67 * pow(Ra, 0.25)/pow(1 + pow(0.437/Pr, 9/16) ,4/9) )/ ${Liquid_width} '
block = 0
[]
[]
[FunctorMaterials]
[functor_constants]
type = ADGenericFunctorMaterial
prop_names = 'cp k k_solid'
prop_values = '${cp} ${k} ${k_solid}'
[]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
temperature = 'T'
rho = ${rho}
block = 0
[]
[ins_fv_solid]
type = INSFVEnthalpyFunctorMaterial
temperature = 'Ts'
rho = ${rho}
cp = ${cp}
h = h_solid
rho_h = rho_h_solid
block = 1
[]
[const_functor]
type = ADGenericFunctorMaterial
prop_names = 'alpha_b'
prop_values = '${alpha_b}'
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -sub_pc_factor_shift_type -ksp_gmres_restart'
petsc_options_value = ' lu NONZERO 200'
line_search = 'none'
[TimeStepper]
type = IterationAdaptiveDT
dt = 0.01
optimal_iterations = 20
iteration_window = 2
[]
nl_max_its = 30
nl_abs_tol = 1e-10
steady_state_detection = true
steady_state_tolerance = 1e-09
[]
[Postprocessors]
[max_T]
type = ADElementExtremeFunctorValue
functor = T
block = 0
[]
[max_Ts]
type = ADElementExtremeFunctorValue
functor = Ts
block = 1
[]
[]
[Outputs]
exodus = false
csv = true
[]
(modules/navier_stokes/examples/laser-welding/2d-fv.i)
period=.2e-4 # s
endtime=${fparse 3 * period} # s
timestep=${fparse period / 100} # s
surfacetemp=2700 # K
bottomtemp=2700 # K
sb=5.67e-8 # W/(m^2 K^4)
advected_interp_method='upwind'
velocity_interp_method='rc'
rho='rho'
mu='mu'
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -.7e-3 # m
xmax = 0.7e-3 # m
ymin = -.35e-3 # m
ymax = 0
nx = 75
ny = 20
displacements = 'disp_x disp_y'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
use_displaced_mesh = true
disp_x = disp_x
disp_y = disp_y
[]
[]
[Problem]
extra_tag_vectors = 'e_time e_advection e_conduction e_laser e_radiation e_mesh_advection'
[]
[AuxVariables]
[mu_out]
type = MooseVariableFVReal
[]
[e_time]
type = MooseVariableFVReal
[]
[e_advection]
type = MooseVariableFVReal
[]
[e_mesh_advection]
type = MooseVariableFVReal
[]
[e_conduction]
type = MooseVariableFVReal
[]
[e_laser]
type = MooseVariableFVReal
[]
[e_radiation]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[mu_out]
type = FunctorAux
functor = mu
variable = mu_out
execute_on = timestep_end
[]
[e_time]
variable = e_time
vector_tag = e_time
v = T
execute_on = 'timestep_end'
type = TagVectorAux
[]
[e_advection]
variable = e_advection
vector_tag = e_advection
v = T
execute_on = 'timestep_end'
type = TagVectorAux
[]
[e_mesh_advection]
variable = e_mesh_advection
vector_tag = e_mesh_advection
v = T
execute_on = 'timestep_end'
type = TagVectorAux
[]
[e_conduction]
variable = e_conduction
vector_tag = e_conduction
v = T
execute_on = 'timestep_end'
type = TagVectorAux
[]
[e_laser]
variable = e_laser
vector_tag = e_laser
v = T
execute_on = 'timestep_end'
type = TagVectorAux
[]
[e_radiation]
variable = e_radiation
vector_tag = e_radiation
v = T
execute_on = 'timestep_end'
type = TagVectorAux
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
[]
[vel_y]
type = INSFVVelocityVariable
[]
[T]
type = INSFVEnergyVariable
[]
[pressure]
type = INSFVPressureVariable
[]
[disp_x]
[]
[disp_y]
[]
[]
[ICs]
[T]
type = FunctionIC
variable = T
function = '${surfacetemp} + ((${surfacetemp} - ${bottomtemp}) / .35e-3) * y'
[]
[]
[Kernels]
[disp_x]
type = MatDiffusion
variable = disp_x
diffusivity = 1e6
[]
[disp_y]
type = MatDiffusion
variable = disp_y
diffusivity = 1e6
[]
[]
[FVKernels]
# pressure equation
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
use_displaced_mesh = true
boundaries_to_force = top
[]
# momentum equations
# u equation
[u_time]
type = INSFVMomentumTimeDerivative
variable = vel_x
rho = ${rho}
momentum_component = 'x'
use_displaced_mesh = true
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
momentum_component = 'x'
use_displaced_mesh = true
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
use_displaced_mesh = true
[]
[u_pressure]
type = INSFVMomentumPressureFlux
variable = vel_x
momentum_component = 'x'
pressure = pressure
use_displaced_mesh = true
[]
[u_mesh_advection_volumetric]
type = INSFVMomentumMeshAdvection
variable = vel_x
momentum_component = 'x'
rho = ${rho}
disp_x = disp_x
disp_y = disp_y
add_to_a = false
use_displaced_mesh = true
[]
# v equation
[v_time]
type = INSFVMomentumTimeDerivative
variable = vel_y
rho = ${rho}
momentum_component = 'y'
use_displaced_mesh = true
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
momentum_component = 'y'
use_displaced_mesh = true
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
use_displaced_mesh = true
[]
[v_pressure]
type = INSFVMomentumPressureFlux
variable = vel_y
momentum_component = 'y'
pressure = pressure
use_displaced_mesh = true
[]
[v_mesh_advection_volumetric]
type = INSFVMomentumMeshAdvection
variable = vel_y
momentum_component = 'y'
rho = ${rho}
disp_x = disp_x
disp_y = disp_y
add_to_a = false
use_displaced_mesh = true
[]
# energy equation
[temperature_time]
type = INSFVEnergyTimeDerivative
variable = T
rho = ${rho}
dh_dt = dh_dt
use_displaced_mesh = true
extra_vector_tags = 'e_time'
[]
[temperature_advection]
type = INSFVEnergyAdvection
variable = T
use_displaced_mesh = true
extra_vector_tags = 'e_advection'
[]
[temperature_conduction]
type = FVDiffusion
coeff = 'k'
variable = T
use_displaced_mesh = true
extra_vector_tags = 'e_conduction'
[]
[temperature_mesh_advection_volumetric]
type = INSFVMeshAdvection
variable = T
rho = ${rho}
disp_x = disp_x
disp_y = disp_y
advected_quantity = 'h'
use_displaced_mesh = true
extra_vector_tags = 'e_mesh_advection'
[]
[]
[FVBCs]
# momentum boundary conditions
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'bottom right left'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'bottom right left'
function = 0
[]
[vapor_recoil_x]
type = INSFVVaporRecoilPressureMomentumFluxBC
variable = vel_x
boundary = 'top'
momentum_component = 'x'
rc_pressure = rc_pressure
use_displaced_mesh = true
[]
[vapor_recoil_y]
type = INSFVVaporRecoilPressureMomentumFluxBC
variable = vel_y
boundary = 'top'
momentum_component = 'y'
rc_pressure = rc_pressure
use_displaced_mesh = true
[]
# energy boundary conditions
[T_cold]
type = FVDirichletBC
variable = T
boundary = 'bottom'
value = '${bottomtemp}'
[]
[radiation_flux]
type = FVFunctorRadiativeBC
variable = T
boundary = 'top'
emissivity = '1'
Tinfinity = 300
stefan_boltzmann_constant = ${sb}
use_displaced_mesh = true
extra_vector_tags = 'e_radiation'
[]
[weld_flux]
type = FVGaussianEnergyFluxBC
variable = T
boundary = 'top'
P0 = 159.96989792079225
R = 1.25e-4
x_beam_coord = '2e-4 * sin(t * 2 * pi / ${period})'
y_beam_coord = 0
z_beam_coord = 0
use_displaced_mesh = true
extra_vector_tags = 'e_laser'
[]
[]
[BCs]
# displacement boundary conditions
[x_no_disp]
type = DirichletBC
variable = disp_x
boundary = 'bottom'
value = 0
[]
[y_no_disp]
type = DirichletBC
variable = disp_y
boundary = 'bottom'
value = 0
[]
[displace_x_top]
type = INSADDisplaceBoundaryBC
boundary = 'top'
variable = 'disp_x'
velocity = 'vel'
component = 0
associated_subdomain = 0
[]
[displace_y_top]
type = INSADDisplaceBoundaryBC
boundary = 'top'
variable = 'disp_y'
velocity = 'vel'
component = 1
associated_subdomain = 0
[]
[displace_x_top_dummy]
type = INSADDummyDisplaceBoundaryIntegratedBC
boundary = 'top'
variable = 'disp_x'
velocity = 'vel'
component = 0
[]
[displace_y_top_dummy]
type = INSADDummyDisplaceBoundaryIntegratedBC
boundary = 'top'
variable = 'disp_y'
velocity = 'vel'
component = 1
[]
[]
[FunctorMaterials]
[steel]
type = AriaLaserWeld304LStainlessSteelFunctorMaterial
temperature = T
beta = 1e7
[]
[disp_vec_value_and_dot]
type = ADGenericVectorFunctorMaterial
prop_names = 'disp_vec'
prop_values = 'disp_x disp_y 0'
[]
[vel]
type = ADGenericVectorFunctorMaterial
prop_names = 'vel'
prop_values = 'vel_x vel_y 0'
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_type -mat_mffd_err'
petsc_options_value = 'lu NONZERO strumpack 1e-6'
[]
[]
[Executioner]
type = Transient
end_time = ${endtime}
dtmin = 1e-8
dtmax = ${timestep}
petsc_options = '-snes_converged_reason -ksp_converged_reason -options_left'
solve_type = 'PJFNK'
line_search = 'none'
nl_max_its = 12
l_max_its = 100
[TimeStepper]
type = IterationAdaptiveDT
optimal_iterations = 7
dt = ${timestep}
linear_iteration_ratio = 1e6
growth_factor = 1.1
[]
[]
[Outputs]
exodus = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[laser_flux]
type = TagVectorSum
vector = 'e_laser'
[]
[volume_rho_cp_dT]
type = TagVectorSum
vector = 'e_time'
[]
[conduction]
type = TagVectorSum
vector = 'e_conduction'
[]
[advection]
type = TagVectorSum
vector = 'e_advection'
[]
[mesh_advection]
type = TagVectorSum
vector = 'e_mesh_advection'
[]
[radiation]
type = TagVectorSum
vector = 'e_radiation'
[]
[total_sum]
type = ParsedPostprocessor
expression = 'laser_flux + volume_rho_cp_dT + advection + mesh_advection + conduction + radiation'
pp_names = 'laser_flux volume_rho_cp_dT advection mesh_advection conduction radiation'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/channel-flow/2d-rc-transient.i)
# Fluid properties
mu = 1.1
rho = 1.1
cp = 1.1
k = 1e-3
# Operating conditions
u_inlet = 1
T_inlet = 200
T_solid = 190
p_outlet = 10
h_fs = 0.01
# Numerical scheme
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 5
ymin = -1
ymax = 1
nx = 50
ny = 20
[]
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = ${u_inlet}
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 1e-12
[]
[pressure]
type = INSFVPressureVariable
[]
[T_fluid]
type = INSFVEnergyVariable
initial_condition = ${T_inlet}
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[u_time]
type = INSFVMomentumTimeDerivative
variable = vel_x
rho = ${rho}
momentum_component = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[v_time]
type = INSFVMomentumTimeDerivative
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[energy_time]
type = INSFVEnergyTimeDerivative
variable = T_fluid
rho = ${rho}
dh_dt = dh_dt
[]
[energy_advection]
type = INSFVEnergyAdvection
variable = T_fluid
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[energy_diffusion]
type = FVDiffusion
variable = T_fluid
coeff = ${k}
[]
[energy_convection]
type = PINSFVEnergyAmbientConvection
variable = T_fluid
is_solid = false
T_fluid = 'T_fluid'
T_solid = 'T_solid'
h_solid_fluid = 'h_cv'
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_x
function = '1'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_y
function = 0
[]
[inlet-T]
type = FVNeumannBC
variable = T_fluid
value = '${fparse u_inlet * rho * cp * T_inlet}'
boundary = 'left'
[]
[no-slip-u]
type = INSFVNoSlipWallBC
boundary = 'top'
variable = vel_x
function = 0
[]
[no-slip-v]
type = INSFVNoSlipWallBC
boundary = 'top'
variable = vel_y
function = 0
[]
[symmetry-u]
type = INSFVSymmetryVelocityBC
boundary = 'bottom'
variable = vel_x
u = vel_x
v = vel_y
mu = ${mu}
momentum_component = 'x'
[]
[symmetry-v]
type = INSFVSymmetryVelocityBC
boundary = 'bottom'
variable = vel_y
u = vel_x
v = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[symmetry-p]
type = INSFVSymmetryPressureBC
boundary = 'bottom'
variable = pressure
[]
[outlet_u]
type = INSFVMomentumAdvectionOutflowBC
variable = vel_x
u = vel_x
v = vel_y
boundary = 'right'
momentum_component = 'x'
rho = ${rho}
[]
[outlet_v]
type = INSFVMomentumAdvectionOutflowBC
variable = vel_y
u = vel_x
v = vel_y
boundary = 'right'
momentum_component = 'y'
rho = ${rho}
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = '${p_outlet}'
[]
[]
[FunctorMaterials]
[constants]
type = ADGenericFunctorMaterial
prop_names = 'h_cv T_solid'
prop_values = '${h_fs} ${T_solid}'
[]
[functor_constants]
type = ADGenericFunctorMaterial
prop_names = 'cp'
prop_values = '${cp}'
[]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
rho = ${rho}
temperature = 'T_fluid'
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
line_search = 'none'
nl_rel_tol = 7e-13
dt = 0.4
end_time = 0.8
[]
[Outputs]
exodus = true
csv = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/iks/flow-around-square/flow-around-square.i)
# Water properties
mu = 1.0E-3
rho = 1000.0
k = 0.598
cp = 4186
# Solid properties
cp_s = 830
rho_s = 1680
k_s = 3.5
# Other parameters
p_outlet = 0
u_inlet = -1e-4
h_conv = 50
[Mesh]
[generated_mesh]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 10
xmin = 0
ymin = 0
ymax = 0.1
xmax = 0.1
[]
[subdomain1]
input = generated_mesh
type = SubdomainBoundingBoxGenerator
block_name = subdomain1
bottom_left = '0.04 0.04 0'
block_id = 1
top_right = '0.06 0.06 0'
[]
[interface]
input = subdomain1
type = SideSetsBetweenSubdomainsGenerator
primary_block = 0
paired_block = 1
new_boundary = interface
[]
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
block = 0
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 1e-4
block = 0
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 1e-4
block = 0
[]
[pressure]
type = INSFVPressureVariable
block = 0
[]
[T]
type = INSFVEnergyVariable
initial_condition = 283.15
scaling = 1e-5
block = 0
[]
[Ts]
type = INSFVEnergyVariable
initial_condition = 333.15
scaling = 1e-5
block = 1
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
rho = ${rho}
block = 0
[]
[u_time]
type = INSFVMomentumTimeDerivative
variable = vel_x
rho = ${rho}
momentum_component = 'x'
block = 0
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = ${rho}
momentum_component = 'x'
block = 0
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
block = 0
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
block = 0
[]
[v_time]
type = INSFVMomentumTimeDerivative
variable = vel_y
rho = ${rho}
momentum_component = 'y'
block = 0
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
block = 0
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
block = 0
[]
[energy_time]
type = INSFVEnergyTimeDerivative
variable = T
rho = ${rho}
dh_dt = dh_dt
block = 0
[]
[temp_conduction]
type = FVDiffusion
coeff = 'k'
variable = T
block = 0
[]
[temp_advection]
type = INSFVEnergyAdvection
variable = T
block = 0
[]
[solid_energy_time]
type = INSFVEnergyTimeDerivative
variable = Ts
rho = ${rho_s}
dh_dt = dh_solid_dt
block = 1
[]
[solid_temp_conduction]
type = FVDiffusion
coeff = 'k_s'
variable = Ts
block = 1
[]
[]
[FVInterfaceKernels]
[convection]
type = FVConvectionCorrelationInterface
variable1 = T
variable2 = Ts
subdomain1 = 0
subdomain2 = 1
boundary = interface
h = ${h_conv}
T_solid = Ts
T_fluid = T
wall_cell_is_bulk = true
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'top'
variable = vel_x
function = 0
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'top'
variable = vel_y
function = ${u_inlet}
[]
[inlet_T]
type = FVDirichletBC
variable = T
boundary = 'top'
value = 283.15
[]
[no-slip-u]
type = INSFVNoSlipWallBC
boundary = 'left right interface'
variable = vel_x
function = 0
[]
[no-slip-v]
type = INSFVNoSlipWallBC
boundary = 'left right interface'
variable = vel_y
function = 0
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'bottom'
variable = pressure
function = '${p_outlet}'
[]
[]
[FunctorMaterials]
[functor_constants]
type = ADGenericFunctorMaterial
prop_names = 'cp k'
prop_values = '${cp} ${k}'
block = 0
[]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
temperature = 'T'
rho = ${rho}
block = 0
[]
[solid_functor_constants]
type = ADGenericFunctorMaterial
prop_names = 'cp_s k_s'
prop_values = '${cp_s} ${k_s}'
block = 1
[]
[solid_ins_fv]
type = INSFVEnthalpyFunctorMaterial
temperature = 'Ts'
rho = ${rho_s}
cp = ${cp_s}
block = 1
h = h_solid
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
petsc_options_value = 'asm 100 lu NONZERO'
line_search = 'none'
nl_rel_tol = 1e-8
dt = 10
end_time = 10
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/solidification/solidification_no_advection.i)
rho_solid = 1.0
rho_liquid = 1.0
k_solid = 0.03
k_liquid = 0.1
cp_solid = 1.0
cp_liquid = 1.0
T_liquidus = 260
T_solidus = 240
L = 1.0
T_hot = 300.0
T_cold = 200.0
N = 10
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = ${N}
ny = ${N}
[]
[]
[AuxVariables]
[fl]
type = MooseVariableFVReal
initial_condition = 1.0
[]
[density]
type = MooseVariableFVReal
[]
[th_cond]
type = MooseVariableFVReal
[]
[cp_var]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[compute_fl]
type = NSLiquidFractionAux
variable = fl
temperature = T
T_liquidus = '${T_liquidus}'
T_solidus = '${T_solidus}'
execute_on = 'TIMESTEP_END'
[]
[rho_out]
type = FunctorAux
functor = 'rho_mixture'
variable = 'density'
[]
[th_cond_out]
type = FunctorAux
functor = 'k_mixture'
variable = 'th_cond'
[]
[cp_out]
type = FunctorAux
functor = 'cp_mixture'
variable = 'cp_var'
[]
[]
[Variables]
[T]
type = INSFVEnergyVariable
initial_condition = '${T_hot}'
[]
[]
[FVKernels]
[T_time]
type = INSFVEnergyTimeDerivative
variable = T
rho = ${rho_liquid}
[]
[energy_diffusion]
type = FVDiffusion
coeff = 'k_mixture'
variable = T
[]
[energy_source]
type = NSFVPhaseChangeSource
variable = T
L = ${L}
liquid_fraction = fl
T_liquidus = ${T_liquidus}
T_solidus = ${T_solidus}
rho = 'rho_mixture'
[]
[]
[FVBCs]
[heated_wall]
type = FVDirichletBC
variable = T
value = '${T_hot}'
boundary = 'top'
[]
[cooled_wall]
type = FVDirichletBC
variable = T
value = '${T_cold}'
boundary = 'bottom'
[]
[]
[FunctorMaterials]
[eff_cp]
type = NSFVMixtureFunctorMaterial
phase_2_names = '${cp_solid} ${k_solid} ${rho_solid}'
phase_1_names = '${cp_liquid} ${k_liquid} ${rho_liquid}'
prop_names = 'cp_mixture k_mixture rho_mixture'
phase_1_fraction = fl
[]
[h]
type = INSFVEnthalpyFunctorMaterial
cp = ${cp_liquid}
temperature = T
rho = ${rho_liquid}
[]
[]
[Executioner]
type = Transient
dt = 0.5
end_time = 50.0
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_abs_tol = 1e-12
nl_max_its = 50
steady_state_detection = true
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/ins/solidification/pipe_solidification.i)
mu = 8.8871e-4
rho_solid = 997.561
rho_liquid = 997.561
k_solid = 0.6203
k_liquid = 0.6203
cp_solid = 4181.72
cp_liquid = 4181.72
L = 3e5
T_liquidus = 285
T_solidus = 280
advected_interp_method = 'average'
velocity_interp_method = 'rc'
U_inlet = '${fparse 0.5 * mu / rho_liquid / 0.5}'
T_inlet = 300.0
T_cold = 200.0
Nx = 30
Ny = 5
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Mesh]
coord_type = 'RZ'
rz_coord_axis = 'X'
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = '${fparse 0.5 * 1.0}'
nx = ${Nx}
ny = ${Ny}
bias_y = '${fparse 1 / 1.2}'
[]
[rename1]
type = RenameBoundaryGenerator
input = gen
old_boundary = 'left'
new_boundary = 'inlet'
[]
[rename2]
type = RenameBoundaryGenerator
input = rename1
old_boundary = 'right'
new_boundary = 'outlet'
[]
[rename3]
type = RenameBoundaryGenerator
input = rename2
old_boundary = 'bottom'
new_boundary = 'symmetry'
[]
[rename4]
type = RenameBoundaryGenerator
input = rename3
old_boundary = 'top'
new_boundary = 'wall'
[]
[rename5]
type = ParsedGenerateSideset
input = rename4
normal = '0 1 0'
combinatorial_geometry = 'x>2.0 & x<8.0 & y>0.49999'
new_sideset_name = 'cooled_wall'
[]
[]
[AuxVariables]
[U]
type = MooseVariableFVReal
[]
[fl]
type = MooseVariableFVReal
initial_condition = 1.0
[]
[density]
type = MooseVariableFVReal
[]
[th_cond]
type = MooseVariableFVReal
[]
[cp_var]
type = MooseVariableFVReal
[]
[darcy_coef]
type = MooseVariableFVReal
[]
[fch_coef]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[mag]
type = VectorMagnitudeAux
variable = U
x = vel_x
y = vel_y
[]
[compute_fl]
type = NSLiquidFractionAux
variable = fl
temperature = T
T_liquidus = '${T_liquidus}'
T_solidus = '${T_solidus}'
execute_on = 'TIMESTEP_END'
[]
[rho_out]
type = FunctorAux
functor = 'rho_mixture'
variable = 'density'
[]
[th_cond_out]
type = FunctorAux
functor = 'k_mixture'
variable = 'th_cond'
[]
[cp_out]
type = FunctorAux
functor = 'cp_mixture'
variable = 'cp_var'
[]
[darcy_out]
type = FunctorAux
functor = 'Darcy_coefficient'
variable = 'darcy_coef'
[]
[fch_out]
type = FunctorAux
functor = 'Forchheimer_coefficient'
variable = 'fch_coef'
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0.0
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0.0
[]
[pressure]
type = INSFVPressureVariable
[]
[T]
type = INSFVEnergyVariable
initial_condition = '${T_inlet}'
scaling = 1.0
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = rho_mixture
[]
[u_time]
type = INSFVMomentumTimeDerivative
variable = vel_x
rho = rho_mixture
momentum_component = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = rho_mixture
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = ${mu}
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[u_friction]
type = PINSFVMomentumFriction
variable = vel_x
momentum_component = 'x'
u = vel_x
v = vel_y
Darcy_name = 'Darcy_coeff'
Forchheimer_name = 'Forchheimer_coeff'
rho = ${rho_liquid}
mu = ${mu}
standard_friction_formulation = false
[]
[v_time]
type = INSFVMomentumTimeDerivative
variable = vel_y
rho = rho_mixture
momentum_component = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = rho_mixture
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[v_friction]
type = PINSFVMomentumFriction
variable = vel_y
momentum_component = 'y'
u = vel_x
v = vel_y
Darcy_name = 'Darcy_coeff'
Forchheimer_name = 'Forchheimer_coeff'
rho = ${rho_liquid}
mu = ${mu}
standard_friction_formulation = false
[]
[T_time]
type = INSFVEnergyTimeDerivative
variable = T
rho = rho_mixture
dh_dt = dh_dt
[]
[energy_advection]
type = INSFVEnergyAdvection
variable = T
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[energy_diffusion]
type = FVDiffusion
coeff = k_mixture
variable = T
[]
[energy_source]
type = NSFVPhaseChangeSource
variable = T
L = ${L}
liquid_fraction = fl
T_liquidus = ${T_liquidus}
T_solidus = ${T_solidus}
rho = 'rho_mixture'
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'inlet'
variable = vel_x
function = '${U_inlet}'
[]
[sym_u]
type = INSFVSymmetryVelocityBC
boundary = 'symmetry'
variable = vel_x
u = vel_x
v = vel_y
mu = ${mu}
momentum_component = 'x'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'inlet'
variable = vel_y
function = 0
[]
[walls-u]
type = INSFVNoSlipWallBC
boundary = 'wall'
variable = vel_x
function = 0
[]
[walls-v]
type = INSFVNoSlipWallBC
boundary = 'wall'
variable = vel_y
function = 0
[]
[sym_v]
type = INSFVSymmetryVelocityBC
boundary = 'symmetry'
variable = vel_y
u = vel_x
v = vel_y
mu = ${mu}
momentum_component = y
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'outlet'
variable = pressure
function = 0
[]
[sym_p]
type = INSFVSymmetryPressureBC
boundary = 'symmetry'
variable = pressure
[]
[sym_T]
type = INSFVSymmetryScalarBC
variable = T
boundary = 'symmetry'
[]
[cooled_wall]
type = FVFunctorDirichletBC
variable = T
functor = '${T_cold}'
boundary = 'cooled_wall'
[]
[]
[FunctorMaterials]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
rho = rho_mixture
cp = cp_mixture
temperature = 'T'
[]
[eff_cp]
type = NSFVMixtureFunctorMaterial
phase_2_names = '${cp_solid} ${k_solid} ${rho_solid}'
phase_1_names = '${cp_liquid} ${k_liquid} ${rho_liquid}'
prop_names = 'cp_mixture k_mixture rho_mixture'
phase_1_fraction = fl
[]
[mushy_zone_resistance]
type = INSFVMushyPorousFrictionFunctorMaterial
liquid_fraction = 'fl'
mu = '${mu}'
rho_l = '${rho_liquid}'
[]
[friction]
type = ADGenericVectorFunctorMaterial
prop_names = 'Darcy_coeff Forchheimer_coeff'
prop_values = 'darcy_coef darcy_coef darcy_coef fch_coef fch_coef fch_coef'
[]
[]
[Executioner]
type = Transient
dt = 5e3
end_time = 1e4
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_abs_tol = 1e-8
nl_max_its = 12
[]
[Postprocessors]
[average_T]
type = ElementAverageValue
variable = T
outputs = csv
execute_on = FINAL
[]
[]
[VectorPostprocessors]
[sat]
type = LineValueSampler
warn_discontinuous_face_values = false
start_point = '0.0 0 0'
end_point = '10.0 0 0'
num_points = '${Nx}'
sort_by = x
variable = 'T'
execute_on = FINAL
[]
[]
[Outputs]
exodus = true
[csv]
type = CSV
execute_on = 'FINAL'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/lid-driven/transient-lid-driven-with-energy.i)
mu = 1
rho = 1
k = .01
cp = 1
velocity_interp_method = 'rc'
advected_interp_method = 'average'
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = u
v = v
pressure = pressure
[]
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 32
ny = 32
[]
[pin]
type = ExtraNodesetGenerator
input = gen
new_boundary = 'pin'
nodes = '0'
[]
[]
[Variables]
[u]
type = INSFVVelocityVariable
[]
[v]
type = INSFVVelocityVariable
[]
[pressure]
type = INSFVPressureVariable
[]
[T]
type = INSFVEnergyVariable
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[ICs]
[T]
type = ConstantIC
variable = T
value = 1
[]
[]
[AuxVariables]
[U]
order = CONSTANT
family = MONOMIAL
fv = true
[]
[]
[AuxKernels]
[mag]
type = VectorMagnitudeAux
variable = U
x = u
y = v
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[mean_zero_pressure]
type = FVIntegralValueConstraint
variable = pressure
lambda = lambda
[]
[u_time]
type = INSFVMomentumTimeDerivative
variable = 'u'
rho = ${rho}
momentum_component = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection
variable = u
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
rho = ${rho}
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = u
mu = ${mu}
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = u
momentum_component = 'x'
pressure = pressure
[]
[v_time]
type = INSFVMomentumTimeDerivative
variable = v
rho = ${rho}
momentum_component = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = v
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = v
mu = ${mu}
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = v
momentum_component = 'y'
pressure = pressure
[]
[temp_time]
type = INSFVEnergyTimeDerivative
variable = T
rho = ${rho}
dh_dt = dh_dt
[]
[temp_conduction]
type = FVDiffusion
coeff = 'k'
variable = T
[]
[temp_advection]
type = INSFVEnergyAdvection
variable = T
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = u
boundary = 'top'
function = 'lid_function'
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = u
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = v
boundary = 'left right top bottom'
function = 0
[]
[T_hot]
type = FVDirichletBC
variable = T
boundary = 'bottom'
value = 1
[]
[T_cold]
type = FVDirichletBC
variable = T
boundary = 'top'
value = 0
[]
[]
[FunctorMaterials]
[functor_constants]
type = ADGenericFunctorMaterial
prop_names = 'cp k'
prop_values = '${cp} ${k}'
[]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
temperature = 'T'
rho = ${rho}
[]
[]
[Functions]
[lid_function]
type = ParsedFunction
expression = '4*x*(1-x)'
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
# Run for 100+ timesteps to reach steady state.
num_steps = 5
dt = .5
dtmin = .5
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
line_search = 'none'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/include/fvkernels/WCNSFVEnergyTimeDerivative.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "INSFVEnergyTimeDerivative.h"
/**
* Computes the energy time derivative for the weakly compressible formulation of the energy
* equation, using functor material properties
*/
class WCNSFVEnergyTimeDerivative : public INSFVEnergyTimeDerivative
{
public:
static InputParameters validParams();
WCNSFVEnergyTimeDerivative(const InputParameters & params);
protected:
ADReal computeQpResidual() override;
/// Functor for the time derivative of density, material property or variable
const Moose::Functor<ADReal> & _rho_dot;
/// The specific enthalpy
const Moose::Functor<ADReal> & _h;
};