- 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
FVFunctorTimeKernel
Residual contribution from time derivative of an AD functor (default is the variable this kernel is acting upon if the 'functor' parameter is not supplied) for the finite volume method.
The user may provide the functor
parameter from which to query the time derivative. If the functor
parameter is not provided, then the variable that this kernel acts on will be the functor used. The time derivative is automatically computed for nonlinear and auxiliary variables based on the time integration scheme selected. Time derivatives of Function/ADFunction
functors are computed using those objects timeDerivative
APIs. Time derivatives of functor material properties are not yet implemented. This class should be used in finite volume simulations which leverage the on-the-fly functor evaluation system, which includes incompressible and weakly compressible Navier-Stokes simulations.
When creating a new time derivative kernel, developers should consider inheriting this class as it provides the matrix/vector time tags. If not, those should be added in the validParams()
routine of the new class.
Example input syntax
In this example, the variable v
is the solution of a simple time-dependent diffusion problem. The time derivative term of the equation is added to the numerical system using a FVFunctorTimeKernel
.
[FVKernels]
[./time]
type = FVFunctorTimeKernel
variable = v
[../]
[diff]
type = FVDiffusion
variable = v
coeff = coeff
[]
[]
(test/tests/fvkernels/fv_simple_diffusion/transient.i)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
- 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
- (test/tests/fvkernels/fv_simple_diffusion/transient.i)
- (modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/lid-driven-two-phase.i)
- (modules/navier_stokes/test/tests/finite_volume/wcns/materials/functorfluidprops.i)
- (modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/channel-drift-flux-w-interface-area.i)
- (modules/navier_stokes/test/tests/finite_volume/wcns/boundary_conditions/flux_bcs_direct.i)
- (modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/channel-drift-flux-transient.i)
- (modules/navier_stokes/test/tests/finite_volume/wcns/boundary_conditions/flux_bcs_velocity.i)
- (modules/navier_stokes/test/tests/finite_volume/wcns/boundary_conditions/flux_bcs_mdot.i)
- (modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_interface_area_model/pressure_driven_growth_transient.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-std-wall-nonlinear.i)
- (modules/navier_stokes/test/tests/finite_volume/wcns/boundary_conditions/flux_bcs_reversal.i)
- (modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/rayleigh-bernard-two-phase.i)
Child Objects
(test/tests/fvkernels/fv_simple_diffusion/transient.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[v]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 7
[]
[]
[Kernels]
[]
[FVKernels]
[./time]
type = FVFunctorTimeKernel
variable = v
[../]
[diff]
type = FVDiffusion
variable = v
coeff = coeff
[]
[]
[FVBCs]
[left]
type = FVDirichletBC
variable = v
boundary = left
value = 7
[]
[right]
type = FVDirichletBC
variable = v
boundary = right
value = 42
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '.2'
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
residual_and_jacobian_together = true
num_steps = 20
dt = 0.1
[]
[Outputs]
exodus = true
[]
(test/tests/fvkernels/fv_simple_diffusion/transient.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[v]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 7
[]
[]
[Kernels]
[]
[FVKernels]
[./time]
type = FVFunctorTimeKernel
variable = v
[../]
[diff]
type = FVDiffusion
variable = v
coeff = coeff
[]
[]
[FVBCs]
[left]
type = FVDirichletBC
variable = v
boundary = left
value = 7
[]
[right]
type = FVDirichletBC
variable = v
boundary = right
value = 42
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '.2'
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
residual_and_jacobian_together = true
num_steps = 20
dt = 0.1
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/lid-driven-two-phase.i)
mu = 1.0
rho = 1.0e3
mu_d = 0.3
rho_d = 1.0
dp = 0.01
U_lid = 0.1
g = -9.81
[GlobalParams]
velocity_interp_method = 'rc'
advected_interp_method = 'upwind'
rhie_chow_user_object = 'rc'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = .1
ymin = 0
ymax = .1
nx = 10
ny = 10
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
[]
[vel_y]
type = INSFVVelocityVariable
[]
[pressure]
type = INSFVPressureVariable
[]
[phase_2]
type = INSFVScalarFieldVariable
[]
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[pin_pressure]
type = NSPressurePin
variable = pressure
pin_type = point-value
point = '0 0 0'
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
rho = 'rho_mixture'
[]
[u_time]
type = INSFVMomentumTimeDerivative
variable = vel_x
rho = 'rho_mixture'
momentum_component = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = 'rho_mixture'
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_mixture'
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[u_buoyant]
type = INSFVMomentumGravity
variable = vel_x
rho = 'rho_mixture'
momentum_component = 'x'
gravity = '0 ${g} 0'
[]
# NOTE: the friction terms for u and v are missing
[v_time]
type = INSFVMomentumTimeDerivative
variable = vel_y
rho = 'rho_mixture'
momentum_component = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = 'rho_mixture'
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_mixture'
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[v_buoyant]
type = INSFVMomentumGravity
variable = vel_y
rho = 'rho_mixture'
momentum_component = 'y'
gravity = '0 ${g} 0'
[]
[phase_2_time]
type = FVFunctorTimeKernel
variable = phase_2
[]
[phase_2_advection]
type = INSFVScalarFieldAdvection
variable = phase_2
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
[]
[phase_2_diffusion]
type = FVDiffusion
variable = phase_2
coeff = 1e-3
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = ${U_lid}
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[bottom_phase_2]
type = FVDirichletBC
variable = phase_2
boundary = 'bottom'
value = 1.0
[]
[top_phase_2]
type = FVDirichletBC
variable = phase_2
boundary = 'top'
value = 0.0
[]
[]
[AuxVariables]
[U]
order = CONSTANT
family = MONOMIAL
fv = true
[]
[drag_coefficient]
type = MooseVariableFVReal
[]
[rho_mixture_var]
type = MooseVariableFVReal
[]
[mu_mixture_var]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[mag]
type = VectorMagnitudeAux
variable = U
x = vel_x
y = vel_y
[]
[populate_cd]
type = FunctorAux
variable = drag_coefficient
functor = 'Darcy_coefficient'
[]
[populate_rho_mixture_var]
type = FunctorAux
variable = rho_mixture_var
functor = 'rho_mixture'
[]
[populate_mu_mixture_var]
type = FunctorAux
variable = mu_mixture_var
functor = 'mu_mixture'
[]
[]
[FunctorMaterials]
[populate_u_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_x'
momentum_component = 'x'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
gravity = '0 ${g} 0'
[]
[populate_v_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_y'
momentum_component = 'y'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
gravity = '0 ${g} 0'
[]
[compute_phase_1]
type = ADParsedFunctorMaterial
property_name = phase_1
functor_names = 'phase_2'
expression = '1 - phase_2'
[]
[CD]
type = NSFVDispersePhaseDragFunctorMaterial
rho = 'rho_mixture'
mu = mu_mixture
u = 'vel_x'
v = 'vel_y'
particle_diameter = ${dp}
[]
[mixing_material]
type = NSFVMixtureFunctorMaterial
phase_1_names = '${rho_d} ${mu_d}'
phase_2_names = '${rho} ${mu}'
prop_names = 'rho_mixture mu_mixture'
phase_1_fraction = 'phase_2'
[]
[]
[Postprocessors]
[average_void]
type = ElementAverageValue
variable = 'phase_2'
[]
[max_y_velocity]
type = ElementExtremeValue
variable = 'vel_y'
value_type = max
[]
[min_y_velocity]
type = ElementExtremeValue
variable = 'vel_y'
value_type = min
[]
[max_x_velocity]
type = ElementExtremeValue
variable = 'vel_x'
value_type = max
[]
[min_x_velocity]
type = ElementExtremeValue
variable = 'vel_x'
value_type = min
[]
[max_x_slip_velocity]
type = ElementExtremeFunctorValue
functor = 'vel_slip_x'
value_type = max
[]
[max_y_slip_velocity]
type = ElementExtremeFunctorValue
functor = 'vel_slip_y'
value_type = max
[]
[max_drag_coefficient]
type = ElementExtremeFunctorValue
functor = 'drag_coefficient'
value_type = max
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[TimeStepper]
type = IterationAdaptiveDT
optimal_iterations = 7
iteration_window = 2
growth_factor = 2.0
cutback_factor = 0.5
dt = 1e-3
[]
nl_max_its = 10
nl_rel_tol = 1e-03
nl_abs_tol = 1e-9
l_max_its = 5
end_time = 1e8
[]
[Outputs]
exodus = false
[CSV]
type = CSV
execute_on = 'FINAL'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/materials/functorfluidprops.i)
# Operating conditions
inlet_temp = 300
outlet_pressure = 1e5
inlet_v = 4
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 2
ymin = 0
ymax = 1
nx = 5
ny = 5
[]
[]
[Variables]
[u]
type = INSFVVelocityVariable
initial_condition = ${inlet_v}
[]
[v]
type = INSFVVelocityVariable
initial_condition = 2
[]
[pressure]
type = INSFVPressureVariable
initial_condition = ${outlet_pressure}
[]
[T]
type = INSFVEnergyVariable
initial_condition = ${inlet_temp}
[]
[]
[FVKernels]
[u_time]
type = FVFunctorTimeKernel
variable = u
[]
[v_time]
type = FVFunctorTimeKernel
variable = v
[]
[p_time]
type = FVFunctorTimeKernel
variable = pressure
[]
[T_time]
type = FVFunctorTimeKernel
variable = T
[]
[]
[FluidProperties]
[fp]
type = FlibeFluidProperties
[]
[]
[FunctorMaterials]
[fluid_props_to_mat_props]
type = GeneralFunctorFluidProps
fp = fp
pressure = 'pressure'
T_fluid = 'T'
speed = 'velocity_norm'
# For porous flow
characteristic_length = 2
porosity = 'porosity'
[]
[]
[AuxVariables]
[velocity_norm]
type = MooseVariableFVReal
[]
[porosity]
type = MooseVariableFVReal
initial_condition = 0.4
[]
[rho_var]
type = MooseVariableFVReal
[]
[drho_dp_var]
type = MooseVariableFVReal
[]
[drho_dT_var]
type = MooseVariableFVReal
[]
[rho_dot_var]
type = MooseVariableFVReal
[]
[cp_var]
type = MooseVariableFVReal
[]
[dcp_dp_var]
type = MooseVariableFVReal
[]
[dcp_dT_var]
type = MooseVariableFVReal
[]
[cp_dot_var]
type = MooseVariableFVReal
[]
[cv_var]
type = MooseVariableFVReal
[]
[mu_var]
type = MooseVariableFVReal
[]
[dmu_dp_var]
type = MooseVariableFVReal
[]
[dmu_dT_var]
type = MooseVariableFVReal
[]
[k_var]
type = MooseVariableFVReal
[]
[dk_dp_var]
type = MooseVariableFVReal
[]
[dk_dT_var]
type = MooseVariableFVReal
[]
[Pr_var]
type = MooseVariableFVReal
[]
[dPr_dp_var]
type = MooseVariableFVReal
[]
[dPr_dT_var]
type = MooseVariableFVReal
[]
[Re_var]
type = MooseVariableFVReal
[]
[dRe_dp_var]
type = MooseVariableFVReal
[]
[dRe_dT_var]
type = MooseVariableFVReal
[]
[Re_h_var]
type = MooseVariableFVReal
[]
[Re_i_var]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[speed]
type = VectorMagnitudeAux
variable = 'velocity_norm'
x = u
y = v
[]
# To output the functor material properties
[rho_out]
type = FunctorAux
functor = 'rho'
variable = 'rho_var'
execute_on = 'timestep_begin'
[]
[drho_dp_out]
type = FunctorAux
functor = 'drho/dpressure'
variable = 'drho_dp_var'
execute_on = 'timestep_begin'
[]
[drho_dT_out]
type = FunctorAux
functor = 'drho/dT_fluid'
variable = 'drho_dT_var'
execute_on = 'timestep_begin'
[]
[drho_dt_out]
type = FunctorAux
functor = 'drho_dt'
variable = 'rho_dot_var'
execute_on = 'timestep_begin'
[]
[cp_out]
type = FunctorAux
functor = 'cp'
variable = 'cp_var'
execute_on = 'timestep_begin'
[]
[dcp_dp_out]
type = FunctorAux
functor = 'dcp/dpressure'
variable = 'dcp_dp_var'
execute_on = 'timestep_begin'
[]
[dcp_dT_out]
type = FunctorAux
functor = 'dcp/dT_fluid'
variable = 'dcp_dT_var'
execute_on = 'timestep_begin'
[]
[dcp_dt_out]
type = FunctorAux
functor = 'dcp_dt'
variable = 'cp_dot_var'
execute_on = 'timestep_begin'
[]
[cv_out]
type = FunctorAux
functor = 'cv'
variable = 'cv_var'
execute_on = 'timestep_begin'
[]
[mu_out]
type = FunctorAux
functor = 'mu'
variable = 'mu_var'
execute_on = 'timestep_begin'
[]
[dmu_dp_out]
type = FunctorAux
functor = 'dmu/dpressure'
variable = 'dmu_dp_var'
execute_on = 'timestep_begin'
[]
[dmu_dT_out]
type = FunctorAux
functor = 'dmu/dT_fluid'
variable = 'dmu_dT_var'
execute_on = 'timestep_begin'
[]
[k_out]
type = FunctorAux
functor = 'k'
variable = 'k_var'
execute_on = 'timestep_begin'
[]
[dk_dp_out]
type = FunctorAux
functor = 'dk/dpressure'
variable = 'dk_dp_var'
execute_on = 'timestep_begin'
[]
[dk_dT_out]
type = FunctorAux
functor = 'dk/dT_fluid'
variable = 'dk_dT_var'
execute_on = 'timestep_begin'
[]
[Pr_out]
type = FunctorAux
functor = 'Pr'
variable = 'Pr_var'
execute_on = 'timestep_begin'
[]
[dPr_dp_out]
type = FunctorAux
functor = 'dPr/dpressure'
variable = 'dPr_dp_var'
execute_on = 'timestep_begin'
[]
[dPr_dT_out]
type = FunctorAux
functor = 'dPr/dT_fluid'
variable = 'dPr_dT_var'
execute_on = 'timestep_begin'
[]
[Re_out]
type = FunctorAux
functor = 'Re'
variable = 'Re_var'
execute_on = 'timestep_begin'
[]
[dRe_dp_out]
type = FunctorAux
functor = 'dRe/dpressure'
variable = 'dRe_dp_var'
execute_on = 'timestep_begin'
[]
[dRe_dT_out]
type = FunctorAux
functor = 'dRe/dT_fluid'
variable = 'dRe_dT_var'
execute_on = 'timestep_begin'
[]
[Re_h_out]
type = FunctorAux
functor = 'Re_h'
variable = 'Re_h_var'
execute_on = 'timestep_begin'
[]
[Re_i_out]
type = FunctorAux
functor = 'Re_i'
variable = 'Re_i_var'
execute_on = 'timestep_begin'
[]
[]
[Executioner]
type = Transient
end_time = 0.1
dt = 0.1
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/channel-drift-flux-w-interface-area.i)
mu = 10.0
rho = 100.0
mu_d = 1.0
rho_d = 1.0
l = 2
U = 1
dp = 0.01
inlet_phase_2 = 0.0
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
mass_exchange_coeff = 0.01
[GlobalParams]
rhie_chow_user_object = 'rc'
density_interp_method = 'average'
mu_interp_method = 'average'
[]
[Problem]
identify_variable_groups_in_nl = false
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = '${fparse l * 5}'
ymin = '${fparse -l / 2}'
ymax = '${fparse l / 2}'
nx = 20
ny = 5
[]
uniform_refine = 0
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0
[]
[pressure]
type = INSFVPressureVariable
[]
[phase_2]
type = INSFVScalarFieldVariable
[]
[interface_area]
type = INSFVScalarFieldVariable
[]
[]
[FVKernels]
inactive = 'u_time v_time phase_2_time interface_area_time'
[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_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_drift]
type = WCNSFV2PMomentumDriftFlux
variable = vel_x
rho_d = ${rho_d}
fd = 'rho_mixture_var'
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_mixture'
limit_interpolation = true
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[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_drift]
type = WCNSFV2PMomentumDriftFlux
variable = vel_y
rho_d = ${rho_d}
fd = 'rho_mixture_var'
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_mixture'
limit_interpolation = true
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[phase_2_time]
type = FVFunctorTimeKernel
variable = phase_2
functor = phase_2
[]
[phase_2_advection]
type = INSFVScalarFieldAdvection
variable = phase_2
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = 'upwind'
[]
[phase_2_diffusion]
type = FVDiffusion
variable = phase_2
coeff = 1.0
[]
[phase_2_src]
type = NSFVMixturePhaseInterface
variable = phase_2
phase_coupled = phase_1
alpha = ${mass_exchange_coeff}
[]
[interface_area_time]
type = FVFunctorTimeKernel
variable = interface_area
functor = interface_area
[]
[interface_area_advection]
type = INSFVScalarFieldAdvection
variable = interface_area
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = 'upwind'
[]
[interface_area_diffusion]
type = FVDiffusion
variable = interface_area
coeff = 0.1
[]
[interface_area_source_sink]
type = WCNSFV2PInterfaceAreaSourceSink
variable = interface_area
u = 'vel_x'
v = 'vel_y'
L = 1.0
rho = 'rho_mixture'
rho_d = ${rho_d}
pressure = 'pressure'
k_c = ${fparse mass_exchange_coeff * 100.0}
fd = 'phase_2'
sigma = 1e-3
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_x
functor = '${U}'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_y
functor = '0'
[]
[walls-u]
type = INSFVNoSlipWallBC
boundary = 'top bottom'
variable = vel_x
function = 0
[]
[walls-v]
type = INSFVNoSlipWallBC
boundary = 'top bottom'
variable = vel_y
function = 0
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = '0'
[]
[inlet_phase_2]
type = FVDirichletBC
boundary = 'left'
variable = phase_2
value = ${inlet_phase_2}
[]
[inlet_interface_area]
type = FVDirichletBC
boundary = 'left'
variable = interface_area
value = 0.0
[]
[]
[AuxVariables]
[drag_coefficient]
type = MooseVariableFVReal
[]
[rho_mixture_var]
type = MooseVariableFVReal
[]
[mu_mixture_var]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[populate_cd]
type = FunctorAux
variable = drag_coefficient
functor = 'Darcy_coefficient'
[]
[populate_rho_mixture_var]
type = FunctorAux
variable = rho_mixture_var
functor = 'rho_mixture'
[]
[populate_mu_mixture_var]
type = FunctorAux
variable = mu_mixture_var
functor = 'mu_mixture'
[]
[]
[FunctorMaterials]
[populate_u_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_x'
momentum_component = 'x'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[populate_v_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_y'
momentum_component = 'y'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[compute_phase_1]
type = ADParsedFunctorMaterial
property_name = phase_1
functor_names = 'phase_2'
expression = '1 - phase_2'
[]
[CD]
type = NSFVDispersePhaseDragFunctorMaterial
rho = 'rho_mixture'
mu = mu_mixture
u = 'vel_x'
v = 'vel_y'
particle_diameter = ${dp}
[]
[mixing_material]
type = NSFVMixtureFunctorMaterial
phase_2_names = '${rho} ${mu}'
phase_1_names = '${rho_d} ${mu_d}'
prop_names = 'rho_mixture mu_mixture'
phase_1_fraction = 'phase_2'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
nl_rel_tol = 1e-10
# dt = 0.1
# end_time = 1.0
# nl_max_its = 10
[]
[Debug]
show_var_residual_norms = true
[]
[Preconditioning]
[SMP]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[]
[]
[Outputs]
exodus = true
[]
[Postprocessors]
[Re]
type = ParsedPostprocessor
function = '${rho} * ${l} * ${U}'
pp_names = ''
[]
[rho_outlet]
type = SideAverageValue
boundary = 'right'
variable = 'rho_mixture_var'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/boundary_conditions/flux_bcs_direct.i)
rho = 'rho'
l = 10
inlet_area = 1
velocity_interp_method = 'rc'
advected_interp_method = 'average'
# Artificial fluid properties
# For a real case, use a GeneralFluidFunctorProperties and a viscosity rampdown
# or initialize very well!
k = 1
cp = 1000
mu = 1e2
# Operating conditions
inlet_temp = 300
outlet_pressure = 1e5
inlet_velocity = 0.001
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = 1
nx = 10
ny = 5
[]
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = u
v = v
pressure = pressure
[]
[]
[Variables]
[u]
type = INSFVVelocityVariable
initial_condition = ${inlet_velocity}
[]
[v]
type = INSFVVelocityVariable
initial_condition = 1e-15
[]
[pressure]
type = INSFVPressureVariable
initial_condition = ${outlet_pressure}
[]
[T]
type = INSFVEnergyVariable
initial_condition = ${inlet_temp}
[]
[scalar]
type = MooseVariableFVReal
initial_condition = 0.1
[]
[]
[AuxVariables]
[power_density]
type = MooseVariableFVReal
initial_condition = 1e4
[]
[]
[FVKernels]
[mass_time]
type = WCNSFVMassTimeDerivative
variable = pressure
drho_dt = drho_dt
[]
[mass]
type = WCNSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[u_time]
type = WCNSFVMomentumTimeDerivative
variable = u
drho_dt = drho_dt
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 = WCNSFVMomentumTimeDerivative
variable = v
drho_dt = drho_dt
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 = WCNSFVEnergyTimeDerivative
variable = T
rho = rho
drho_dt = drho_dt
h = h
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}
[]
[heat_source]
type = FVCoupledForce
variable = T
v = power_density
[]
# Scalar concentration equation
[scalar_time]
type = FVFunctorTimeKernel
variable = scalar
[]
[scalar_advection]
type = INSFVScalarFieldAdvection
variable = scalar
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[scalar_diffusion]
type = FVDiffusion
variable = scalar
coeff = 1.1
[]
[scalar_source]
type = FVBodyForce
variable = scalar
function = 2.1
[]
[]
[FVBCs]
# Inlet
[inlet_mass]
type = WCNSFVMassFluxBC
variable = pressure
boundary = 'left'
mdot_pp = 'inlet_mdot'
area_pp = 'surface_inlet'
vel_x = u
vel_y = v
rho = 'rho'
[]
[inlet_u]
type = WCNSFVMomentumFluxBC
variable = u
boundary = 'left'
mdot_pp = 'inlet_mdot'
area_pp = 'surface_inlet'
rho = 'rho'
momentum_component = 'x'
vel_x = u
vel_y = v
[]
[inlet_v]
type = WCNSFVMomentumFluxBC
variable = v
boundary = 'left'
mdot_pp = 0
area_pp = 'surface_inlet'
rho = 'rho'
momentum_component = 'y'
vel_x = u
vel_y = v
[]
[inlet_T]
type = WCNSFVEnergyFluxBC
variable = T
T_fluid = T
boundary = 'left'
energy_pp = 'inlet_Edot'
area_pp = 'surface_inlet'
vel_x = u
vel_y = v
rho = 'rho'
cp = cp
[]
[inlet_scalar]
type = WCNSFVScalarFluxBC
variable = scalar
boundary = 'left'
scalar_flux_pp = 'inlet_scalar_flux'
area_pp = 'surface_inlet'
vel_x = u
vel_y = v
rho = 'rho'
passive_scalar = scalar
[]
[outlet_p]
type = INSFVOutletPressureBC
variable = pressure
boundary = 'right'
function = ${outlet_pressure}
[]
# Walls
[no_slip_x]
type = INSFVNoSlipWallBC
variable = u
boundary = 'top bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = v
boundary = 'top bottom'
function = 0
[]
[]
# used for the boundary conditions in this example
[Postprocessors]
[inlet_mdot]
type = Receiver
default = ${fparse 1980 * inlet_velocity * inlet_area}
[]
[surface_inlet]
type = AreaPostprocessor
boundary = 'left'
execute_on = 'INITIAL'
[]
[inlet_Edot]
type = Receiver
default = ${fparse 1980 * inlet_velocity * 2530 * inlet_temp * inlet_area}
[]
[inlet_scalar_flux]
type = Receiver
default = ${fparse inlet_velocity * 0.2 * inlet_area}
[]
[]
[FluidProperties]
[fp]
type = SimpleFluidProperties
density0 = 1980
cp = 2530
[]
[]
[FunctorMaterials]
[const_functor]
type = ADGenericFunctorMaterial
prop_names = 'cp k'
prop_values = '${cp} ${k}'
[]
[rho]
type = RhoFromPTFunctorMaterial
fp = fp
temperature = T
pressure = pressure
[]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
temperature = 'T'
rho = ${rho}
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e-2
optimal_iterations = 6
[]
end_time = 1
nl_abs_tol = 1e-9
nl_max_its = 50
line_search = 'none'
automatic_scaling = true
[]
[Outputs]
exodus = true
execute_on = FINAL
[]
(modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/channel-drift-flux-transient.i)
mu = 1.0
rho = 10.0
mu_d = 0.1
rho_d = 1.0
l = 2
U = 1
dp = 0.01
inlet_phase_2 = 0.1
advected_interp_method = 'average'
velocity_interp_method = 'rc'
[GlobalParams]
rhie_chow_user_object = 'rc'
density_interp_method = 'average'
mu_interp_method = 'average'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = '${fparse l * 5}'
ymin = '${fparse -l / 2}'
ymax = '${fparse l / 2}'
nx = 10
ny = 4
[]
uniform_refine = 0
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0
[]
[pressure]
type = INSFVPressureVariable
[]
[phase_2]
type = INSFVScalarFieldVariable
[]
[]
[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_drift]
type = WCNSFV2PMomentumDriftFlux
variable = vel_x
rho_d = ${rho_d}
fd = 'phase_2'
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_mixture'
limit_interpolation = true
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[u_friction]
type = PINSFVMomentumFriction
Darcy_name = Darcy_coefficient_vec
is_porous_medium = false
momentum_component = x
mu = mu_mixture
rho = rho_mixture
variable = vel_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_drift]
type = WCNSFV2PMomentumDriftFlux
variable = vel_y
rho_d = ${rho_d}
fd = 'phase_2'
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_mixture'
limit_interpolation = true
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[v_friction]
type = PINSFVMomentumFriction
Darcy_name = Darcy_coefficient_vec
is_porous_medium = false
momentum_component = y
mu = mu_mixture
rho = rho_mixture
variable = vel_y
[]
[phase_2_time]
type = FVFunctorTimeKernel
variable = phase_2
functor = phase_2
[]
[phase_2_advection]
type = INSFVScalarFieldAdvection
variable = phase_2
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = 'upwind'
[]
[phase_2_src]
type = NSFVMixturePhaseInterface
variable = phase_2
phase_coupled = phase_1
alpha = 0.1
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_x
functor = '${U}'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_y
functor = '0'
[]
[walls-u]
type = INSFVNoSlipWallBC
boundary = 'top bottom'
variable = vel_x
function = 0
[]
[walls-v]
type = INSFVNoSlipWallBC
boundary = 'top bottom'
variable = vel_y
function = 0
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = '0'
[]
[inlet_phase_2]
type = FVDirichletBC
boundary = 'left'
variable = phase_2
value = ${inlet_phase_2}
[]
[]
[AuxVariables]
[drag_coefficient]
type = MooseVariableFVReal
[]
[rho_mixture_var]
type = MooseVariableFVReal
[]
[mu_mixture_var]
type = MooseVariableFVReal
[]
[vel_slip_x_var]
type = MooseVariableFVReal
[]
[vel_slip_y_var]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[populate_cd]
type = FunctorAux
variable = drag_coefficient
functor = 'Darcy_coefficient'
[]
[populate_rho_mixture_var]
type = FunctorAux
variable = rho_mixture_var
functor = 'rho_mixture'
[]
[populate_mu_mixture_var]
type = FunctorAux
variable = mu_mixture_var
functor = 'mu_mixture'
[]
[populate_vx_slip_var]
type = FunctorAux
variable = vel_slip_x_var
functor = 'vel_slip_x'
[]
[populate_vy_slip_var]
type = FunctorAux
variable = vel_slip_y_var
functor = 'vel_slip_y'
[]
[]
[FunctorMaterials]
[populate_u_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_x'
momentum_component = 'x'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[populate_v_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_y'
momentum_component = 'y'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[compute_phase_1]
type = ADParsedFunctorMaterial
property_name = phase_1
functor_names = 'phase_2'
expression = '1 - phase_2'
[]
[CD]
type = NSFVDispersePhaseDragFunctorMaterial
rho = 'rho_mixture'
mu = mu_mixture
u = 'vel_x'
v = 'vel_y'
particle_diameter = ${dp}
[]
[mixing_material]
type = NSFVMixtureFunctorMaterial
phase_2_names = '${rho} ${mu}'
phase_1_names = '${rho_d} ${mu_d}'
prop_names = 'rho_mixture mu_mixture'
phase_1_fraction = 'phase_2'
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
nl_rel_tol = 1e-10
dt = 0.1
end_time = 1.0
[]
[Preconditioning]
[SMP]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[]
[]
[Outputs]
exodus = false
[CSV]
type = CSV
execute_on = 'TIMESTEP_END'
[]
[]
[Postprocessors]
[Re]
type = ParsedPostprocessor
expression = '${rho} * ${l} * ${U}'
[]
[rho_outlet]
type = SideAverageValue
boundary = 'right'
variable = 'rho_mixture_var'
[]
[vslip_x]
type = SideExtremeValue
boundary = 'left'
variable = 'vel_slip_x_var'
[]
[vslip_y]
type = SideExtremeValue
boundary = 'left'
variable = 'vel_slip_y_var'
[]
[vslip_value]
type = ParsedPostprocessor
expression = 'sqrt(vslip_x*vslip_x + vslip_y*vslip_y)*vslip_x/abs(vslip_x)'
pp_names = 'vslip_x vslip_y'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/boundary_conditions/flux_bcs_velocity.i)
rho = 'rho'
l = 10
velocity_interp_method = 'rc'
advected_interp_method = 'average'
# Artificial fluid properties
# For a real case, use a GeneralFluidFunctorProperties and a viscosity rampdown
# or initialize very well!
k = 1
cp = 1000
mu = 1e2
# Operating conditions
inlet_temp = 300
outlet_pressure = 1e5
inlet_velocity = 0.001
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = 1
nx = 10
ny = 5
[]
[]
[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 = ${inlet_velocity}
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 1e-15
[]
[pressure]
type = INSFVPressureVariable
initial_condition = ${outlet_pressure}
[]
[T_fluid]
type = INSFVEnergyVariable
initial_condition = ${inlet_temp}
[]
[scalar]
type = MooseVariableFVReal
initial_condition = 0.1
[]
[]
[AuxVariables]
[power_density]
type = MooseVariableFVReal
initial_condition = 1e4
[]
[]
[FVKernels]
[mass_time]
type = WCNSFVMassTimeDerivative
variable = pressure
drho_dt = drho_dt
[]
[mass]
type = WCNSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[u_time]
type = WCNSFVMomentumTimeDerivative
variable = vel_x
drho_dt = drho_dt
rho = rho
momentum_component = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_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 = WCNSFVMomentumTimeDerivative
variable = vel_y
drho_dt = drho_dt
rho = rho
momentum_component = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_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
[]
[temp_time]
type = WCNSFVEnergyTimeDerivative
variable = T_fluid
rho = rho
drho_dt = drho_dt
h = h
dh_dt = dh_dt
[]
[temp_conduction]
type = FVDiffusion
coeff = 'k'
variable = T_fluid
[]
[temp_advection]
type = INSFVEnergyAdvection
variable = T_fluid
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[heat_source]
type = FVCoupledForce
variable = T_fluid
v = power_density
[]
# Scalar concentration equation
[scalar_time]
type = FVFunctorTimeKernel
variable = scalar
[]
[scalar_advection]
type = INSFVScalarFieldAdvection
variable = scalar
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[scalar_diffusion]
type = FVDiffusion
variable = scalar
coeff = 1.1
[]
[scalar_source]
type = FVBodyForce
variable = scalar
function = 2.1
[]
[]
[FVBCs]
# Inlet
[inlet_mass]
type = WCNSFVMassFluxBC
variable = pressure
boundary = 'left'
velocity_pp = 'inlet_u'
rho = 'rho'
vel_x = vel_x
vel_y = vel_y
[]
[inlet_u]
type = WCNSFVMomentumFluxBC
variable = vel_x
boundary = 'left'
velocity_pp = 'inlet_u'
rho = 'rho'
momentum_component = 'x'
vel_x = vel_x
vel_y = vel_y
[]
[inlet_v]
type = WCNSFVMomentumFluxBC
variable = vel_y
boundary = 'left'
velocity_pp = 0
rho = 'rho'
momentum_component = 'y'
vel_x = vel_x
vel_y = vel_y
[]
[inlet_T]
type = WCNSFVEnergyFluxBC
variable = T_fluid
T_fluid = T_fluid
boundary = 'left'
velocity_pp = 'inlet_u'
temperature_pp = 'inlet_T'
rho = 'rho'
cp = 'cp'
vel_x = vel_x
vel_y = vel_y
[]
[inlet_scalar]
type = WCNSFVScalarFluxBC
variable = scalar
boundary = 'left'
scalar_value_pp = 'inlet_scalar_value'
velocity_pp = 'inlet_u'
vel_x = vel_x
vel_y = vel_y
rho = rho
passive_scalar = scalar
[]
[outlet_p]
type = INSFVOutletPressureBC
variable = pressure
boundary = 'right'
function = ${outlet_pressure}
[]
# Walls
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'top bottom'
function = 0
[]
[]
# used for the boundary conditions in this example
[Postprocessors]
[inlet_u]
type = Receiver
default = ${inlet_velocity}
[]
[area_pp_left]
type = AreaPostprocessor
boundary = 'left'
execute_on = 'INITIAL'
[]
[inlet_T]
type = Receiver
default = ${inlet_temp}
[]
[inlet_scalar_value]
type = Receiver
default = 0.2
[]
[]
[FluidProperties]
[fp]
type = FlibeFluidProperties
[]
[]
[FunctorMaterials]
[const_functor]
type = ADGenericFunctorMaterial
prop_names = 'cp k'
prop_values = '${cp} ${k}'
[]
[rho]
type = RhoFromPTFunctorMaterial
fp = fp
temperature = T_fluid
pressure = pressure
[]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
temperature = 'T_fluid'
rho = ${rho}
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e-2
optimal_iterations = 6
[]
end_time = 1
nl_abs_tol = 1e-9
nl_max_its = 50
line_search = 'none'
automatic_scaling = true
[]
[Outputs]
exodus = true
execute_on = FINAL
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/boundary_conditions/flux_bcs_mdot.i)
rho = 'rho'
l = 10
inlet_area = 1
velocity_interp_method = 'rc'
advected_interp_method = 'average'
# Artificial fluid properties
# For a real case, use a GeneralFluidFunctorProperties and a viscosity rampdown
# or initialize very well!
k = 1
cp = 1000
mu = 1e2
# Operating conditions
inlet_temp = 300
outlet_pressure = 1e5
inlet_velocity = 0.001
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = 1
nx = 10
ny = 5
[]
[]
[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 = ${inlet_velocity}
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 1e-15
[]
[pressure]
type = INSFVPressureVariable
initial_condition = ${outlet_pressure}
[]
[T_fluid]
type = INSFVEnergyVariable
initial_condition = ${inlet_temp}
[]
[scalar]
type = MooseVariableFVReal
initial_condition = 0.1
[]
[]
[AuxVariables]
[power_density]
type = MooseVariableFVReal
initial_condition = 1e4
[]
[]
[FVKernels]
# Mass equation
[mass_time]
type = WCNSFVMassTimeDerivative
variable = pressure
drho_dt = drho_dt
[]
[mass]
type = WCNSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
# X component momentum equation
[u_time]
type = WCNSFVMomentumTimeDerivative
variable = vel_x
drho_dt = drho_dt
rho = rho
momentum_component = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_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
[]
# Y component momentum equation
[v_time]
type = WCNSFVMomentumTimeDerivative
variable = vel_y
drho_dt = drho_dt
rho = rho
momentum_component = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_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 equation
[temp_time]
type = WCNSFVEnergyTimeDerivative
variable = T_fluid
rho = rho
drho_dt = drho_dt
h = h
dh_dt = dh_dt
[]
[temp_conduction]
type = FVDiffusion
coeff = 'k'
variable = T_fluid
[]
[temp_advection]
type = INSFVEnergyAdvection
variable = T_fluid
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[heat_source]
type = FVCoupledForce
variable = T_fluid
v = power_density
[]
# Scalar concentration equation
[scalar_time]
type = FVFunctorTimeKernel
variable = scalar
[]
[scalar_advection]
type = INSFVScalarFieldAdvection
variable = scalar
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[scalar_diffusion]
type = FVDiffusion
variable = scalar
coeff = 1.1
[]
[scalar_source]
type = FVBodyForce
variable = scalar
function = 2.1
[]
[]
[FVBCs]
# Inlet
[inlet_mass]
type = WCNSFVMassFluxBC
variable = pressure
boundary = 'left'
mdot_pp = 'inlet_mdot'
area_pp = 'area_pp_left'
rho = 'rho'
vel_x = vel_x
vel_y = vel_y
[]
[inlet_u]
type = WCNSFVMomentumFluxBC
variable = vel_x
boundary = 'left'
mdot_pp = 'inlet_mdot'
area_pp = 'area_pp_left'
rho = 'rho'
momentum_component = 'x'
vel_x = vel_x
vel_y = vel_y
[]
[inlet_v]
type = WCNSFVMomentumFluxBC
variable = vel_y
boundary = 'left'
mdot_pp = 0
area_pp = 'area_pp_left'
rho = 'rho'
momentum_component = 'y'
vel_x = vel_x
vel_y = vel_y
[]
[inlet_T]
type = WCNSFVEnergyFluxBC
variable = T_fluid
T_fluid = T_fluid
boundary = 'left'
temperature_pp = 'inlet_T'
mdot_pp = 'inlet_mdot'
area_pp = 'area_pp_left'
rho = 'rho'
cp = 'cp'
vel_x = vel_x
vel_y = vel_y
[]
[inlet_scalar]
type = WCNSFVScalarFluxBC
variable = scalar
boundary = 'left'
scalar_value_pp = 'inlet_scalar_value'
mdot_pp = 'inlet_mdot'
area_pp = 'area_pp_left'
rho = 'rho'
vel_x = vel_x
vel_y = vel_y
passive_scalar = scalar
[]
[outlet_p]
type = INSFVOutletPressureBC
variable = pressure
boundary = 'right'
function = ${outlet_pressure}
[]
# Walls
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'top bottom'
function = 0
[]
[]
# used for the boundary conditions in this example
[Postprocessors]
[inlet_mdot]
type = Receiver
default = ${fparse 1980 * inlet_velocity * inlet_area}
[]
[area_pp_left]
type = AreaPostprocessor
boundary = 'left'
execute_on = 'INITIAL'
[]
[inlet_T]
type = Receiver
default = ${inlet_temp}
[]
[inlet_scalar_value]
type = Receiver
default = 0.2
[]
[]
[FluidProperties]
[fp]
type = FlibeFluidProperties
[]
[]
[FunctorMaterials]
[const_functor]
type = ADGenericFunctorMaterial
prop_names = 'cp k'
prop_values = '${cp} ${k}'
[]
[rho]
type = RhoFromPTFunctorMaterial
fp = fp
temperature = T_fluid
pressure = pressure
[]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
temperature = 'T_fluid'
rho = ${rho}
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e-2
optimal_iterations = 6
[]
end_time = 1
nl_abs_tol = 1e-9
nl_max_its = 50
line_search = 'none'
automatic_scaling = true
[]
[Outputs]
exodus = true
execute_on = FINAL
[]
(modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_interface_area_model/pressure_driven_growth_transient.i)
###############################################################################
# Validation test based on Hibiki and Ishii experiment [1] reported in Figure 3
# [1] Hibiki, T., & Ishii, M. (2000). One-group interfacial area transport of bubbly flows in vertical round tubes.
# International Journal of Heat and Mass Transfer, 43(15), 2711-2726.
###############################################################################
mu = 1.0
rho = 1000.0
mu_d = 1.0
rho_d = 1.0
l = ${fparse 50.8/1000.0}
U = 0.491230114
dp = 0.001
inlet_phase_2 = 0.049
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
mass_exchange_coeff = 0.0
inlet_interface_area = ${fparse 6.0*inlet_phase_2/dp}
outlet_pressure = 1e6
[GlobalParams]
rhie_chow_user_object = 'rc'
density_interp_method = 'average'
mu_interp_method = 'average'
[]
[Problem]
identify_variable_groups_in_nl = false
[]
[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 = '${fparse l * 60}'
ymin = 0
ymax = '${fparse l / 2}'
nx = 20
ny = 5
[]
uniform_refine = 0
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 0
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 0
[]
[pressure]
type = INSFVPressureVariable
[]
[phase_2]
type = INSFVScalarFieldVariable
initial_condition = ${inlet_phase_2}
[]
[interface_area]
type = INSFVScalarFieldVariable
initial_condition = ${inlet_interface_area}
[]
[]
[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_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_drift]
type = WCNSFV2PMomentumDriftFlux
variable = vel_x
rho_d = ${rho_d}
fd = 'rho_mixture_var'
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_mixture'
limit_interpolation = true
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[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_drift]
type = WCNSFV2PMomentumDriftFlux
variable = vel_y
rho_d = ${rho_d}
fd = 'rho_mixture_var'
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_mixture'
limit_interpolation = true
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[phase_2_time]
type = FVFunctorTimeKernel
variable = phase_2
functor = phase_2
[]
[phase_2_advection]
type = INSFVScalarFieldAdvection
variable = phase_2
u_slip = 'vel_x'
v_slip = 'vel_y'
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = 'upwind'
[]
[phase_2_diffusion]
type = FVDiffusion
variable = phase_2
coeff = 1.0
[]
[phase_2_src]
type = NSFVMixturePhaseInterface
variable = phase_2
phase_coupled = phase_1
alpha = ${mass_exchange_coeff}
[]
[interface_area_time]
type = FVFunctorTimeKernel
variable = interface_area
functor = interface_area
[]
[interface_area_advection]
type = INSFVScalarFieldAdvection
variable = interface_area
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = 'upwind'
[]
[interface_area_diffusion]
type = FVDiffusion
variable = interface_area
coeff = 0.1
[]
[interface_area_source_sink]
type = WCNSFV2PInterfaceAreaSourceSink
variable = interface_area
u = 'vel_x'
v = 'vel_y'
L = ${fparse l/2}
rho = 'rho_mixture'
rho_d = 'rho'
pressure = 'pressure'
k_c = '${fparse mass_exchange_coeff}'
fd = 'phase_2'
sigma = 1e-3
cutoff_fraction = 0.0
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_x
functor = '${U}'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = vel_y
functor = '0'
[]
[walls-u]
type = INSFVNoSlipWallBC
boundary = 'top'
variable = vel_x
function = 0
[]
[walls-v]
type = INSFVNoSlipWallBC
boundary = 'top'
variable = vel_y
function = 0
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = '${outlet_pressure}'
[]
[inlet_phase_2]
type = FVDirichletBC
boundary = 'left'
variable = phase_2
value = ${inlet_phase_2}
[]
[inlet_interface_area]
type = FVDirichletBC
boundary = 'left'
variable = interface_area
value = ${inlet_interface_area}
[]
[symmetry-u]
type = PINSFVSymmetryVelocityBC
boundary = 'bottom'
variable = vel_x
u = vel_x
v = vel_y
mu = 'mu_mixture'
momentum_component = 'x'
[]
[symmetry-v]
type = PINSFVSymmetryVelocityBC
boundary = 'bottom'
variable = vel_y
u = vel_x
v = vel_y
mu = 'mu_mixture'
momentum_component = 'y'
[]
[symmetry-p]
type = INSFVSymmetryPressureBC
boundary = 'bottom'
variable = pressure
[]
[symmetry-phase-2]
type = INSFVSymmetryScalarBC
boundary = 'bottom'
variable = phase_2
[]
[symmetry-interface-area]
type = INSFVSymmetryScalarBC
boundary = 'bottom'
variable = interface_area
[]
[]
[AuxVariables]
[drag_coefficient]
type = MooseVariableFVReal
[]
[rho_mixture_var]
type = MooseVariableFVReal
[]
[mu_mixture_var]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[populate_cd]
type = FunctorAux
variable = drag_coefficient
functor = 'Darcy_coefficient'
[]
[populate_rho_mixture_var]
type = FunctorAux
variable = rho_mixture_var
functor = 'rho_mixture'
[]
[populate_mu_mixture_var]
type = FunctorAux
variable = mu_mixture_var
functor = 'mu_mixture'
[]
[]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[FunctorMaterials]
[bubble_properties]
type = GeneralFunctorFluidProps
fp = 'fp'
pressure = 'pressure'
T_fluid = 300.0
speed = 1.0
characteristic_length = 1.0
porosity = 1.0
output_properties = 'rho'
outputs = 'out'
[]
[populate_u_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_x'
momentum_component = 'x'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[populate_v_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_y'
momentum_component = 'y'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[compute_phase_1]
type = ADParsedFunctorMaterial
property_name = phase_1
functor_names = 'phase_2'
expression = '1 - phase_2'
[]
[CD]
type = NSFVDispersePhaseDragFunctorMaterial
rho = 'rho_mixture'
mu = mu_mixture
u = 'vel_x'
v = 'vel_y'
particle_diameter = ${dp}
[]
[mixing_material]
type = NSFVMixtureFunctorMaterial
phase_2_names = '${rho} ${mu}'
phase_1_names = 'rho ${mu_d}'
prop_names = 'rho_mixture mu_mixture'
phase_1_fraction = 'phase_2'
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
nl_abs_tol = 1e-7
dt = 0.1
end_time = 1.0
nl_max_its = 10
line_search = 'none'
[]
[Debug]
show_var_residual_norms = true
[]
[Preconditioning]
[SMP]
type = SMP
full = true
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[]
[]
[Outputs]
[out]
type = Exodus
[]
[]
[Postprocessors]
[Re]
type = ParsedPostprocessor
expression = '${rho} * ${l} * ${U}'
pp_names = ''
[]
[rho_outlet]
type = SideAverageValue
boundary = 'right'
variable = 'rho_mixture_var'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/turbulence/lid-driven/lid-driven-turb-std-wall-nonlinear.i)
##########################################################
# Lid-driven cavity test
# Reynolds: 5,000
# Author: Dr. Mauricio Tano
# Last Update: November, 2023
# Turbulent model using:
# k-epsilon model
# No wall functions
# Newton Solve
##########################################################
### Thermophysical Properties ###
mu = 2e-5
rho = 1.0
### Operation Conditions ###
lid_velocity = 1.0
side_length = 0.1
### Initial Conditions ###
intensity = 0.01
k_init = '${fparse 1.5*(intensity * lid_velocity)^2}'
eps_init = '${fparse C_mu^0.75 * k_init^1.5 / side_length}'
### k-epsilon Closure Parameters ###
sigma_k = 1.0
sigma_eps = 1.3
C1_eps = 1.44
C2_eps = 1.92
C_mu = 0.09
### Modeling parameters ###
walls = ''
linearized_model = false
[GlobalParams]
rhie_chow_user_object = 'rc'
advected_interp_method = 'upwind'
velocity_interp_method = 'rc'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${side_length}
ymin = 0
ymax = ${side_length}
nx = 10
ny = 10
[]
[]
[Problem]
previous_nl_solution_required = true
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
initial_condition = 1e-10
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 1e-10
[]
[pressure]
type = INSFVPressureVariable
initial_condition = 0.2
[]
[TKE]
type = INSFVEnergyVariable
initial_condition = ${k_init}
[]
[TKED]
type = INSFVEnergyVariable
initial_condition = ${eps_init}
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
rho = ${rho}
[]
[mean_zero_pressure]
type = FVIntegralValueConstraint
variable = pressure
lambda = lambda
[]
[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'
mu_interp_method = average
[]
[u_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_t'
momentum_component = 'x'
complete_expansion = true
u = vel_x
v = vel_y
mu_interp_method = average
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_time]
type = INSFVMomentumTimeDerivative
variable = vel_y
rho = ${rho}
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = ${mu}
momentum_component = 'y'
mu_interp_method = average
[]
[v_viscosity_turbulent]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_t'
momentum_component = 'y'
complete_expansion = true
u = vel_x
v = vel_y
mu_interp_method = average
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[TKE_time]
type = FVFunctorTimeKernel
variable = TKE
[]
[TKE_advection]
type = INSFVTurbulentAdvection
variable = TKE
rho = ${rho}
[]
[TKE_diffusion]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = ${mu}
[]
[TKE_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKE
coeff = 'mu_t'
scaling_coef = ${sigma_k}
coeff_interp_method = average
[]
[TKE_source_sink]
type = INSFVTKESourceSink
variable = TKE
u = vel_x
v = vel_y
epsilon = TKED
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
walls = ${walls}
linearized_model = ${linearized_model}
[]
[TKED_time]
type = FVFunctorTimeKernel
variable = TKED
[]
[TKED_advection]
type = INSFVTurbulentAdvection
variable = TKED
rho = ${rho}
walls = ${walls}
[]
[TKED_diffusion]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = ${mu}
walls = ${walls}
[]
[TKED_diffusion_turbulent]
type = INSFVTurbulentDiffusion
variable = TKED
coeff = 'mu_t'
scaling_coef = ${sigma_eps}
walls = ${walls}
coeff_interp_method = average
[]
[TKED_source_sink]
type = INSFVTKEDSourceSink
variable = TKED
u = vel_x
v = vel_y
tke = TKE
rho = ${rho}
mu = ${mu}
mu_t = 'mu_t'
C1_eps = ${C1_eps}
C2_eps = ${C2_eps}
walls = ${walls}
linearized_model = ${linearized_model}
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = ${lid_velocity}
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[walls_TKE]
type = FVDirichletBC
boundary = 'left right top bottom'
variable = TKE
value = ${k_init}
[]
[walls_TKED]
type = FVDirichletBC
boundary = 'left right top bottom'
variable = TKED
value = ${eps_init}
[]
[]
[FunctorMaterials]
[mu_t_material]
type = INSFVkEpsilonViscosityFunctorMaterial
tke = TKE
epsilon = TKED
rho = ${rho}
[]
[]
[Executioner]
type = Transient
end_time = 200
dt = 0.01
steady_state_detection = true
steady_state_tolerance = 1e-3
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type -snes_linesearch_damping'
petsc_options_value = 'lu NONZERO 0.5'
nl_abs_tol = 1e-8
nl_rel_tol = 1e-8
nl_max_its = 50
line_search = none
[]
[Outputs]
exodus = true
csv = false
perf_graph = false
print_nonlinear_residuals = true
print_linear_residuals = true
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/boundary_conditions/flux_bcs_reversal.i)
rho = 'rho'
l = 10
inlet_area = 1
velocity_interp_method = 'rc'
advected_interp_method = 'average'
# Artificial fluid properties
# For a real case, use a GeneralFluidFunctorProperties and a viscosity rampdown
# or initialize very well!
k = 1
cp = 1000
mu = 1e2
# Operating conditions
inlet_temp = 300
outlet_pressure = 1e5
inlet_velocity = 0.1
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = 1
nx = 6
ny = 3
[]
[]
[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 = ${inlet_velocity}
[]
[vel_y]
type = INSFVVelocityVariable
initial_condition = 1e-15
[]
[pressure]
type = INSFVPressureVariable
initial_condition = ${outlet_pressure}
[]
[T_fluid]
type = INSFVEnergyVariable
initial_condition = ${inlet_temp}
[]
[scalar]
type = MooseVariableFVReal
initial_condition = 0.1
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[AuxVariables]
[power_density]
type = MooseVariableFVReal
initial_condition = 1e6
[]
[]
[FVKernels]
# Mass equation
[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
phi0 = 0.0
[]
# X component momentum equation
[u_time]
type = WCNSFVMomentumTimeDerivative
variable = vel_x
drho_dt = drho_dt
rho = rho
momentum_component = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_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
[]
# Y component momentum equation
[v_time]
type = WCNSFVMomentumTimeDerivative
variable = vel_y
drho_dt = drho_dt
rho = rho
momentum_component = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_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 equation
[temp_time]
type = WCNSFVEnergyTimeDerivative
variable = T_fluid
rho = rho
drho_dt = drho_dt
dh_dt = dh_dt
h = h
[]
[temp_conduction]
type = FVDiffusion
coeff = 'k'
variable = T_fluid
[]
[temp_advection]
type = INSFVEnergyAdvection
variable = T_fluid
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[heat_source]
type = FVCoupledForce
variable = T_fluid
v = power_density
[]
# Scalar concentration equation
[scalar_time]
type = FVFunctorTimeKernel
variable = scalar
[]
[scalar_advection]
type = INSFVScalarFieldAdvection
variable = scalar
velocity_interp_method = ${velocity_interp_method}
advected_interp_method = ${advected_interp_method}
[]
[scalar_diffusion]
type = FVDiffusion
variable = scalar
coeff = 1.1
[]
[scalar_source]
type = FVBodyForce
variable = scalar
function = 2.1
[]
[]
[FVBCs]
# Inlet
[inlet_mass]
type = WCNSFVMassFluxBC
variable = pressure
boundary = 'left'
mdot_pp = 'inlet_mdot'
area_pp = 'area_pp_left'
rho = 'rho'
vel_x = vel_x
vel_y = vel_y
[]
[inlet_u]
type = WCNSFVMomentumFluxBC
variable = vel_x
boundary = 'left'
mdot_pp = 'inlet_mdot'
area_pp = 'area_pp_left'
rho = 'rho'
momentum_component = 'x'
vel_x = vel_x
vel_y = vel_y
[]
[inlet_v]
type = WCNSFVMomentumFluxBC
variable = vel_y
boundary = 'left'
mdot_pp = 0
area_pp = 'area_pp_left'
rho = 'rho'
momentum_component = 'y'
vel_x = vel_x
vel_y = vel_y
[]
[inlet_T]
type = WCNSFVEnergyFluxBC
variable = T_fluid
T_fluid = T_fluid
boundary = 'left'
temperature_pp = 'inlet_T'
mdot_pp = 'inlet_mdot'
area_pp = 'area_pp_left'
rho = 'rho'
cp = 'cp'
vel_x = vel_x
vel_y = vel_y
[]
[inlet_scalar]
type = WCNSFVScalarFluxBC
variable = scalar
boundary = 'left'
scalar_value_pp = 'inlet_scalar_value'
mdot_pp = 'inlet_mdot'
area_pp = 'area_pp_left'
rho = 'rho'
vel_x = vel_x
vel_y = vel_y
passive_scalar = scalar
[]
[outlet_mass]
type = WCNSFVMassFluxBC
variable = pressure
boundary = 'right'
mdot_pp = 'outlet_mdot'
area_pp = 'area_pp_left'
rho = 'rho'
vel_x = vel_x
vel_y = vel_y
[]
[outlet_u]
type = WCNSFVMomentumFluxBC
variable = vel_x
boundary = 'right'
mdot_pp = 'outlet_mdot'
area_pp = 'area_pp_left'
rho = 'rho'
momentum_component = 'x'
vel_x = vel_x
vel_y = vel_y
[]
[outlet_v]
type = WCNSFVMomentumFluxBC
variable = vel_y
boundary = 'right'
mdot_pp = 0
area_pp = 'area_pp_left'
rho = 'rho'
momentum_component = 'y'
vel_x = vel_x
vel_y = vel_y
[]
[outlet_T]
type = WCNSFVEnergyFluxBC
variable = T_fluid
T_fluid = T_fluid
boundary = 'right'
temperature_pp = 'inlet_T'
mdot_pp = 'outlet_mdot'
area_pp = 'area_pp_left'
rho = 'rho'
cp = 'cp'
vel_x = vel_x
vel_y = vel_y
[]
[outlet_scalar]
type = WCNSFVScalarFluxBC
variable = scalar
boundary = 'right'
scalar_value_pp = 'inlet_scalar_value'
mdot_pp = 'outlet_mdot'
area_pp = 'area_pp_left'
rho = 'rho'
vel_x = vel_x
vel_y = vel_y
passive_scalar = scalar
[]
# Walls
[no_slip_x]
type = INSFVNaturalFreeSlipBC
variable = vel_x
momentum_component = x
boundary = 'top bottom'
[]
[no_slip_y]
type = INSFVNaturalFreeSlipBC
variable = vel_y
momentum_component = y
boundary = 'top bottom'
[]
[]
# used for the boundary conditions in this example
[Postprocessors]
[inlet_mdot]
type = Receiver
default = ${fparse 1980 * inlet_velocity * inlet_area}
#outputs = none
[]
[outlet_mdot]
type = Receiver
default = ${fparse -1980 * inlet_velocity * inlet_area}
outputs = none
[]
[area_pp_left]
type = AreaPostprocessor
boundary = 'left'
execute_on = 'INITIAL'
outputs = none
[]
[inlet_T]
type = Receiver
default = ${inlet_temp}
outputs = none
[]
[inlet_scalar_value]
type = Receiver
default = 0.2
outputs = none
[]
[left_mdot]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
advected_quantity = rho
boundary = left
#advected_interp_method = ${advected_interp_method}
[]
[right_mdot]
type = VolumetricFlowRate
vel_x = vel_x
vel_y = vel_y
advected_quantity = rho
boundary = right
advected_interp_method = upwind #${advected_interp_method}
[]
[]
[FluidProperties]
[fp]
type = FlibeFluidProperties
[]
[]
[FunctorMaterials]
[const_functor]
type = ADGenericFunctorMaterial
prop_names = 'cp k rho'
prop_values = '${cp} ${k} 1980'
[]
#[rho]
# type = RhoFromPTFunctorMaterial
# fp = fp
# temperature = T_fluid
# pressure = pressure
#[]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial
temperature = 'T_fluid'
rho = ${rho}
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[TimeStepper]
type = IterationAdaptiveDT
dt = 1e-1
optimal_iterations = 6
growth_factor = 4
[]
end_time = 500000
nl_abs_tol = 1e-7
nl_max_its = 50
line_search = 'none'
automatic_scaling = true
[]
[Outputs]
exodus = true
execute_on = FINAL
[]
(modules/navier_stokes/test/tests/finite_volume/two_phase/mixture_model/rayleigh-bernard-two-phase.i)
mu = 1.0
rho = 1e3
mu_d = 0.3
rho_d = 1.0
dp = 0.01
U_lid = 0.0
g = -9.81
[GlobalParams]
velocity_interp_method = 'rc'
advected_interp_method = 'upwind'
rhie_chow_user_object = 'rc'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = .1
ymin = 0
ymax = .1
nx = 11
ny = 11
[]
[]
[Variables]
[vel_x]
type = INSFVVelocityVariable
[]
[vel_y]
type = INSFVVelocityVariable
[]
[pressure]
type = INSFVPressureVariable
[]
[phase_2]
type = INSFVScalarFieldVariable
[]
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = vel_x
v = vel_y
pressure = pressure
[]
[pin_pressure]
type = NSPressurePin
variable = pressure
pin_type = point-value
point = '0 0 0'
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
rho = 'rho_mixture'
[]
[u_time]
type = INSFVMomentumTimeDerivative
variable = vel_x
rho = 'rho_mixture'
momentum_component = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection
variable = vel_x
rho = 'rho_mixture'
momentum_component = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion
variable = vel_x
mu = 'mu_mixture'
momentum_component = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure
variable = vel_x
momentum_component = 'x'
pressure = pressure
[]
[u_buoyant]
type = INSFVMomentumGravity
variable = vel_x
rho = 'rho_mixture'
momentum_component = 'x'
gravity = '0 ${g} 0'
[]
[v_time]
type = INSFVMomentumTimeDerivative
variable = vel_y
rho = 'rho_mixture'
momentum_component = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = vel_y
rho = 'rho_mixture'
momentum_component = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion
variable = vel_y
mu = 'mu_mixture'
momentum_component = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure
variable = vel_y
momentum_component = 'y'
pressure = pressure
[]
[v_buoyant]
type = INSFVMomentumGravity
variable = vel_y
rho = 'rho_mixture'
momentum_component = 'y'
gravity = '0 ${g} 0'
[]
[phase_2_time]
type = FVFunctorTimeKernel
variable = phase_2
[]
[phase_2_advection]
type = INSFVScalarFieldAdvection
variable = phase_2
u_slip = 'vel_slip_x'
v_slip = 'vel_slip_y'
[]
[phase_2_diffusion]
type = FVDiffusion
variable = phase_2
coeff = 1e-3
[]
[]
[FVBCs]
[top_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'top'
function = ${U_lid}
[]
[no_slip_x]
type = INSFVNoSlipWallBC
variable = vel_x
boundary = 'left right bottom'
function = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC
variable = vel_y
boundary = 'left right top bottom'
function = 0
[]
[bottom_phase_2]
type = FVDirichletBC
variable = phase_2
boundary = 'bottom'
value = 1.0
[]
[top_phase_2]
type = FVDirichletBC
variable = phase_2
boundary = 'top'
value = 0.0
[]
[]
[AuxVariables]
[U]
order = CONSTANT
family = MONOMIAL
fv = true
[]
[drag_coefficient]
type = MooseVariableFVReal
[]
[rho_mixture_var]
type = MooseVariableFVReal
[]
[mu_mixture_var]
type = MooseVariableFVReal
[]
[phase_1]
type = MooseVariableFVReal
[]
[]
[AuxKernels]
[mag]
type = VectorMagnitudeAux
variable = U
x = vel_x
y = vel_y
[]
[populate_cd]
type = FunctorAux
variable = drag_coefficient
functor = 'Darcy_coefficient'
[]
[populate_rho_mixture_var]
type = FunctorAux
variable = rho_mixture_var
functor = 'rho_mixture'
[]
[populate_mu_mixture_var]
type = FunctorAux
variable = mu_mixture_var
functor = 'mu_mixture'
[]
[compute_phase_1]
type = ParsedAux
variable = phase_1
coupled_variables = 'phase_2'
expression = '1 - phase_2'
[]
[]
[FunctorMaterials]
[CD]
type = NSFVDispersePhaseDragFunctorMaterial
rho = 'rho_mixture'
mu = mu_mixture
u = 'vel_x'
v = 'vel_y'
particle_diameter = ${dp}
[]
[mixing_material]
type = NSFVMixtureFunctorMaterial
phase_1_names = '${rho_d} ${mu_d}'
phase_2_names = '${rho} ${mu}'
prop_names = 'rho_mixture mu_mixture'
phase_1_fraction = 'phase_2'
[]
[populate_u_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_x'
momentum_component = 'x'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[populate_v_slip]
type = WCNSFV2PSlipVelocityFunctorMaterial
slip_velocity_name = 'vel_slip_y'
momentum_component = 'y'
u = 'vel_x'
v = 'vel_y'
rho = ${rho}
mu = 'mu_mixture'
rho_d = ${rho_d}
particle_diameter = ${dp}
linear_coef_name = 'Darcy_coefficient'
[]
[]
[Postprocessors]
[average_void]
type = ElementAverageValue
variable = 'phase_2'
[]
[max_y_velocity]
type = ElementExtremeValue
variable = 'vel_y'
value_type = max
[]
[min_y_velocity]
type = ElementExtremeValue
variable = 'vel_y'
value_type = min
[]
[max_x_velocity]
type = ElementExtremeValue
variable = 'vel_x'
value_type = max
[]
[min_x_velocity]
type = ElementExtremeValue
variable = 'vel_x'
value_type = min
[]
[max_x_slip_velocity]
type = ElementExtremeFunctorValue
functor = 'vel_slip_x'
value_type = max
[]
[max_y_slip_velocity]
type = ElementExtremeFunctorValue
functor = 'vel_slip_y'
value_type = max
[]
[max_drag_coefficient]
type = ElementExtremeFunctorValue
functor = 'drag_coefficient'
value_type = max
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[TimeStepper]
type = IterationAdaptiveDT
optimal_iterations = 10
iteration_window = 2
growth_factor = 2
cutback_factor = 0.5
dt = 1e-3
[]
nl_max_its = 20
nl_rel_tol = 1e-03
nl_abs_tol = 1e-9
l_max_its = 5
end_time = 1e8
[]
[Outputs]
exodus = false
[CSV]
type = CSV
execute_on = 'FINAL'
[]
[]
(modules/navier_stokes/include/fvkernels/WCNSFVMassTimeDerivative.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 "FVFunctorTimeKernel.h"
/**
* Computes the mass time derivative for the weakly compressible formulation of the mass
* equation, using functor material properties
*/
class WCNSFVMassTimeDerivative : public FVFunctorTimeKernel
{
public:
static InputParameters validParams();
WCNSFVMassTimeDerivative(const InputParameters & params);
protected:
ADReal computeQpResidual() override;
const Moose::Functor<ADReal> & _rho_dot;
};
(modules/navier_stokes/include/fvkernels/INSFVEnergyTimeDerivative.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 "FVFunctorTimeKernel.h"
class INSFVEnergyTimeDerivative : public FVFunctorTimeKernel
{
public:
static InputParameters validParams();
INSFVEnergyTimeDerivative(const InputParameters & params);
protected:
ADReal computeQpResidual() override;
/// the density
const Moose::Functor<ADReal> & _rho;
/// The time derivative of the specific enthalpy
const Moose::Functor<ADReal> & _h_dot;
};
(modules/navier_stokes/include/fvkernels/INSFVTimeKernel.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 "FVFunctorTimeKernel.h"
#include "INSFVMomentumResidualObject.h"
/**
* All navier-stokes momentum time derivative terms should inherit from this class
*/
class INSFVTimeKernel : public FVFunctorTimeKernel, public INSFVMomentumResidualObject
{
public:
static InputParameters validParams();
INSFVTimeKernel(const InputParameters & params);
using INSFVMomentumResidualObject::gatherRCData;
void gatherRCData(const FaceInfo &) override final {}
virtual ~INSFVTimeKernel() = default;
void computeResidual() override final {}
void computeJacobian() override final {}
using FVFunctorTimeKernel::computeOffDiagJacobian;
void computeOffDiagJacobian() override final {}
void computeResidualAndJacobian() override final {}
protected:
ADReal computeQpResidual() override final
{
mooseError("INSFVTimeKernels must implement gatherRCData and not computeQpResidual");
}
/**
* Process into either the system residual or Jacobian
*/
void addResidualAndJacobian(const ADReal & residual, dof_id_type dof);
/// Whether to contribute to RC coefficients
const bool _contribute_to_rc_coeffs;
private:
using FVFunctorTimeKernel::_current_elem;
};
(modules/navier_stokes/include/fvkernels/PINSFVEnergyTimeDerivative.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 "FVFunctorTimeKernel.h"
class PINSFVEnergyTimeDerivative : public FVFunctorTimeKernel
{
public:
static InputParameters validParams();
PINSFVEnergyTimeDerivative(const InputParameters & params);
protected:
ADReal computeQpResidual() override;
/// the density
const Moose::Functor<ADReal> & _rho;
/// the time derivative of the density
const Moose::Functor<ADReal> * const _rho_dot;
/// the specific heat or isobaric heat capacity
const Moose::Functor<ADReal> * const _cp;
/// the specific enthalpy
const Moose::Functor<ADReal> * const _h;
/// the time derivative of the specific enthalpy
const Moose::Functor<ADReal> * const _h_dot;
/// the porosity
const Moose::Functor<ADReal> & _eps;
/// whether this kernel is being used for a solid or a fluid temperature
const bool _is_solid;
/// scales the value of the kernel, used for faster steady state during pseudo transient
const Real _scaling;
/// whether a zero scaling factor has been specifed
const bool _zero_scaling;
};