Weakly Compressible Finite Volume Navier Stokes
The weakly compressible implementation of the Navier Stokes equations in finite volume is largely based of the incompressible version, except that constant densities were replaced by functor material properties. The incompressible kernels were adapted to use the functor material properties, so most kernels can be adapted directly. The only kernels that different between the two formulations are the time derivative kernels, as they include the density time derivative.
In addition, the NavierStokesFV action syntax can also be used to set up weakly-compressible simulations.
Lid Driven Cavity Heated Flow
We show below an example input file for a cavity with a moving lid, with different heat sources on both sides of the cavity.
mu = 1
rho = 'rho'
k = 1
cp = 1
l = 10
velocity_interp_method = 'rc'
advected_interp_method = 'average'
cold_temp=300
hot_temp=310
[GlobalParams<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
two_term_boundary_expansion = true
rhie_chow_user_object = 'rc'
[]
[UserObjects<<<{"href": "../../syntax/UserObjects/index.html"}>>>]
[rc]
type = INSFVRhieChowInterpolator<<<{"description": "Computes the Rhie-Chow velocity based on gathered 'a' coefficient data.", "href": "../../source/userobjects/INSFVRhieChowInterpolator.html"}>>>
u<<<{"description": "The x-component of velocity"}>>> = u
v<<<{"description": "The y-component of velocity"}>>> = v
pressure<<<{"description": "The pressure variable."}>>> = pressure
[]
[]
[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
[gen]
type = GeneratedMeshGenerator<<<{"description": "Create a line, square, or cube mesh with uniformly spaced or biased elements.", "href": "../../source/meshgenerators/GeneratedMeshGenerator.html"}>>>
dim<<<{"description": "The dimension of the mesh to be generated"}>>> = 2
xmin<<<{"description": "Lower X Coordinate of the generated mesh"}>>> = 0
xmax<<<{"description": "Upper X Coordinate of the generated mesh"}>>> = ${l}
ymin<<<{"description": "Lower Y Coordinate of the generated mesh"}>>> = 0
ymax<<<{"description": "Upper Y Coordinate of the generated mesh"}>>> = ${l}
nx<<<{"description": "Number of elements in the X direction"}>>> = 16
ny<<<{"description": "Number of elements in the Y direction"}>>> = 16
[]
[]
[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
[u]
type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e-15
[]
[v]
type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e-15
[]
[pressure]
type = INSFVPressureVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVPressureVariable.html"}>>>
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e5
[]
[T]
type = INSFVEnergyVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVEnergyVariable.html"}>>>
scaling<<<{"description": "Specifies a scaling factor to apply to this variable"}>>> = 1e-4
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${cold_temp}
[]
[]
[AuxVariables<<<{"href": "../../syntax/AuxVariables/index.html"}>>>]
[U]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
fv = true
[]
[vel_x]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[vel_y]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[viz_T]
order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
[]
[]
[AuxKernels<<<{"href": "../../syntax/AuxKernels/index.html"}>>>]
[mag]
type = VectorMagnitudeAux<<<{"description": "Creates a field representing the magnitude of three coupled variables using an Euclidean norm.", "href": "../../source/auxkernels/VectorMagnitudeAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = U
x<<<{"description": "x-component of the vector"}>>> = u
y<<<{"description": "y-component of the vector"}>>> = v
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
[]
[vel_x]
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_x
expression<<<{"description": "Parsed function expression to compute"}>>> = 'u'
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'u'
[]
[vel_y]
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = vel_y
expression<<<{"description": "Parsed function expression to compute"}>>> = 'v'
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'v'
[]
[viz_T]
type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../../source/auxkernels/ParsedAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = viz_T
expression<<<{"description": "Parsed function expression to compute"}>>> = 'T'
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
coupled_variables<<<{"description": "Vector of coupled variable names"}>>> = 'T'
[]
[]
[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
[mass_time]
type = WCNSFVMassTimeDerivative<<<{"description": "Adds the time derivative term to the weakly-compressible Navier-Stokes continuity equation.", "href": "../../source/fvkernels/WCNSFVMassTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
drho_dt<<<{"description": "The time derivative of the density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = drho_dt
[]
[mass]
type = WCNSFVMassAdvection<<<{"description": "Object for advecting mass, e.g. rho", "href": "../../source/fvkernels/WCNSFVMassAdvection.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
[]
[u_time]
type = WCNSFVMomentumTimeDerivative<<<{"description": "Adds the time derivative term to the incompressible Navier-Stokes momentum equation.", "href": "../../source/fvkernels/WCNSFVMomentumTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
drho_dt<<<{"description": "The time derivative of the density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = drho_dt
rho<<<{"description": "The density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = rho
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
[]
[u_gravity]
type = INSFVMomentumGravity<<<{"description": "Computes a body force due to gravity in Rhie-Chow based simulations.", "href": "../../source/fvkernels/INSFVMomentumGravity.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
gravity<<<{"description": "Direction of the gravity vector"}>>> = '0 -1 0'
rho<<<{"description": "The value for the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
[]
[v_time]
type = WCNSFVMomentumTimeDerivative<<<{"description": "Adds the time derivative term to the incompressible Navier-Stokes momentum equation.", "href": "../../source/fvkernels/WCNSFVMomentumTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = v
drho_dt<<<{"description": "The time derivative of the density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = drho_dt
rho<<<{"description": "The density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = rho
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = v
velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = v
mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
[]
[v_pressure]
type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = v
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
[]
[v_gravity]
type = INSFVMomentumGravity<<<{"description": "Computes a body force due to gravity in Rhie-Chow based simulations.", "href": "../../source/fvkernels/INSFVMomentumGravity.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = v
gravity<<<{"description": "Direction of the gravity vector"}>>> = '0 -1 0'
rho<<<{"description": "The value for the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
[]
[temp_conduction]
type = FVDiffusion<<<{"description": "Computes residual for diffusion operator for finite volume method.", "href": "../../source/fvkernels/FVDiffusion.html"}>>>
coeff<<<{"description": "diffusion coefficient. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'k'
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
[]
[temp_advection]
type = INSFVEnergyAdvection<<<{"description": "Advects energy, e.g. rho*cp*T. A user may still override what quantity is advected, but the default is rho*cp*T", "href": "../../source/fvkernels/INSFVEnergyAdvection.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T
velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
[]
[]
[FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
[no_slip_x]
type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = u
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left right top bottom'
function<<<{"description": "The exact solution function."}>>> = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = v
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left right top bottom'
function<<<{"description": "The exact solution function."}>>> = 0
[]
[T_hot]
type = FVDirichletBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../source/fvbcs/FVDirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = T
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = left
value<<<{"description": "value to enforce at the boundary face"}>>> = ${hot_temp}
[]
[T_cold]
type = FVDirichletBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../source/fvbcs/FVDirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = T
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = right
value<<<{"description": "value to enforce at the boundary face"}>>> = ${cold_temp}
[]
[]
[FluidProperties<<<{"href": "../../syntax/FluidProperties/index.html"}>>>]
[fp]
type = IdealGasFluidProperties<<<{"description": "Fluid properties for an ideal gas", "href": "../../source/fluidproperties/IdealGasFluidProperties.html"}>>>
[]
[]
[FunctorMaterials<<<{"href": "../../syntax/FunctorMaterials/index.html"}>>>]
[const_functor]
type = ADGenericFunctorMaterial<<<{"description": "FunctorMaterial object for declaring properties that are populated by evaluation of a Functor (a constant, variable, function or functor material property) objects.", "href": "../../source/functormaterials/GenericFunctorMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = 'cp k'
prop_values<<<{"description": "The corresponding names of the functors that are going to provide the values for the variables. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '${cp} ${k}'
[]
[rho]
type = RhoFromPTFunctorMaterial<<<{"description": "Computes the density from coupled pressure and temperature functors (variables, functions, functor material properties", "href": "../../source/functormaterials/RhoFromPTFunctorMaterial.html"}>>>
fp<<<{"description": "fluid userobject"}>>> = fp
temperature<<<{"description": "temperature functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = T
pressure<<<{"description": "pressure functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
[]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial<<<{"description": "This is the material class used to compute enthalpy for the incompressible/weakly-compressible finite-volume implementation of the Navier-Stokes equations.", "href": "../../source/functormaterials/INSFVEnthalpyFunctorMaterial.html"}>>>
temperature<<<{"description": "the temperature. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'T'
rho<<<{"description": "The value for the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
[]
[]
[Functions<<<{"href": "../../syntax/Functions/index.html"}>>>]
[lid_function]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = '4*x*(1-x)'
[]
[]
[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
steady_state_detection = true
[TimeStepper<<<{"href": "../../syntax/Executioner/TimeStepper/index.html"}>>>]
type = IterationAdaptiveDT
dt = 1e-5
optimal_iterations = 6
[]
nl_abs_tol = 1e-9
normalize_solution_diff_norm_by_dt = false
nl_max_its = 10
[]
[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
[out]
type = Exodus<<<{"description": "Object for output data in the Exodus format", "href": "../../source/outputs/Exodus.html"}>>>
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/boussinesq/transient-wcnsfv.i)Heated Channel Flow
We show below an example input file for a heated channel. The fluid progressively heats up through the channel, with a volumetric heat source.
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_v = 0.001
[Mesh<<<{"href": "../../syntax/Mesh/index.html"}>>>]
[gen]
type = GeneratedMeshGenerator<<<{"description": "Create a line, square, or cube mesh with uniformly spaced or biased elements.", "href": "../../source/meshgenerators/GeneratedMeshGenerator.html"}>>>
dim<<<{"description": "The dimension of the mesh to be generated"}>>> = 2
xmin<<<{"description": "Lower X Coordinate of the generated mesh"}>>> = 0
xmax<<<{"description": "Upper X Coordinate of the generated mesh"}>>> = ${l}
ymin<<<{"description": "Lower Y Coordinate of the generated mesh"}>>> = 0
ymax<<<{"description": "Upper Y Coordinate of the generated mesh"}>>> = 1
nx<<<{"description": "Number of elements in the X direction"}>>> = 20
ny<<<{"description": "Number of elements in the Y direction"}>>> = 10
[]
[]
[GlobalParams<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
rhie_chow_user_object = 'rc'
[]
[UserObjects<<<{"href": "../../syntax/UserObjects/index.html"}>>>]
[rc]
type = INSFVRhieChowInterpolator<<<{"description": "Computes the Rhie-Chow velocity based on gathered 'a' coefficient data.", "href": "../../source/userobjects/INSFVRhieChowInterpolator.html"}>>>
u<<<{"description": "The x-component of velocity"}>>> = vel_x
v<<<{"description": "The y-component of velocity"}>>> = vel_y
pressure<<<{"description": "The pressure variable."}>>> = pressure
[]
[]
[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
[vel_x]
type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${inlet_v}
[]
[vel_y]
type = INSFVVelocityVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVVelocityVariable.html"}>>>
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e-15
[]
[pressure]
type = INSFVPressureVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVPressureVariable.html"}>>>
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${outlet_pressure}
[]
[T_fluid]
type = INSFVEnergyVariable<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/INSFVEnergyVariable.html"}>>>
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = ${inlet_temp}
[]
[]
[AuxVariables<<<{"href": "../../syntax/AuxVariables/index.html"}>>>]
[mixing_length]
type = MooseVariableFVReal<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/MooseVariableFV.html"}>>>
[]
[power_density]
type = MooseVariableFVReal<<<{"description": "Base class for Moose variables. This should never be the terminal object type", "href": "../../source/variables/MooseVariableFV.html"}>>>
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 1e4
[]
[]
[FVKernels<<<{"href": "../../syntax/FVKernels/index.html"}>>>]
inactive<<<{"description": "If specified blocks matching these identifiers will be skipped."}>>> = 'u_turb v_turb temp_turb'
[mass_time]
type = WCNSFVMassTimeDerivative<<<{"description": "Adds the time derivative term to the weakly-compressible Navier-Stokes continuity equation.", "href": "../../source/fvkernels/WCNSFVMassTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
drho_dt<<<{"description": "The time derivative of the density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = drho_dt
[]
[mass]
type = WCNSFVMassAdvection<<<{"description": "Object for advecting mass, e.g. rho", "href": "../../source/fvkernels/WCNSFVMassAdvection.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = pressure
advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
[]
[u_time]
type = WCNSFVMomentumTimeDerivative<<<{"description": "Adds the time derivative term to the incompressible Navier-Stokes momentum equation.", "href": "../../source/fvkernels/WCNSFVMomentumTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
drho_dt<<<{"description": "The time derivative of the density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = drho_dt
rho<<<{"description": "The density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = rho
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
[]
[u_advection]
type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
[]
[u_viscosity]
type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
[]
[u_pressure]
type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
[]
[u_turb]
type = INSFVMixingLengthReynoldsStress<<<{"description": "Computes the force due to the Reynolds stress term in the incompressible Reynolds-averaged Navier-Stokes equations.", "href": "../../source/fvkernels/INSFVMixingLengthReynoldsStress.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_x
rho<<<{"description": "fluid density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
mixing_length<<<{"description": "Turbulent eddy mixing length. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'mixing_length'
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'x'
u<<<{"description": "The velocity in the x direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = vel_x
v<<<{"description": "The velocity in the y direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = vel_y
[]
[v_time]
type = WCNSFVMomentumTimeDerivative<<<{"description": "Adds the time derivative term to the incompressible Navier-Stokes momentum equation.", "href": "../../source/fvkernels/WCNSFVMomentumTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
drho_dt<<<{"description": "The time derivative of the density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = drho_dt
rho<<<{"description": "The density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = rho
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
[]
[v_advection]
type = INSFVMomentumAdvection<<<{"description": "Object for advecting momentum, e.g. rho*u", "href": "../../source/fvkernels/INSFVMomentumAdvection.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
rho<<<{"description": "Density functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
[]
[v_viscosity]
type = INSFVMomentumDiffusion<<<{"description": "Implements the Laplace form of the viscous stress in the Navier-Stokes equation.", "href": "../../source/fvkernels/INSFVMomentumDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
mu<<<{"description": "The viscosity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${mu}
[]
[v_pressure]
type = INSFVMomentumPressure<<<{"description": "Introduces the coupled pressure term into the Navier-Stokes momentum equation.", "href": "../../source/fvkernels/INSFVMomentumPressure.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
pressure<<<{"description": "The pressure. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
[]
[v_turb]
type = INSFVMixingLengthReynoldsStress<<<{"description": "Computes the force due to the Reynolds stress term in the incompressible Reynolds-averaged Navier-Stokes equations.", "href": "../../source/fvkernels/INSFVMixingLengthReynoldsStress.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = vel_y
rho<<<{"description": "fluid density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
mixing_length<<<{"description": "Turbulent eddy mixing length. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'mixing_length'
momentum_component<<<{"description": "The component of the momentum equation that this kernel applies to."}>>> = 'y'
u<<<{"description": "The velocity in the x direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = vel_x
v<<<{"description": "The velocity in the y direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = vel_y
[]
[temp_time]
type = WCNSFVEnergyTimeDerivative<<<{"description": "Adds the time derivative term to the incompressible Navier-Stokes momentum equation.", "href": "../../source/fvkernels/WCNSFVEnergyTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T_fluid
rho<<<{"description": "Density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = rho
drho_dt<<<{"description": "The time derivative of the density material property. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = drho_dt
h<<<{"description": "The specific enthalpy. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = h
dh_dt<<<{"description": "The time derivative of the specific enthalpy. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = dh_dt
[]
[temp_conduction]
type = FVDiffusion<<<{"description": "Computes residual for diffusion operator for finite volume method.", "href": "../../source/fvkernels/FVDiffusion.html"}>>>
coeff<<<{"description": "diffusion coefficient. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'k'
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T_fluid
[]
[temp_advection]
type = INSFVEnergyAdvection<<<{"description": "Advects energy, e.g. rho*cp*T. A user may still override what quantity is advected, but the default is rho*cp*T", "href": "../../source/fvkernels/INSFVEnergyAdvection.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T_fluid
velocity_interp_method<<<{"description": "The interpolation to use for the velocity. Options are 'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow."}>>> = ${velocity_interp_method}
advected_interp_method<<<{"description": "The interpolation to use for the advected quantity. Options are 'upwind', 'average', 'sou' (for second-order upwind), 'min_mod', 'vanLeer', 'quick', 'venkatakrishnan', and 'skewness-corrected' with the default being 'upwind'."}>>> = ${advected_interp_method}
[]
[heat_source]
type = FVCoupledForce<<<{"description": "Implements a source term proportional to the value of a coupled variable.", "href": "../../source/fvkernels/FVCoupledForce.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T_fluid
v<<<{"description": "The coupled functor which provides the force. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = power_density
[]
[temp_turb]
type = WCNSFVMixingLengthEnergyDiffusion<<<{"description": "Computes the turbulent diffusive flux that appears in Reynolds-averaged fluid energy conservation equations.", "href": "../../source/fvkernels/WCNSFVMixingLengthEnergyDiffusion.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = T_fluid
rho<<<{"description": "Density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = rho
cp<<<{"description": "Specific heat capacity. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = cp
mixing_length<<<{"description": "The turbulent mixing length. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'mixing_length'
schmidt_number<<<{"description": "The turbulent Schmidt number (or turbulent Prandtl number if the passive scalar is energy) that relates the turbulent scalar diffusivity to the turbulent momentum diffusivity."}>>> = 1
u<<<{"description": "The velocity in the x direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = vel_x
v<<<{"description": "The velocity in the y direction. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = vel_y
[]
[]
[FVBCs<<<{"href": "../../syntax/FVBCs/index.html"}>>>]
[no_slip_x]
type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_x
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top bottom'
function<<<{"description": "The exact solution function."}>>> = 0
[]
[no_slip_y]
type = INSFVNoSlipWallBC<<<{"description": "Implements a no slip boundary condition.", "href": "../../source/fvbcs/INSFVNoSlipWallBC.html"}>>>
variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_y
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'top bottom'
function<<<{"description": "The exact solution function."}>>> = 0
[]
# Inlet
[inlet_u]
type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_x
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
function = ${inlet_v}
[]
[inlet_v]
type = INSFVInletVelocityBC<<<{"description": "Uses the value of a functor to set a Dirichlet boundary value.", "href": "../../source/fvbcs/INSFVInletVelocityBC.html"}>>>
variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = vel_y
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
function = 0
[]
[inlet_T]
type = FVDirichletBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../source/fvbcs/FVDirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = T_fluid
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
value<<<{"description": "value to enforce at the boundary face"}>>> = ${inlet_temp}
[]
[outlet_p]
type = INSFVOutletPressureBC<<<{"description": "Defines a Dirichlet boundary condition for finite volume method.", "href": "../../source/fvbcs/INSFVOutletPressureBC.html"}>>>
variable<<<{"description": "The name of the variable that this boundary condition applies to"}>>> = pressure
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'right'
function<<<{"description": "The boundary pressure as a regular function"}>>> = ${outlet_pressure}
[]
[]
[FluidProperties<<<{"href": "../../syntax/FluidProperties/index.html"}>>>]
[fp]
type = FlibeFluidProperties<<<{"description": "Fluid properties for flibe", "href": "../../source/fluidproperties/FlibeFluidProperties.html"}>>>
[]
[]
[FunctorMaterials<<<{"href": "../../syntax/FunctorMaterials/index.html"}>>>]
[const_functor]
type = ADGenericFunctorMaterial<<<{"description": "FunctorMaterial object for declaring properties that are populated by evaluation of a Functor (a constant, variable, function or functor material property) objects.", "href": "../../source/functormaterials/GenericFunctorMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = 'cp k'
prop_values<<<{"description": "The corresponding names of the functors that are going to provide the values for the variables. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = '${cp} ${k}'
[]
[rho]
type = RhoFromPTFunctorMaterial<<<{"description": "Computes the density from coupled pressure and temperature functors (variables, functions, functor material properties", "href": "../../source/functormaterials/RhoFromPTFunctorMaterial.html"}>>>
fp<<<{"description": "fluid userobject"}>>> = fp
temperature<<<{"description": "temperature functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = T_fluid
pressure<<<{"description": "pressure functor. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = pressure
[]
[ins_fv]
type = INSFVEnthalpyFunctorMaterial<<<{"description": "This is the material class used to compute enthalpy for the incompressible/weakly-compressible finite-volume implementation of the Navier-Stokes equations.", "href": "../../source/functormaterials/INSFVEnthalpyFunctorMaterial.html"}>>>
temperature<<<{"description": "the temperature. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 'T_fluid'
rho<<<{"description": "The value for the density. A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = ${rho}
[]
[]
[AuxKernels<<<{"href": "../../syntax/AuxKernels/index.html"}>>>]
inactive<<<{"description": "If specified blocks matching these identifiers will be skipped."}>>> = 'mixing_len'
[mixing_len]
type = WallDistanceMixingLengthAux<<<{"description": "Computes the turbulent mixing length by assuming that it is proportional to the distance from the nearest wall. The mixinglength is capped at a distance proportional to inputted parameter delta.", "href": "../../source/auxkernels/WallDistanceMixingLengthAux.html"}>>>
walls<<<{"description": "Boundaries that correspond to solid walls."}>>> = 'top'
variable<<<{"description": "The name of the variable that this object applies to"}>>> = mixing_length
execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial'
delta<<<{"description": ". A functor is any of the following: a variable, a functor material property, a function, a postprocessor or a number."}>>> = 0.5
[]
[]
[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
[TimeStepper<<<{"href": "../../syntax/Executioner/TimeStepper/index.html"}>>>]
type = IterationAdaptiveDT
dt = 1e-3
optimal_iterations = 6
[]
end_time = 15
nl_abs_tol = 1e-9
nl_max_its = 50
line_search = 'none'
automatic_scaling = true
off_diagonals_in_auto_scaling = true
compute_scaling_once = false
[]
[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
[]
(modules/navier_stokes/test/tests/finite_volume/wcns/channel-flow/2d-transient.i)