- boundaryThe list of boundary IDs from the mesh where this boundary condition applies
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:The list of boundary IDs from the mesh where this boundary condition applies
- functionThe exact solution function.
C++ Type:FunctionName
Controllable:No
Description:The exact solution function.
- variableThe name of the variable that this boundary condition applies to
C++ Type:NonlinearVariableName
Controllable:No
Description:The name of the variable that this boundary condition applies to
FVFunctionDirichletBC
Imposes the essential boundary condition , where is a (possibly) time and space-dependent MOOSE Function.
Overview
FVFunctionDirichletBC
is very similar to FVDirichletBC except the parameter value
is replaced by the function
parameter, where the latter is a FunctionName
or alternatively a direct input of a parsed function. FVFunctionDirichletBC
is generally useful; it's critical for implementing MMS studies.
Example Input File Syntax
[FVBCs]
[boundary]
type = FVFunctionDirichletBC
boundary = 'left right'
function = 'exact'
variable = v
[]
[]
(test/tests/fvkernels/mms/advection-diffusion.i)Input Parameters
- displacementsThe displacements
C++ Type:std::vector<VariableName>
Controllable:No
Description:The displacements
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:Yes
Description:Set the enabled status of the MooseObject.
- implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Controllable:No
Description:Determines whether this object is calculated using an implicit or explicit form
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
Description:The extra tags for the matrices this Kernel should fill
- extra_vector_tagsThe extra tags for the vectors this Kernel should fill
C++ Type:std::vector<TagName>
Controllable:No
Description:The extra tags for the vectors this Kernel should fill
- matrix_tagssystemThe tag for the matrices this Kernel should fill
Default:system
C++ Type:MultiMooseEnum
Options:nontime, system
Controllable:No
Description:The tag for the matrices this Kernel should fill
- vector_tagsnontimeThe tag for the vectors this Kernel should fill
Default:nontime
C++ Type:MultiMooseEnum
Options:nontime, time
Controllable:No
Description:The tag for the vectors this Kernel should fill
Tagging Parameters
Input Files
- (test/tests/fvkernels/mms/advective-outflow/advection-diffusion.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/varying-eps-basic-kt-primitive.i)
- (test/tests/fvkernels/mms/mat-advection-diffusion.i)
- (test/tests/postprocessors/interface_value/interface_fv_variable_value_postprocessor.i)
- (test/tests/fvkernels/mms/grad-reconstruction/mat-rz.i)
- (test/tests/fvkernels/mms/non-orthogonal/advection-diffusion-reaction.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/basic-primitive-pcnsfv-kt.i)
- (test/tests/fvkernels/mms/skewness-correction/adv-diff-react/skewed.i)
- (test/tests/fvkernels/mms/diffusion.i)
- (test/tests/fvkernels/mms/cylindrical/advection-diffusion-reaction.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/mms/channel-flow/cylindrical/2d-average-with-temp.i)
- (test/tests/fvkernels/mms/grad-reconstruction/mat-cartesian.i)
- (test/tests/postprocessors/element_integral_var_pps/pps_old_value_fv.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/varying-eps-basic-kt-mixed.i)
- (test/tests/fvkernels/mms/non-orthogonal/extended-adr.i)
- (modules/navier_stokes/test/tests/finite_volume/pins/mms/1d-rc-no-diffusion-strong-bc.i)
- (test/tests/fvkernels/mms/advective-outflow/advection.i)
- (test/tests/fvkernels/mms/advective-outflow/advection-outflow.i)
- (test/tests/fvkernels/mms/skewness-correction/two_term_extrapol/advection-outflow.i)
- (test/tests/fvkernels/mms/skewness-correction/diffusion/skewed.i)
- (test/tests/fvkernels/mms/grad-reconstruction/rz.i)
- (test/tests/fvkernels/mms/grad-reconstruction/cartesian.i)
- (test/tests/fvkernels/mms/advective-outflow/kt-limited-advection.i)
- (test/tests/fvkernels/mms/advection-diffusion.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/mms/1d/hllc-mms.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/basic-conserved-pcnsfv-kt.i)
- (modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/pwcnsfv.i)
- (test/tests/postprocessors/interface_diffusive_flux/interface_diffusive_flux_fv.i)
- (test/tests/fvkernels/mms/cylindrical/diffusion.i)
- (test/tests/fvkernels/mms/mass-mom-mat-advection-diffusion/input.i)
- (modules/navier_stokes/test/tests/finite_volume/ins/mms/channel-flow/2d-average-with-temp.i)
Child Objects
(test/tests/fvkernels/mms/advection-diffusion.i)
diff=1.1
a=1.1
[GlobalParams]
advected_interp_method = 'average'
[]
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 1
xmin = -0.6
xmax = 0.6
nx = 64
[../]
[]
[Variables]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
[../]
[]
[FVKernels]
[./advection]
type = FVAdvection
variable = v
velocity = '${a} 0 0'
[../]
[./diffusion]
type = FVDiffusion
variable = v
coeff = coeff
[../]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[boundary]
type = FVFunctionDirichletBC
boundary = 'left right'
function = 'exact'
variable = v
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '${diff}'
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = '3*x^2 + 2*x + 1'
[]
[forcing]
type = ParsedFunction
value = '-${diff}*6 + ${a} * (6*x + 2)'
# value = '-${diff}*6'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/advective-outflow/advection-diffusion.i)
diff=1
a=1
[GlobalParams]
advected_interp_method = 'average'
[]
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 1
xmin = -1
xmax = 0
nx = 2
[../]
[]
[Variables]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
[../]
[]
[FVKernels]
[./advection]
type = FVAdvection
variable = v
velocity = '${a} 0 0'
force_boundary_execution = true
[../]
[./diffusion]
type = FVDiffusion
variable = v
coeff = coeff
[../]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[left]
type = FVFunctionDirichletBC
boundary = 'left'
function = 'exact'
variable = v
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '${diff}'
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'cos(x)'
[]
[forcing]
type = ParsedFunction
value = 'cos(x) - sin(x)'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/varying-eps-basic-kt-primitive.i)
[GlobalParams]
fp = fp
limiter = 'central_difference'
two_term_boundary_expansion = true
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = .1
xmax = .6
nx = 2
[]
[]
[Modules]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
[]
[sup_vel_x]
type = MooseVariableFVReal
[]
[T_fluid]
type = MooseVariableFVReal
[]
[]
[ICs]
[pressure]
type = FunctionIC
variable = pressure
function = 'exact_p'
[]
[sup_vel_x]
type = FunctionIC
variable = sup_vel_x
function = 'exact_sup_vel_x'
[]
[T_fluid]
type = FunctionIC
variable = T_fluid
function = 'exact_T'
[]
[]
[FVKernels]
[mass_advection]
type = PCNSFVKT
variable = pressure
eqn = "mass"
[]
[mass_fn]
type = FVBodyForce
variable = pressure
function = 'forcing_rho'
[]
[momentum_x_advection]
type = PCNSFVKT
variable = sup_vel_x
momentum_component = x
eqn = "momentum"
[]
[eps_grad]
type = PNSFVPGradEpsilon
variable = sup_vel_x
momentum_component = 'x'
epsilon_function = 'eps'
[]
[momentum_fn]
type = FVBodyForce
variable = sup_vel_x
function = 'forcing_rho_ud'
[]
[fluid_energy_advection]
type = PCNSFVKT
variable = T_fluid
eqn = "energy"
[]
[energy_fn]
type = FVBodyForce
variable = T_fluid
function = 'forcing_rho_et'
[]
[]
[FVBCs]
[mass_left]
variable = pressure
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'mass'
[]
[momentum_left]
variable = sup_vel_x
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'momentum'
momentum_component = 'x'
[]
[energy_left]
variable = T_fluid
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'energy'
[]
[mass_right]
variable = pressure
type = PCNSFVStrongBC
boundary = right
eqn = 'mass'
pressure = 'exact_p'
[]
[momentum_right]
variable = sup_vel_x
type = PCNSFVStrongBC
boundary = right
eqn = 'momentum'
momentum_component = 'x'
pressure = 'exact_p'
[]
[energy_right]
variable = T_fluid
type = PCNSFVStrongBC
boundary = right
eqn = 'energy'
pressure = 'exact_p'
[]
# help gradient reconstruction
[pressure_right]
type = FVFunctionDirichletBC
variable = pressure
function = exact_p
boundary = 'right'
[]
[sup_vel_x_left]
type = FVFunctionDirichletBC
variable = sup_vel_x
function = exact_sup_vel_x
boundary = 'left'
[]
[T_fluid_left]
type = FVFunctionDirichletBC
variable = T_fluid
function = exact_T
boundary = 'left'
[]
[]
[Materials]
[var_mat]
type = PorousPrimitiveVarMaterial
pressure = pressure
superficial_vel_x = sup_vel_x
T_fluid = T_fluid
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[]
[Functions]
[exact_rho]
type = ParsedFunction
value = '3.48788261470924*cos(x)'
[]
[forcing_rho]
type = ParsedFunction
value = '-3.83667087618017*sin(1.1*x)*cos(1.3*x) - 4.53424739912202*sin(1.3*x)*cos(1.1*x)'
[]
[exact_rho_ud]
type = ParsedFunction
value = '3.48788261470924*cos(1.1*x)*cos(1.3*x)'
[]
[forcing_rho_ud]
type = ParsedFunction
value = '(-(10.6975765229419*cos(1.5*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.5*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 16.0463647844128*sin(1.5*x)/cos(x))*cos(x))*cos(1.3*x) + 3.48788261470924*sin(x)*cos(1.1*x)^2*cos(1.3*x)/cos(x)^2 - 7.67334175236034*sin(1.1*x)*cos(1.1*x)*cos(1.3*x)/cos(x) - 4.53424739912202*sin(1.3*x)*cos(1.1*x)^2/cos(x)'
[]
[exact_rho_et]
type = ParsedFunction
value = '26.7439413073546*cos(1.5*x)'
[]
[forcing_rho_et]
type = ParsedFunction
value = '1.0*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(x)*cos(1.1*x)*cos(1.3*x)/cos(x)^2 - 1.1*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(1.1*x)*cos(1.3*x)/cos(x) - 1.3*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(1.3*x)*cos(1.1*x)/cos(x) + 1.0*(-(10.6975765229419*cos(1.5*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.5*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 16.0463647844128*sin(1.5*x)/cos(x))*cos(x) - 40.1159119610319*sin(1.5*x))*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[exact_T]
type = ParsedFunction
value = '0.0106975765229418*cos(1.5*x)/cos(x) - 0.000697576522941848*cos(1.1*x)^2/cos(x)^2'
[]
[exact_eps_p]
type = ParsedFunction
value = '3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)*cos(1.3*x)'
[]
[exact_p]
type = ParsedFunction
value = '3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_sup_vel_x]
type = ParsedFunction
value = '1.0*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[eps]
type = ParsedFunction
value = 'cos(1.3*x)'
[]
[exact_superficial_velocity]
type = ParsedVectorFunction
value_x = '1.0*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[]
[Executioner]
solve_type = NEWTON
type = Transient
num_steps = 1
dtmin = 1
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_max_its = 50
line_search = bt
nl_rel_tol = 1e-12
nl_abs_tol = 1e-12
[]
[Outputs]
exodus = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2pressure]
type = ElementL2Error
variable = pressure
function = exact_p
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2sup_vel_x]
variable = sup_vel_x
function = exact_sup_vel_x
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2T_fluid]
variable = T_fluid
function = exact_T
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/mat-advection-diffusion.i)
diff=1.1
a=1.1
[GlobalParams]
advected_interp_method = 'average'
[]
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 1
xmin = -0.6
xmax = 0.6
nx = 64
[../]
[]
[Variables]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
[../]
[]
[FVKernels]
[./advection]
type = FVMatAdvection
variable = v
vel = 'fv_velocity'
[../]
[./diffusion]
type = FVDiffusion
variable = v
coeff = coeff
[../]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[boundary]
type = FVFunctionDirichletBC
boundary = 'left right'
function = 'exact'
variable = v
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '${diff}'
[]
[adv_material]
type = ADCoupledVelocityMaterial
vel_x = '${a}'
rho = 'v'
velocity = 'fv_velocity'
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = '3*x^2 + 2*x + 1'
[]
[forcing]
type = ParsedFunction
value = '-${diff}*6 + ${a} * (6*x + 2)'
# value = '-${diff}*6'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/postprocessors/interface_value/interface_fv_variable_value_postprocessor.i)
postprocessor_type = InterfaceAverageVariableValuePostprocessor
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 6
xmax = 3
ny = 9
ymax = 3
elem_type = QUAD4
[]
[./subdomain_id]
input = gen
type = SubdomainBoundingBoxGenerator
bottom_left = '0 0 0'
top_right = '2 1 0'
block_id = 1
[../]
[./interface]
type = SideSetsBetweenSubdomainsGenerator
input = subdomain_id
primary_block = '0'
paired_block = '1'
new_boundary = 'interface'
[../]
[]
[Functions]
[./fn_exact]
type = ParsedFunction
value = 'x*x+y*y'
[../]
[./ffn]
type = ParsedFunction
value = -4
[../]
[]
[Variables]
[./u]
family = MONOMIAL
order = CONSTANT
fv = true
[../]
[]
[FVKernels]
[./diff]
type = FVDiffusion
variable = u
coeff = 1
[../]
[./ffn]
type = FVBodyForce
variable = u
function = ffn
[../]
[]
[FVBCs]
[./all]
type = FVFunctionDirichletBC
variable = u
boundary = '0 1 2 3'
function = fn_exact
[../]
[]
[Materials]
[./stateful1]
type = GenericConstantMaterial
block = 0
prop_names = 'diffusivity'
prop_values = 10
[../]
[./stateful2]
type = GenericConstantMaterial
block = 1
prop_names = 'diffusivity'
prop_values = 4
[../]
[]
[AuxKernels]
[./diffusivity_1]
type = MaterialRealAux
property = diffusivity
variable = diffusivity_1
[]
[./diffusivity_2]
type = MaterialRealAux
property = diffusivity
variable = diffusivity_2
[]
[]
[AuxVariables]
[./diffusivity_1]
family = MONOMIAL
order = CONSTANT
[]
[./diffusivity_2]
family = MONOMIAL
order = CONSTANT
[]
[]
[Postprocessors]
[./diffusivity_average]
type = ${postprocessor_type}
interface_value_type = average
variable = diffusivity_1
neighbor_variable = diffusivity_2
execute_on = TIMESTEP_END
boundary = 'interface'
[../]
[./diffusivity_jump_primary_secondary]
type = ${postprocessor_type}
interface_value_type = jump_primary_minus_secondary
variable = diffusivity_1
neighbor_variable = diffusivity_2
execute_on = TIMESTEP_END
boundary = 'interface'
[../]
[./diffusivity_jump_secondary_primary]
type = ${postprocessor_type}
interface_value_type = jump_secondary_minus_primary
variable = diffusivity_1
neighbor_variable = diffusivity_2
execute_on = TIMESTEP_END
boundary = 'interface'
[../]
[./diffusivity_jump_abs]
type = ${postprocessor_type}
interface_value_type = jump_abs
variable = diffusivity_1
neighbor_variable = diffusivity_2
execute_on = TIMESTEP_END
boundary = 'interface'
[../]
[./diffusivity_primary]
type = ${postprocessor_type}
interface_value_type = primary
variable = diffusivity_1
neighbor_variable = diffusivity_2
execute_on = TIMESTEP_END
boundary = 'interface'
[../]
[./diffusivity_secondary]
type = ${postprocessor_type}
interface_value_type = secondary
variable = diffusivity_1
neighbor_variable = diffusivity_2
execute_on = TIMESTEP_END
boundary = 'interface'
[../]
[./diffusivity_single_variable]
type = ${postprocessor_type}
interface_value_type = primary
variable = diffusivity_1
execute_on = TIMESTEP_END
boundary = 'interface'
[../]
[]
[Executioner]
type = Steady
solve_type = NEWTON
[]
[Outputs]
file_base = ${raw ${postprocessor_type} _fv}
exodus = true
[]
(test/tests/fvkernels/mms/grad-reconstruction/mat-rz.i)
a=1.1
diff=1.1
[Mesh]
[gen_mesh]
type = GeneratedMeshGenerator
dim = 2
xmin = 2
xmax = 3
ymin = 0
ymax = 1
nx = 2
ny = 2
[]
[]
[Problem]
coord_type = 'RZ'
[]
[Variables]
[v]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 1
[]
[]
[FVKernels]
[advection]
type = FVElementalAdvection
variable = v
velocity = '${a} ${a} 0'
advected_quantity = 'mat_u'
grad_advected_quantity = 'mat_grad_u'
[]
[reaction]
type = FVReaction
variable = v
[]
[diff_v]
type = FVDiffusion
variable = v
coeff = ${diff}
[]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[diri]
type = FVFunctionDirichletBC
boundary = 'left right top bottom'
function = 'exact'
variable = v
[]
[]
[Materials]
[mat]
type = ADCoupledGradientMaterial
mat_prop = 'mat_u'
grad_mat_prop = 'mat_grad_u'
u = v
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'sin(x)*cos(y)'
[]
[forcing]
type = ParsedFunction
value = '-a*sin(x)*sin(y) + diff*sin(x)*cos(y) + sin(x)*cos(y) + (x*a*cos(x)*cos(y) + a*sin(x)*cos(y))/x - (-x*diff*sin(x)*cos(y) + diff*cos(x)*cos(y))/x'
vars = 'a diff'
vals = '${a} ${diff}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_shift_type -sub_pc_type'
petsc_options_value = 'asm NONZERO lu'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/non-orthogonal/advection-diffusion-reaction.i)
a=1.1
diff=1.1
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 2
xmin = 2
xmax = 3
ymin = 0
ymax = 1
nx = 2
ny = 2
elem_type = TRI3
[../]
[]
[Variables]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 1
[../]
[]
[FVKernels]
[./advection]
type = FVAdvection
variable = v
velocity = '${a} ${fparse 2*a} 0'
advected_interp_method = 'average'
[../]
[reaction]
type = FVReaction
variable = v
[]
[diff_v]
type = FVDiffusion
variable = v
coeff = ${diff}
[]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[exact]
type = FVFunctionDirichletBC
boundary = 'left right top bottom'
function = 'exact'
variable = v
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'sin(x)*cos(y)'
[]
[forcing]
type = ParsedFunction
value = '-2*a*sin(x)*sin(y) + a*cos(x)*cos(y) + 2*diff*sin(x)*cos(y) + sin(x)*cos(y)'
vars = 'a diff'
vals = '${a} ${diff}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type'
petsc_options_value = 'hypre'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/basic-primitive-pcnsfv-kt.i)
[GlobalParams]
fp = fp
limiter = 'central_difference'
two_term_boundary_expansion = true
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = .1
xmax = .6
nx = 2
[]
[]
[Modules]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
[]
[sup_vel_x]
type = MooseVariableFVReal
[]
[T_fluid]
type = MooseVariableFVReal
[]
[]
[ICs]
[pressure]
type = FunctionIC
variable = pressure
function = 'exact_p'
[]
[sup_vel_x]
type = FunctionIC
variable = sup_vel_x
function = 'exact_sup_vel_x'
[]
[T_fluid]
type = FunctionIC
variable = T_fluid
function = 'exact_T'
[]
[]
[FVKernels]
[mass_advection]
type = PCNSFVKT
variable = pressure
eqn = "mass"
[]
[mass_fn]
type = FVBodyForce
variable = pressure
function = 'forcing_rho'
[]
[momentum_x_advection]
type = PCNSFVKT
variable = sup_vel_x
momentum_component = x
eqn = "momentum"
[]
[momentum_fn]
type = FVBodyForce
variable = sup_vel_x
function = 'forcing_rho_ud'
[]
[fluid_energy_advection]
type = PCNSFVKT
variable = T_fluid
eqn = "energy"
[]
[energy_fn]
type = FVBodyForce
variable = T_fluid
function = 'forcing_rho_et'
[]
[]
[FVBCs]
[mass_left]
variable = pressure
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'mass'
[]
[momentum_left]
variable = sup_vel_x
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'momentum'
momentum_component = 'x'
[]
[energy_left]
variable = T_fluid
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'energy'
[]
[mass_right]
variable = pressure
type = PCNSFVStrongBC
boundary = right
eqn = 'mass'
pressure = 'exact_p'
[]
[momentum_right]
variable = sup_vel_x
type = PCNSFVStrongBC
boundary = right
eqn = 'momentum'
momentum_component = 'x'
pressure = 'exact_p'
[]
[energy_right]
variable = T_fluid
type = PCNSFVStrongBC
boundary = right
eqn = 'energy'
pressure = 'exact_p'
[]
# help gradient reconstruction
[pressure_right]
type = FVFunctionDirichletBC
variable = pressure
function = exact_p
boundary = 'right'
[]
[sup_vel_x_left]
type = FVFunctionDirichletBC
variable = sup_vel_x
function = exact_sup_vel_x
boundary = 'left'
[]
[T_fluid_left]
type = FVFunctionDirichletBC
variable = T_fluid
function = exact_T
boundary = 'left'
[]
[]
[Materials]
[var_mat]
type = PorousPrimitiveVarMaterial
pressure = pressure
superficial_vel_x = sup_vel_x
T_fluid = T_fluid
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[]
[Functions]
[exact_rho]
type = ParsedFunction
value = '3.48788261470924*cos(x)'
[]
[forcing_rho]
type = ParsedFunction
value = '-3.45300378856215*sin(1.1*x)'
[]
[exact_rho_ud]
type = ParsedFunction
value = '3.13909435323832*cos(1.1*x)'
[]
[forcing_rho_ud]
type = ParsedFunction
value = '-0.9*(10.6975765229419*cos(1.2*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + 0.9*(10.6975765229419*sin(x)*cos(1.2*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 12.8370918275302*sin(1.2*x)/cos(x))*cos(x) + 3.13909435323832*sin(x)*cos(1.1*x)^2/cos(x)^2 - 6.9060075771243*sin(1.1*x)*cos(1.1*x)/cos(x)'
[]
[exact_rho_et]
type = ParsedFunction
value = '26.7439413073546*cos(1.2*x)'
[]
[forcing_rho_et]
type = ParsedFunction
value = '0.9*(3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.2*x))*sin(x)*cos(1.1*x)/cos(x)^2 - 0.99*(3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.2*x))*sin(1.1*x)/cos(x) + 0.9*(-(10.6975765229419*cos(1.2*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.2*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 12.8370918275302*sin(1.2*x)/cos(x))*cos(x) - 32.0927295688256*sin(1.2*x))*cos(1.1*x)/cos(x)'
[]
[exact_T]
type = ParsedFunction
value = '0.0106975765229418*cos(1.2*x)/cos(x) - 0.000697576522941848*cos(1.1*x)^2/cos(x)^2'
[]
[exact_eps_p]
type = ParsedFunction
value = '3.13909435323832*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_p]
type = ParsedFunction
value = '3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_sup_vel_x]
type = ParsedFunction
value = '0.9*cos(1.1*x)/cos(x)'
[]
[exact_superficial_velocity]
type = ParsedVectorFunction
value_x = '0.9*cos(1.1*x)/cos(x)'
[]
[eps]
type = ParsedFunction
value = '0.9'
[]
[]
[Executioner]
solve_type = NEWTON
type = Transient
num_steps = 1
dtmin = 1
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_max_its = 50
line_search = bt
nl_rel_tol = 1e-12
nl_abs_tol = 1e-12
[]
[Outputs]
exodus = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2pressure]
type = ElementL2Error
variable = pressure
function = exact_p
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2sup_vel_x]
variable = sup_vel_x
function = exact_sup_vel_x
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2T_fluid]
variable = T_fluid
function = exact_T
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/skewness-correction/adv-diff-react/skewed.i)
a=1.1
diff=1.1
[Mesh]
[gen_mesh]
type = FileMeshGenerator
file = skewed.msh
[]
[]
[Variables]
[v]
initial_condition = 1
type = MooseVariableFVReal
face_interp_method = 'skewness-corrected'
cache_face_gradients = false
cache_face_values = true
[]
[]
[FVKernels]
[diff_v]
type = FVDiffusion
variable = v
coeff = ${diff}
[]
[advection]
type = FVAdvection
variable = v
velocity = '${a} ${fparse 2*a} 0'
advected_interp_method = 'average'
[]
[reaction]
type = FVReaction
variable = v
[]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[exact]
type = FVFunctionDirichletBC
boundary = 'left right top bottom'
function = 'exact'
variable = v
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'sin(x)*cos(y)'
[]
[forcing]
type = ParsedFunction
value = '-2*a*sin(x)*sin(y) + a*cos(x)*cos(y) + 2*diff*sin(x)*cos(y) + sin(x)*cos(y)'
vars = 'a diff'
vals = '${a} ${diff}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
[]
[]
(test/tests/fvkernels/mms/diffusion.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 2
[]
[Variables]
# [u]
# []
[v]
family = MONOMIAL
order = CONSTANT
fv = true
[]
[]
[FVKernels]
[diff_v]
type = FVDiffusion
variable = v
coeff = coeff
[]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[boundary]
type = FVFunctionDirichletBC
boundary = 'left right'
function = 'exact'
variable = v
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '1'
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = '3*x^2 + 2*x + 1'
[]
[forcing]
type = ParsedFunction
value = '-6'
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
# [./L2u]
# type = ElementL2Error
# variable = u
# function = exact
# outputs = 'console'
# execute_on = 'timestep_end'
# [../]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/cylindrical/advection-diffusion-reaction.i)
a=1.1
diff=1.1
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 2
xmin = 2
xmax = 3
ymin = 0
ymax = 1
nx = 2
ny = 2
[../]
[]
[Problem]
coord_type = 'RZ'
[]
[Variables]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 1
[../]
[]
[FVKernels]
[./advection]
type = FVAdvection
variable = v
velocity = '${a} ${a} 0'
advected_interp_method = 'average'
[../]
[reaction]
type = FVReaction
variable = v
[]
[diff_v]
type = FVDiffusion
variable = v
coeff = ${diff}
[]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[exact]
type = FVFunctionDirichletBC
boundary = 'left right top bottom'
function = 'exact'
variable = v
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'sin(x)*cos(y)'
[]
[forcing]
type = ParsedFunction
value = '-a*sin(x)*sin(y) + diff*sin(x)*cos(y) + sin(x)*cos(y) + (x*a*cos(x)*cos(y) + a*sin(x)*cos(y))/x - (-x*diff*sin(x)*cos(y) + diff*cos(x)*cos(y))/x'
vars = 'a diff'
vals = '${a} ${diff}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_shift_type -sub_pc_type'
petsc_options_value = 'asm NONZERO lu'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/mms/channel-flow/cylindrical/2d-average-with-temp.i)
mu=1.1
rho=1.1
k=1.1
cp=1.1
advected_interp_method='average'
velocity_interp_method='average'
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 2
ny = 2
[]
[]
[Problem]
fv_bcs_integrity_check = true
coord_type = 'RZ'
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = u
v = v
pressure = pressure
[]
[]
[Variables]
[u]
type = INSFVVelocityVariable
initial_condition = 1
two_term_boundary_expansion = false
[]
[v]
type = INSFVVelocityVariable
initial_condition = 1
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
two_term_boundary_expansion = false
[]
[temperature]
type = INSFVEnergyVariable
two_term_boundary_expansion = false
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[mass_forcing]
type = FVBodyForce
variable = pressure
function = forcing_p
[]
[u_advection]
type = INSFVMomentumAdvection
variable = u
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_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
[]
[u_forcing]
type = INSFVBodyForce
variable = u
functor = forcing_u
momentum_component = 'x'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = v
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_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
[]
[v_forcing]
type = INSFVBodyForce
variable = v
functor = forcing_v
momentum_component = 'y'
[]
[temp_conduction]
type = FVDiffusion
coeff = 'k'
variable = temperature
[]
[temp_advection]
type = INSFVEnergyAdvection
variable = temperature
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
[]
[temp_forcing]
type = FVBodyForce
variable = temperature
function = forcing_t
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'bottom'
variable = u
function = 'exact_u'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'bottom'
variable = v
function = 'exact_v'
[]
[no-slip-wall-u]
type = INSFVNoSlipWallBC
boundary = 'right'
variable = u
function = 'exact_u'
[]
[no-slip-wall-v]
type = INSFVNoSlipWallBC
boundary = 'right'
variable = v
function = 'exact_v'
[]
[outlet-p]
type = INSFVOutletPressureBC
boundary = 'top'
variable = pressure
function = 'exact_p'
[]
[axis-u]
type = INSFVSymmetryVelocityBC
boundary = 'left'
variable = u
u = u
v = v
mu = ${mu}
momentum_component = x
[]
[axis-v]
type = INSFVSymmetryVelocityBC
boundary = 'left'
variable = v
u = u
v = v
mu = ${mu}
momentum_component = y
[]
[axis-p]
type = INSFVSymmetryPressureBC
boundary = 'left'
variable = pressure
[]
[axis-inlet-wall-t]
type = FVFunctionDirichletBC
boundary = 'left bottom right'
variable = temperature
function = 'exact_t'
[]
[]
[Materials]
[const_functor]
type = ADGenericFunctorMaterial
prop_names = 'cp k'
prop_values = '${cp} ${k}'
[]
[ins_fv]
type = INSFVEnthalpyMaterial
temperature = 'temperature'
rho = ${rho}
[]
[]
[Functions]
[exact_u]
type = ParsedFunction
value = 'sin(x*pi)^2*sin((1/2)*y*pi)'
[]
[exact_rhou]
type = ParsedFunction
value = 'rho*sin(x*pi)^2*sin((1/2)*y*pi)'
vars = 'rho'
vals = '${rho}'
[]
[forcing_u]
type = ADParsedFunction
value = '(1/4)*pi^2*mu*sin(x*pi)^2*sin((1/2)*y*pi) - pi*sin(x*pi)*cos((1/2)*y*pi) + (4*x*pi*rho*sin(x*pi)^3*sin((1/2)*y*pi)^2*cos(x*pi) + rho*sin(x*pi)^4*sin((1/2)*y*pi)^2)/x + (-x*pi*rho*sin(x*pi)^2*sin((1/2)*y*pi)*sin(y*pi)*cos(x*pi) + (1/2)*x*pi*rho*sin(x*pi)^2*cos(x*pi)*cos((1/2)*y*pi)*cos(y*pi))/x - (-2*x*pi^2*mu*sin(x*pi)^2*sin((1/2)*y*pi) + 2*x*pi^2*mu*sin((1/2)*y*pi)*cos(x*pi)^2 + 2*pi*mu*sin(x*pi)*sin((1/2)*y*pi)*cos(x*pi))/x'
vars = 'mu rho'
vals = '${mu} ${rho}'
[]
[exact_v]
type = ParsedFunction
value = 'cos(x*pi)*cos(y*pi)'
[]
[exact_rhov]
type = ParsedFunction
value = 'rho*cos(x*pi)*cos(y*pi)'
vars = 'rho'
vals = '${rho}'
[]
[forcing_v]
type = ADParsedFunction
value = 'pi^2*mu*cos(x*pi)*cos(y*pi) - 2*pi*rho*sin(y*pi)*cos(x*pi)^2*cos(y*pi) - 1/2*pi*sin((1/2)*y*pi)*cos(x*pi) - (-x*pi^2*mu*cos(x*pi)*cos(y*pi) - pi*mu*sin(x*pi)*cos(y*pi))/x + (-x*pi*rho*sin(x*pi)^3*sin((1/2)*y*pi)*cos(y*pi) + 2*x*pi*rho*sin(x*pi)*sin((1/2)*y*pi)*cos(x*pi)^2*cos(y*pi) + rho*sin(x*pi)^2*sin((1/2)*y*pi)*cos(x*pi)*cos(y*pi))/x'
vars = 'mu rho'
vals = '${mu} ${rho}'
[]
[exact_p]
type = ParsedFunction
value = 'cos(x*pi)*cos((1/2)*y*pi)'
[]
[forcing_p]
type = ParsedFunction
value = '-pi*rho*sin(y*pi)*cos(x*pi) + (2*x*pi*rho*sin(x*pi)*sin((1/2)*y*pi)*cos(x*pi) + rho*sin(x*pi)^2*sin((1/2)*y*pi))/x'
vars = 'rho'
vals = '${rho}'
[]
[exact_t]
type = ParsedFunction
value = 'sin(x*pi)*sin((1/2)*y*pi)'
[]
[forcing_t]
type = ParsedFunction
value = '(1/4)*pi^2*k*sin(x*pi)*sin((1/2)*y*pi) - (-x*pi^2*k*sin(x*pi)*sin((1/2)*y*pi) + pi*k*sin((1/2)*y*pi)*cos(x*pi))/x + (3*x*pi*cp*rho*sin(x*pi)^2*sin((1/2)*y*pi)^2*cos(x*pi) + cp*rho*sin(x*pi)^3*sin((1/2)*y*pi)^2)/x + (-x*pi*cp*rho*sin(x*pi)*sin((1/2)*y*pi)*sin(y*pi)*cos(x*pi) + (1/2)*x*pi*cp*rho*sin(x*pi)*cos(x*pi)*cos((1/2)*y*pi)*cos(y*pi))/x'
vars = 'k rho cp'
vals = '${k} ${rho} ${cp}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
petsc_options_value = 'asm 100 lu NONZERO'
line_search = 'none'
[]
[Outputs]
exodus = true
csv = true
[dof]
type = DOFMap
execute_on = 'initial'
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[./L2u]
type = ElementL2Error
variable = u
function = exact_u
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2v]
type = ElementL2Error
variable = v
function = exact_v
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2p]
variable = pressure
function = exact_p
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2t]
variable = temperature
function = exact_t
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[]
(test/tests/fvkernels/mms/grad-reconstruction/mat-cartesian.i)
a=1.1
diff=1.1
[Mesh]
[gen_mesh]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 2
ny = 2
[]
[]
[Variables]
[v]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 1
[]
[]
[FVKernels]
[advection]
type = FVElementalAdvection
variable = v
velocity = '${a} ${fparse 2 * a} 0'
advected_quantity = 'mat_u'
grad_advected_quantity = 'mat_grad_u'
[]
[reaction]
type = FVReaction
variable = v
[]
[diff_v]
type = FVDiffusion
variable = v
coeff = ${diff}
[]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[diri]
type = FVFunctionDirichletBC
boundary = 'left right top bottom'
function = 'exact'
variable = v
[]
[]
[Materials]
[mat]
type = ADCoupledGradientMaterial
mat_prop = 'mat_u'
grad_mat_prop = 'mat_grad_u'
u = v
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'sin(x)*cos(y)'
[]
[forcing]
type = ParsedFunction
value = '-2*a*sin(x)*sin(y) + a*cos(x)*cos(y) + 2*diff*sin(x)*cos(y) + sin(x)*cos(y)'
vars = 'a diff'
vals = '${a} ${diff}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_shift_type -sub_pc_type'
petsc_options_value = 'asm NONZERO lu'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/postprocessors/element_integral_var_pps/pps_old_value_fv.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 4
ny = 4
elem_type = QUAD4
[]
[Variables]
[./u]
order = CONSTANT
family = MONOMIAL
fv = true
initial_condition = 1
[../]
[]
[Functions]
[./force_fn]
type = ParsedFunction
value = '1'
[../]
[./exact_fn]
type = ParsedFunction
value = 't'
[../]
[]
[FVKernels]
[./diff_u]
type = FVDiffusion
variable = u
coeff = '1'
block = '0'
[../]
[./ffn_u]
type = FVBodyForce
variable = u
function = force_fn
[../]
[]
[FVBCs]
[./all_u]
type = FVFunctionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
[../]
[]
[Postprocessors]
[./a]
type = ElementIntegralVariablePostprocessor
variable = u
execute_on = 'initial timestep_end'
[../]
[./total_a]
type = TimeIntegratedPostprocessor
value = a
execute_on = 'initial timestep_end'
[../]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Outputs]
execute_on = 'timestep_end'
exodus = true
[]
(modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/varying-eps-basic-kt-mixed.i)
[GlobalParams]
fp = fp
limiter = 'central_difference'
two_term_boundary_expansion = true
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = .1
xmax = .6
nx = 2
[]
[]
[Modules]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = MooseVariableFVReal
[]
[sup_mom_x]
type = MooseVariableFVReal
[]
[T_fluid]
type = MooseVariableFVReal
[]
[]
[ICs]
[pressure]
type = FunctionIC
variable = pressure
function = 'exact_p'
[]
[sup_mom_x]
type = FunctionIC
variable = sup_mom_x
function = 'exact_rho_ud'
[]
[T_fluid]
type = FunctionIC
variable = T_fluid
function = 'exact_T'
[]
[]
[FVKernels]
[mass_advection]
type = PCNSFVKT
variable = pressure
eqn = "mass"
[]
[mass_fn]
type = FVBodyForce
variable = pressure
function = 'forcing_rho'
[]
[momentum_x_advection]
type = PCNSFVKT
variable = sup_mom_x
momentum_component = x
eqn = "momentum"
[]
[eps_grad]
type = PNSFVPGradEpsilon
variable = sup_mom_x
momentum_component = 'x'
epsilon_function = 'eps'
[]
[momentum_fn]
type = FVBodyForce
variable = sup_mom_x
function = 'forcing_rho_ud'
[]
[fluid_energy_advection]
type = PCNSFVKT
variable = T_fluid
eqn = "energy"
[]
[energy_fn]
type = FVBodyForce
variable = T_fluid
function = 'forcing_rho_et'
[]
[]
[FVBCs]
[mass_left]
variable = pressure
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'mass'
[]
[momentum_left]
variable = sup_mom_x
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'momentum'
momentum_component = 'x'
[]
[energy_left]
variable = T_fluid
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'energy'
[]
[mass_right]
variable = pressure
type = PCNSFVStrongBC
boundary = right
eqn = 'mass'
pressure = 'exact_p'
[]
[momentum_right]
variable = sup_mom_x
type = PCNSFVStrongBC
boundary = right
eqn = 'momentum'
momentum_component = 'x'
pressure = 'exact_p'
[]
[energy_right]
variable = T_fluid
type = PCNSFVStrongBC
boundary = right
eqn = 'energy'
pressure = 'exact_p'
[]
# help gradient reconstruction
[pressure_right]
type = FVFunctionDirichletBC
variable = pressure
function = exact_p
boundary = 'right'
[]
[sup_mom_x_left]
type = FVFunctionDirichletBC
variable = sup_mom_x
function = exact_rho_ud
boundary = 'left'
[]
[T_fluid_left]
type = FVFunctionDirichletBC
variable = T_fluid
function = exact_T
boundary = 'left'
[]
[]
[Materials]
[var_mat]
type = PorousMixedVarMaterial
pressure = pressure
superficial_rhou = sup_mom_x
T_fluid = T_fluid
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[]
[Functions]
[exact_rho]
type = ParsedFunction
value = '3.48788261470924*cos(x)'
[]
[forcing_rho]
type = ParsedFunction
value = '-3.83667087618017*sin(1.1*x)*cos(1.3*x) - 4.53424739912202*sin(1.3*x)*cos(1.1*x)'
[]
[exact_rho_ud]
type = ParsedFunction
value = '3.48788261470924*cos(1.1*x)*cos(1.3*x)'
[]
[forcing_rho_ud]
type = ParsedFunction
value = '(-(10.6975765229419*cos(1.5*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.5*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 16.0463647844128*sin(1.5*x)/cos(x))*cos(x))*cos(1.3*x) + 3.48788261470924*sin(x)*cos(1.1*x)^2*cos(1.3*x)/cos(x)^2 - 7.67334175236034*sin(1.1*x)*cos(1.1*x)*cos(1.3*x)/cos(x) - 4.53424739912202*sin(1.3*x)*cos(1.1*x)^2/cos(x)'
[]
[exact_rho_et]
type = ParsedFunction
value = '26.7439413073546*cos(1.5*x)'
[]
[forcing_rho_et]
type = ParsedFunction
value = '1.0*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(x)*cos(1.1*x)*cos(1.3*x)/cos(x)^2 - 1.1*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(1.1*x)*cos(1.3*x)/cos(x) - 1.3*(3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.5*x))*sin(1.3*x)*cos(1.1*x)/cos(x) + 1.0*(-(10.6975765229419*cos(1.5*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.5*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 16.0463647844128*sin(1.5*x)/cos(x))*cos(x) - 40.1159119610319*sin(1.5*x))*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[exact_T]
type = ParsedFunction
value = '0.0106975765229418*cos(1.5*x)/cos(x) - 0.000697576522941848*cos(1.1*x)^2/cos(x)^2'
[]
[exact_eps_p]
type = ParsedFunction
value = '3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)*cos(1.3*x)'
[]
[exact_p]
type = ParsedFunction
value = '3.48788261470924*(3.06706896551724*cos(1.5*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_sup_vel_x]
type = ParsedFunction
value = '1.0*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[eps]
type = ParsedFunction
value = 'cos(1.3*x)'
[]
[exact_superficial_velocity]
type = ParsedVectorFunction
value_x = '1.0*cos(1.1*x)*cos(1.3*x)/cos(x)'
[]
[]
[Executioner]
solve_type = NEWTON
type = Transient
num_steps = 1
dtmin = 1
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_max_its = 50
line_search = bt
nl_rel_tol = 1e-12
nl_abs_tol = 1e-12
[]
[Outputs]
exodus = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2pressure]
type = ElementL2Error
variable = pressure
function = exact_p
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2sup_mom_x]
variable = sup_mom_x
function = exact_rho_ud
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2T_fluid]
variable = T_fluid
function = exact_T
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/non-orthogonal/extended-adr.i)
a=1.1
diff=1.1
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 2
xmin = 2
xmax = 3
ymin = 0
ymax = 1
nx = 2
ny = 2
elem_type = TRI3
[../]
[]
[Variables]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 1
type = MooseVariableFVReal
face_interp_method = vertex-based
[../]
[]
[FVKernels]
[./advection]
type = FVAdvection
variable = v
velocity = '${a} ${fparse 2*a} 0'
advected_interp_method = 'average'
[../]
[reaction]
type = FVReaction
variable = v
[]
[diff_v]
type = FVDiffusion
variable = v
coeff = ${diff}
use_point_neighbors = true
[]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[exact]
type = FVFunctionDirichletBC
boundary = 'left right top bottom'
function = 'exact'
variable = v
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'sin(x)*cos(y)'
[]
[forcing]
type = ParsedFunction
value = '-2*a*sin(x)*sin(y) + a*cos(x)*cos(y) + 2*diff*sin(x)*cos(y) + sin(x)*cos(y)'
vars = 'a diff'
vals = '${a} ${diff}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type'
petsc_options_value = 'hypre'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/pins/mms/1d-rc-no-diffusion-strong-bc.i)
mu=1e-15
rho=1.1
advected_interp_method='upwind'
velocity_interp_method='rc'
[Mesh]
[mesh]
type = GeneratedMeshGenerator
dim = 1
nx = 2
xmax = 0.5
[]
[]
[GlobalParams]
two_term_boundary_expansion = true
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = PINSFVRhieChowInterpolator
u = u
pressure = pressure
porosity = porosity
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[u]
type = PINSFVSuperficialVelocityVariable
initial_condition = .1
[]
[pressure]
type = INSFVPressureVariable
[]
[]
[AuxVariables]
[porosity]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 0.8
[]
[]
[Problem]
error_on_jacobian_nonzero_reallocation = true
[]
[Functions]
[exact_u]
type = ParsedFunction
value = 'cos((1/2)*x*pi)'
[]
[forcing_u]
type = ADParsedFunction
value = '-1.25*pi*rho*sin((1/2)*x*pi)*cos((1/2)*x*pi) + 0.8*cos(x)'
vars = 'mu rho'
vals = '${mu} ${rho}'
[]
[exact_p]
type = ParsedFunction
value = 'sin(x)'
[]
[forcing_p]
type = ParsedFunction
value = '-1/2*pi*rho*sin((1/2)*x*pi)'
vars = 'rho'
vals = '${rho}'
[]
[]
[FVKernels]
[mass]
type = PINSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[mass_forcing]
type = FVBodyForce
variable = pressure
function = forcing_p
[]
[u_advection]
type = PINSFVMomentumAdvection
variable = u
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
porosity = porosity
momentum_component = 'x'
[]
[u_pressure]
type = PINSFVMomentumPressureFlux
variable = u
pressure = pressure
porosity = porosity
momentum_component = 'x'
force_boundary_execution = false
[]
[u_forcing]
type = INSFVBodyForce
variable = u
functor = forcing_u
momentum_component = 'x'
[]
[]
[FVBCs]
[mass]
variable = pressure
type = PINSFVFunctorBC
boundary = 'left right'
superficial_vel_x = u
pressure = pressure
eqn = 'mass'
porosity = porosity
[]
[momentum]
variable = u
type = PINSFVFunctorBC
boundary = 'left right'
superficial_vel_x = u
pressure = pressure
eqn = 'momentum'
momentum_component = 'x'
porosity = porosity
[]
[inlet-u]
type = FVFunctionDirichletBC
boundary = 'left'
variable = u
function = 'exact_u'
[]
[outlet_p]
type = FVFunctionDirichletBC
boundary = 'right'
variable = pressure
function = 'exact_p'
[]
[]
[Materials]
[const]
type = ADGenericFunctorMaterial
prop_names = 'rho'
prop_values = '${rho}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
petsc_options_value = 'asm 100 lu NONZERO'
line_search = 'bt'
[]
[Postprocessors]
[inlet_p]
type = SideAverageValue
variable = 'pressure'
boundary = 'left'
[]
[outlet-u]
type = SideIntegralVariablePostprocessor
variable = u
boundary = 'right'
[]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2u]
type = ElementL2Error
variable = u
function = exact_u
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2p]
variable = pressure
function = exact_p
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
[Outputs]
exodus = true
csv = true
[]
(test/tests/fvkernels/mms/advective-outflow/advection.i)
a=1.1
[GlobalParams]
advected_interp_method = 'average'
[]
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 1
xmin = 0.1
xmax = 1.1
nx = 2
[../]
[]
[Variables]
[./u]
family = MONOMIAL
order = CONSTANT
fv = true
two_term_boundary_expansion = false
type = MooseVariableFVReal
[../]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
two_term_boundary_expansion = true
type = MooseVariableFVReal
[../]
[]
[FVKernels]
[./advection_u]
type = FVAdvection
variable = u
velocity = '${a} 0 0'
force_boundary_execution = true
[../]
[body_u]
type = FVBodyForce
variable = u
function = 'forcing'
[]
[./advection_v]
type = FVAdvection
variable = v
velocity = '${a} 0 0'
force_boundary_execution = true
[../]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[left_u]
type = FVFunctionDirichletBC
boundary = 'left'
function = 'exact'
variable = u
[]
[left_v]
type = FVFunctionDirichletBC
boundary = 'left'
function = 'exact'
variable = v
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'cos(x)'
[]
[forcing]
type = ParsedFunction
value = '-${a} * sin(x)'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./L2u]
type = ElementL2Error
variable = u
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2v]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/advective-outflow/advection-outflow.i)
a=1.1
[GlobalParams]
advected_interp_method = 'average'
[]
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 1
xmin = 0
xmax = 1.1
nx = 2
[../]
[]
[Variables]
[./u]
family = MONOMIAL
order = CONSTANT
fv = true
two_term_boundary_expansion = false
type = MooseVariableFVReal
[../]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
two_term_boundary_expansion = true
type = MooseVariableFVReal
[../]
[]
[FVKernels]
[./advection_u]
type = FVAdvection
variable = u
velocity = '${a} 0 0'
[../]
[body_u]
type = FVBodyForce
variable = u
function = 'forcing'
[]
[./advection_v]
type = FVAdvection
variable = v
velocity = '${a} 0 0'
[../]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[left_u]
type = FVFunctionDirichletBC
boundary = 'left'
function = 'exact'
variable = u
[]
[right_u]
type = FVConstantScalarOutflowBC
variable = u
velocity = '${a} 0 0'
boundary = 'right'
[]
[left_v]
type = FVFunctionDirichletBC
boundary = 'left'
function = 'exact'
variable = v
[]
[right_v]
type = FVConstantScalarOutflowBC
variable = v
velocity = '${a} 0 0'
boundary = 'right'
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'cos(x)'
[]
[forcing]
type = ParsedFunction
value = '-${a} * sin(x)'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./L2u]
type = ElementL2Error
variable = u
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2v]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/skewness-correction/two_term_extrapol/advection-outflow.i)
diff=1
a=1
[GlobalParams]
advected_interp_method = 'average'
[]
[Mesh]
[./gen_mesh]
type = FileMeshGenerator
file = skewed.msh
[../]
[]
[Variables]
[./v]
type = MooseVariableFVReal
face_interp_method = 'skewness-corrected'
cache_face_gradients = false
cache_face_values = true
[../]
[]
[FVKernels]
[./advection]
type = FVAdvection
variable = v
velocity = '${a} 0 0'
[../]
[./diffusion]
type = FVDiffusion
variable = v
coeff = coeff
[../]
[./body]
type = FVBodyForce
variable = v
function = 'forcing'
[../]
[]
[FVBCs]
[left]
type = FVFunctionDirichletBC
boundary = 'left'
function = 'exact'
variable = v
[]
[top]
type = FVNeumannBC
boundary = 'top'
value = 0
variable = v
[]
[bottom]
type = FVNeumannBC
boundary = 'bottom'
value = 0
variable = v
[]
[right]
type = FVConstantScalarOutflowBC
variable = v
velocity = '${a} 0 0'
boundary = 'right'
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '${diff}'
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'cos(x)'
[]
[forcing]
type = ParsedFunction
value = 'cos(x) - sin(x)'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_hypre_type -snes_linesearch_minlambda'
petsc_options_value = 'hypre boomeramg 1e-9'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/skewness-correction/diffusion/skewed.i)
a=1.1
diff=1.1
[Mesh]
[./gen_mesh]
type = FileMeshGenerator
file = skewed.msh
[../]
[]
[Variables]
[./v]
initial_condition = 1
type = MooseVariableFVReal
face_interp_method = 'skewness-corrected'
cache_face_gradients = false
cache_face_values = true
[../]
[]
[FVKernels]
[diff_v]
type = FVDiffusion
variable = v
coeff = ${diff}
[]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[exact]
type = FVFunctionDirichletBC
boundary = 'left right top bottom'
function = 'exact'
variable = v
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'sin(x)*cos(y)'
[]
[forcing]
type = ParsedFunction
value = '2*diff*sin(x)*cos(y)'
vars = 'a diff'
vals = '${a} ${diff}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
[]
[]
(test/tests/fvkernels/mms/grad-reconstruction/rz.i)
a=1.1
diff=1.1
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 2
xmin = 2
xmax = 3
ymin = 0
ymax = 1
nx = 2
ny = 2
[../]
[]
[Problem]
coord_type = 'RZ'
[]
[Variables]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 1
[../]
[]
[FVKernels]
[./advection]
type = FVElementalAdvection
variable = v
velocity = '${a} ${a} 0'
[../]
[reaction]
type = FVReaction
variable = v
[]
[diff_v]
type = FVDiffusion
variable = v
coeff = ${diff}
[]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[diri]
type = FVFunctionDirichletBC
boundary = 'left right top bottom'
function = 'exact'
variable = v
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'sin(x)*cos(y)'
[]
[forcing]
type = ParsedFunction
value = '-a*sin(x)*sin(y) + diff*sin(x)*cos(y) + sin(x)*cos(y) + (x*a*cos(x)*cos(y) + a*sin(x)*cos(y))/x - (-x*diff*sin(x)*cos(y) + diff*cos(x)*cos(y))/x'
vars = 'a diff'
vals = '${a} ${diff}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_shift_type -sub_pc_type'
petsc_options_value = 'asm NONZERO lu'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/grad-reconstruction/cartesian.i)
a=1.1
diff=1.1
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 2
ny = 2
[../]
[]
[Variables]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
initial_condition = 1
[../]
[]
[FVKernels]
[./advection]
type = FVElementalAdvection
variable = v
velocity = '${a} ${fparse 2 * a} 0'
[../]
[reaction]
type = FVReaction
variable = v
[]
[diff_v]
type = FVDiffusion
variable = v
coeff = ${diff}
[]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[diri]
type = FVFunctionDirichletBC
boundary = 'left right top bottom'
function = 'exact'
variable = v
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'sin(x)*cos(y)'
[]
[forcing]
type = ParsedFunction
value = '-2*a*sin(x)*sin(y) + a*cos(x)*cos(y) + 2*diff*sin(x)*cos(y) + sin(x)*cos(y)'
vars = 'a diff'
vals = '${a} ${diff}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -sub_pc_factor_shift_type -sub_pc_type'
petsc_options_value = 'asm NONZERO lu'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/advective-outflow/kt-limited-advection.i)
a=1.1
c=343
max_abs_eig=${fparse c + a}
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 1
xmin = 0.1
xmax = 1.1
nx = 2
[../]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[ICs]
[u]
type = FunctionIC
variable = u
function = exact
[]
[]
[Variables]
[./u]
two_term_boundary_expansion = true
type = MooseVariableFVReal
[../]
[]
[FVKernels]
[./advection_u]
type = FVKTLimitedAdvection
variable = u
velocity = '${a} 0 0'
limiter = 'vanLeer'
max_abs_eig = ${max_abs_eig}
add_artificial_diff = true
[../]
[body_u]
type = FVBodyForce
variable = u
function = 'forcing'
[]
[]
[FVBCs]
[left_u]
type = FVFunctionNeumannBC
boundary = 'left'
function = 'advection'
variable = u
[]
[diri_left]
type = FVFunctionDirichletBC
boundary = 'left'
function = 'exact'
variable = u
[]
[right]
type = FVConstantScalarOutflowBC
variable = u
velocity = '${a} 0 0'
boundary = 'right'
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 'cos(x)'
[]
[advection]
type = ParsedFunction
value = '${a} * cos(x)'
[]
[forcing]
type = ParsedFunction
value = '-${a} * sin(x)'
[]
[]
[Executioner]
type = Steady
petsc_options_iname = '-snes_linesearch_minlambda'
petsc_options_value = '1e-3'
nl_abs_tol = 1e-9
[]
[Outputs]
file_base = 'kt-limited-advection_out'
[csv]
type = CSV
execute_on = 'final'
[]
[exo]
type = Exodus
execute_on = 'final'
[]
[]
[Postprocessors]
[./L2u]
type = ElementL2Error
variable = u
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/advection-diffusion.i)
diff=1.1
a=1.1
[GlobalParams]
advected_interp_method = 'average'
[]
[Mesh]
[./gen_mesh]
type = GeneratedMeshGenerator
dim = 1
xmin = -0.6
xmax = 0.6
nx = 64
[../]
[]
[Variables]
[./v]
family = MONOMIAL
order = CONSTANT
fv = true
[../]
[]
[FVKernels]
[./advection]
type = FVAdvection
variable = v
velocity = '${a} 0 0'
[../]
[./diffusion]
type = FVDiffusion
variable = v
coeff = coeff
[../]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[boundary]
type = FVFunctionDirichletBC
boundary = 'left right'
function = 'exact'
variable = v
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '${diff}'
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = '3*x^2 + 2*x + 1'
[]
[forcing]
type = ParsedFunction
value = '-${diff}*6 + ${a} * (6*x + 2)'
# value = '-${diff}*6'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/cns/mms/1d/hllc-mms.i)
[GlobalParams]
fp = fp
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 1
xmin = .1
xmax = 1
nx = 2
[]
[]
[Variables]
[rho]
type = MooseVariableFVReal
[]
[rho_u]
type = MooseVariableFVReal
[]
[rho_et]
type = MooseVariableFVReal
[]
[]
[ICs]
[rho]
type = FunctionIC
variable = rho
function = 'cos(1.1*x)'
[]
[rho_u]
type = FunctionIC
variable = rho_u
function = '2*sin(1.1*x)'
[]
[rho_et]
type = FunctionIC
variable = rho_et
function = '3*cos(1.1*x)'
[]
[]
[FVKernels]
[mass_advection]
type = CNSFVMassHLLC
variable = rho
[]
[fake_diffusivity]
type = FVDiffusion
variable = rho
coeff = 1
[]
[mass_fn]
type = FVBodyForce
variable = rho
function = 'forcing_rho'
[]
[momentum_advection]
type = CNSFVMomentumHLLC
variable = rho_u
momentum_component = x
[]
[viscosity]
type = FVDiffusion
variable = rho_u
coeff = 1
[]
[momentum_fn]
type = FVBodyForce
variable = rho_u
function = 'forcing_rho_u'
[]
[fluid_energy_advection]
type = CNSFVFluidEnergyHLLC
variable = rho_et
[]
[fake_conduction]
type = FVDiffusion
variable = rho_et
coeff = 1
[]
[energy_fn]
type = FVBodyForce
variable = rho_et
function = 'forcing_rho_et'
[]
[]
[FVBCs]
[rho]
type = FVFunctionDirichletBC
boundary = 'left right'
variable = rho
function = 'exact_rho'
[]
[rho_u]
type = FVFunctionDirichletBC
boundary = 'left right'
variable = rho_u
function = 'exact_rho_u'
[]
[rho_et]
type = FVFunctionDirichletBC
boundary = 'left right'
variable = rho_et
function = 'exact_rho_et'
[]
[]
[Materials]
[var_mat]
type = ConservedVarValuesMaterial
rho = rho
rhou = rho_u
rho_et = rho_et
[]
[]
[Modules]
[FluidProperties]
[fp]
type = TestConservedVarFluidProperties
[]
[]
[]
[Functions]
[exact_rho]
type = ParsedFunction
value = 'cos(x)'
[]
[forcing_rho]
type = ParsedFunction
value = '3*cos(x)'
[]
[exact_rho_u]
type = ParsedFunction
value = '2*sin(x)'
[]
[forcing_rho_u]
type = ParsedFunction
value = '4*sin(x)^3/cos(x)^2 + 9*sin(x)'
[]
[exact_rho_et]
type = ParsedFunction
value = '3*cos(x)'
[]
[forcing_rho_et]
type = ParsedFunction
value = '11*cos(x)'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2rho]
type = ElementL2Error
variable = rho
function = exact_rho
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2rho_u]
variable = rho_u
function = exact_rho_u
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2rho_et]
variable = rho_et
function = exact_rho_et
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/basic-conserved-pcnsfv-kt.i)
[GlobalParams]
fp = fp
limiter = 'central_difference'
two_term_boundary_expansion = true
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = .1
xmax = .6
nx = 2
[]
[]
[Modules]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[rho]
type = MooseVariableFVReal
[]
[rho_ud]
type = MooseVariableFVReal
[]
[rho_et]
type = MooseVariableFVReal
[]
[]
[ICs]
[pressure]
type = FunctionIC
variable = rho
function = 'exact_rho'
[]
[sup_vel_x]
type = FunctionIC
variable = rho_ud
function = 'exact_rho_ud'
[]
[T_fluid]
type = FunctionIC
variable = rho_et
function = 'exact_rho_et'
[]
[]
[FVKernels]
[mass_advection]
type = PCNSFVKT
variable = rho
eqn = "mass"
[]
[mass_fn]
type = FVBodyForce
variable = rho
function = 'forcing_rho'
[]
[momentum_x_advection]
type = PCNSFVKT
variable = rho_ud
momentum_component = x
eqn = "momentum"
[]
[momentum_fn]
type = FVBodyForce
variable = rho_ud
function = 'forcing_rho_ud'
[]
[fluid_energy_advection]
type = PCNSFVKT
variable = rho_et
eqn = "energy"
[]
[energy_fn]
type = FVBodyForce
variable = rho_et
function = 'forcing_rho_et'
[]
[]
[FVBCs]
[mass_left]
variable = rho
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'mass'
[]
[momentum_left]
variable = rho_ud
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'momentum'
momentum_component = 'x'
[]
[energy_left]
variable = rho_et
type = PCNSFVStrongBC
boundary = left
T_fluid = 'exact_T'
superficial_velocity = 'exact_superficial_velocity'
eqn = 'energy'
[]
[mass_right]
variable = rho
type = PCNSFVStrongBC
boundary = right
eqn = 'mass'
pressure = 'exact_p'
[]
[momentum_right]
variable = rho_ud
type = PCNSFVStrongBC
boundary = right
eqn = 'momentum'
momentum_component = 'x'
pressure = 'exact_p'
[]
[energy_right]
variable = rho_et
type = PCNSFVStrongBC
boundary = right
eqn = 'energy'
pressure = 'exact_p'
[]
# help gradient reconstruction
[rho_right]
type = FVFunctionDirichletBC
variable = rho
function = exact_rho
boundary = 'right'
[]
[rho_ud_left]
type = FVFunctionDirichletBC
variable = rho_ud
function = exact_rho_ud
boundary = 'left'
[]
[rho_et_left]
type = FVFunctionDirichletBC
variable = rho_et
function = exact_rho_et
boundary = 'left'
[]
[]
[Materials]
[var_mat]
type = PorousConservedVarMaterial
rho = rho
superficial_rhou = rho_ud
rho_et = rho_et
porosity = porosity
[]
[porosity]
type = GenericFunctionMaterial
prop_names = 'porosity'
prop_values = 'eps'
[]
[]
[Functions]
[exact_rho]
type = ParsedFunction
value = '3.48788261470924*cos(x)'
[]
[forcing_rho]
type = ParsedFunction
value = '-3.45300378856215*sin(1.1*x)'
[]
[exact_rho_ud]
type = ParsedFunction
value = '3.13909435323832*cos(1.1*x)'
[]
[forcing_rho_ud]
type = ParsedFunction
value = '-0.9*(10.6975765229419*cos(1.2*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + 0.9*(10.6975765229419*sin(x)*cos(1.2*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 12.8370918275302*sin(1.2*x)/cos(x))*cos(x) + 3.13909435323832*sin(x)*cos(1.1*x)^2/cos(x)^2 - 6.9060075771243*sin(1.1*x)*cos(1.1*x)/cos(x)'
[]
[exact_rho_et]
type = ParsedFunction
value = '26.7439413073546*cos(1.2*x)'
[]
[forcing_rho_et]
type = ParsedFunction
value = '0.9*(3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.2*x))*sin(x)*cos(1.1*x)/cos(x)^2 - 0.99*(3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x) + 26.7439413073546*cos(1.2*x))*sin(1.1*x)/cos(x) + 0.9*(-(10.6975765229419*cos(1.2*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + (10.6975765229419*sin(x)*cos(1.2*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 12.8370918275302*sin(1.2*x)/cos(x))*cos(x) - 32.0927295688256*sin(1.2*x))*cos(1.1*x)/cos(x)'
[]
[exact_T]
type = ParsedFunction
value = '0.0106975765229418*cos(1.2*x)/cos(x) - 0.000697576522941848*cos(1.1*x)^2/cos(x)^2'
[]
[exact_eps_p]
type = ParsedFunction
value = '3.13909435323832*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_p]
type = ParsedFunction
value = '3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_sup_vel_x]
type = ParsedFunction
value = '0.9*cos(1.1*x)/cos(x)'
[]
[exact_superficial_velocity]
type = ParsedVectorFunction
value_x = '0.9*cos(1.1*x)/cos(x)'
[]
[eps]
type = ParsedFunction
value = '0.9'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
solve_type = NEWTON
type = Transient
num_steps = 1
dtmin = 1
petsc_options = '-snes_linesearch_monitor'
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_max_its = 50
line_search = bt
[]
[Outputs]
exodus = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2rho]
type = ElementL2Error
variable = rho
function = exact_rho
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2rho_ud]
variable = rho_ud
function = exact_rho_ud
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2rho_et]
variable = rho_et
function = exact_rho_et
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/cns/mms/1d-with-bcs/pwcnsfv.i)
rho='rho'
advected_interp_method='upwind'
velocity_interp_method='rc'
gamma=1.4
R=8.3145
molar_mass=29.0e-3
R_specific=${fparse R/molar_mass}
cp=${fparse gamma*R_specific/(gamma-1)}
[GlobalParams]
two_term_boundary_expansion = true
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = PINSFVRhieChowInterpolator
u = sup_vel_x
pressure = pressure
porosity = porosity
[]
[]
[Mesh]
[cartesian]
type = GeneratedMeshGenerator
dim = 1
xmin = .1
xmax = .6
nx = 2
[]
[]
[Modules]
[FluidProperties]
[fp]
type = IdealGasFluidProperties
[]
[]
[]
[Problem]
fv_bcs_integrity_check = false
[]
[Variables]
[pressure]
type = INSFVPressureVariable
[]
[sup_vel_x]
type = PINSFVSuperficialVelocityVariable
[]
[]
[AuxVariables]
[porosity]
type = MooseVariableFVReal
[]
[T_fluid]
type = INSFVEnergyVariable
[]
[]
[ICs]
[pressure]
type = FunctionIC
variable = pressure
function = 'exact_p'
[]
[sup_vel_x]
type = FunctionIC
variable = sup_vel_x
function = 'exact_sup_vel_x'
[]
[T_fluid]
type = FunctionIC
variable = T_fluid
function = 'exact_T'
[]
[eps]
type = FunctionIC
variable = porosity
function = 'eps'
[]
[]
[FVKernels]
[mass_advection]
type = PINSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[mass_fn]
type = FVBodyForce
variable = pressure
function = 'forcing_rho'
[]
[u_advection]
type = PINSFVMomentumAdvection
variable = sup_vel_x
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
porosity = porosity
momentum_component = 'x'
[]
[u_pressure]
type = PINSFVMomentumPressureFlux
variable = sup_vel_x
pressure = pressure
porosity = porosity
momentum_component = 'x'
force_boundary_execution = false
[]
[momentum_fn]
type = INSFVBodyForce
variable = sup_vel_x
functor = 'forcing_rho_ud'
momentum_component = 'x'
[]
[]
[FVBCs]
[mass]
variable = pressure
type = PINSFVFunctorBC
boundary = 'left right'
superficial_vel_x = sup_vel_x
pressure = pressure
eqn = 'mass'
porosity = porosity
[]
[momentum]
variable = sup_vel_x
type = PINSFVFunctorBC
boundary = 'left right'
superficial_vel_x = sup_vel_x
pressure = pressure
eqn = 'momentum'
momentum_component = 'x'
porosity = porosity
[]
# help gradient reconstruction *and* create Dirichlet values for use in PINSFVFunctorBC
[pressure_right]
type = FVFunctionDirichletBC
variable = pressure
function = exact_p
boundary = 'right'
[]
[sup_vel_x_left]
type = FVFunctionDirichletBC
variable = sup_vel_x
function = exact_sup_vel_x
boundary = 'left'
[]
[T_fluid_left]
type = FVFunctionDirichletBC
variable = T_fluid
function = exact_T
boundary = 'left'
[]
[]
[Materials]
[const_functor]
type = ADGenericFunctorMaterial
prop_names = 'cp'
prop_values = '${cp}'
[]
[rho]
type = RhoFromPTFunctorMaterial
fp = fp
temperature = T_fluid
pressure = pressure
[]
[ins_fv]
type = INSFVEnthalpyMaterial
temperature = T_fluid
rho = ${rho}
[]
[]
[Functions]
[forcing_rho]
type = ParsedFunction
value = '-3.45300378856215*sin(1.1*x)'
[]
[forcing_rho_ud]
type = ADParsedFunction
value = '-0.9*(10.6975765229419*cos(1.2*x)/cos(x) - 0.697576522941849*cos(1.1*x)^2/cos(x)^2)*sin(x) + 0.9*(10.6975765229419*sin(x)*cos(1.2*x)/cos(x)^2 - 1.3951530458837*sin(x)*cos(1.1*x)^2/cos(x)^3 + 1.53466835047207*sin(1.1*x)*cos(1.1*x)/cos(x)^2 - 12.8370918275302*sin(1.2*x)/cos(x))*cos(x) + 3.13909435323832*sin(x)*cos(1.1*x)^2/cos(x)^2 - 6.9060075771243*sin(1.1*x)*cos(1.1*x)/cos(x)'
[]
[exact_T]
type = ParsedFunction
value = '0.0106975765229418*cos(1.2*x)/cos(x) - 0.000697576522941848*cos(1.1*x)^2/cos(x)^2'
[]
[exact_p]
type = ParsedFunction
value = '3.48788261470924*(3.06706896551724*cos(1.2*x)/cos(x) - 0.2*cos(1.1*x)^2/cos(x)^2)*cos(x)'
[]
[exact_sup_vel_x]
type = ParsedFunction
value = '0.9*cos(1.1*x)/cos(x)'
[]
[eps]
type = ParsedFunction
value = '0.9'
[]
[]
[Executioner]
solve_type = NEWTON
type = Transient
num_steps = 1
dtmin = 1
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
nl_max_its = 50
line_search = bt
nl_rel_tol = 1e-12
nl_abs_tol = 1e-12
[]
[Outputs]
exodus = true
csv = true
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2pressure]
type = ElementL2Error
variable = pressure
function = exact_p
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2sup_vel_x]
variable = sup_vel_x
function = exact_sup_vel_x
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/postprocessors/interface_diffusive_flux/interface_diffusive_flux_fv.i)
postprocessor_type = InterfaceDiffusiveFluxAverage
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 6
xmax = 3
ny = 9
ymax = 3
elem_type = QUAD4
[]
[subdomain_id]
input = gen
type = SubdomainBoundingBoxGenerator
bottom_left = '0 0 0'
top_right = '2 1 0'
block_id = 1
[]
[interface]
input = subdomain_id
type = SideSetsBetweenSubdomainsGenerator
primary_block = '1'
paired_block = '0'
new_boundary = 'interface'
[]
[]
[Functions]
[fn_exact]
type = ParsedFunction
value = 'x*x+y*y'
[]
[]
[Variables]
[u]
type = MooseVariableFVReal
block = 0
[]
[v]
type = MooseVariableFVReal
block = 1
[]
[]
[FVKernels]
[diff_u]
type = FVDiffusion
variable = u
coeff = 1
[]
[body_u]
type = FVBodyForce
variable = u
function = 1
[]
[diff_v]
type = FVDiffusion
variable = v
coeff = 1
[]
[body_v]
type = FVBodyForce
variable = v
function = -1
[]
[]
[FVInterfaceKernels]
[reaction]
type = FVDiffusionInterface
variable1 = u
variable2 = v
coeff1 = 1
coeff2 = 2
boundary = 'interface'
subdomain1 = '0'
subdomain2 = '1'
[]
[]
[FVBCs]
[all]
type = FVFunctionDirichletBC
variable = u
boundary = '0 1 2 3'
function = fn_exact
[]
[]
[Postprocessors]
[diffusive_flux]
type = ${postprocessor_type}
variable = v
neighbor_variable = u
diffusivity = 1
execute_on = TIMESTEP_END
boundary = 'interface'
[]
[]
[Executioner]
type = Steady
solve_type = NEWTON
[]
[Outputs]
file_base = ${raw ${postprocessor_type} _fv}
exodus = true
[]
(test/tests/fvkernels/mms/cylindrical/diffusion.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
[]
[Variables]
[v]
family = MONOMIAL
order = CONSTANT
fv = true
[]
[]
[Problem]
coord_type = 'RZ'
[]
[FVKernels]
[diff_v]
type = FVDiffusion
variable = v
coeff = coeff
[]
[body_v]
type = FVBodyForce
variable = v
function = 'forcing'
[]
[]
[FVBCs]
[boundary]
type = FVFunctionDirichletBC
boundary = 'left right top bottom'
function = 'exact'
variable = v
[]
[]
[Materials]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '1'
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = '1.1*sin(0.9*x)*cos(1.2*y)'
[]
[forcing]
type = ParsedFunction
value = '1.584*sin(0.9*x)*cos(1.2*y) - (-0.891*x*sin(0.9*x)*cos(1.2*y) + 0.99*cos(0.9*x)*cos(1.2*y))/x'
[]
[]
[Executioner]
type = Steady
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[./error]
type = ElementL2Error
variable = v
function = exact
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/fvkernels/mms/mass-mom-mat-advection-diffusion/input.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 1
nx = 2
xmin = -.6
xmax = .6
[]
[]
[GlobalParams]
advected_interp_method = 'average'
[]
[Variables]
[fv_rho]
order = CONSTANT
family = MONOMIAL
fv = true
initial_condition = 2
[]
[fv_vel]
order = CONSTANT
family = MONOMIAL
fv = true
initial_condition = 2
[]
[]
[FVKernels]
[adv_rho]
type = FVMatAdvection
variable = fv_rho
vel = 'fv_velocity'
[]
[diff_rho]
type = FVDiffusion
variable = fv_rho
coeff = coeff
[]
[forcing_rho]
type = FVBodyForce
variable = fv_rho
function = 'forcing_rho'
[]
[adv_rho_u]
type = FVMatAdvection
variable = fv_vel
vel = 'fv_velocity'
advected_quantity = 'rho_u'
[]
[diff_vel]
type = FVDiffusion
variable = fv_vel
coeff = coeff
[]
[forcing_vel]
type = FVBodyForce
variable = fv_vel
function = 'forcing_vel'
[]
[]
[FVBCs]
[boundary_rho]
type = FVFunctionDirichletBC
boundary = 'left right'
function = 'exact_rho'
variable = fv_rho
[]
[boundary_vel]
type = FVFunctionDirichletBC
boundary = 'left right'
function = 'exact_vel'
variable = fv_vel
[]
[]
[Materials]
[euler_material]
type = ADCoupledVelocityMaterial
vel_x = fv_vel
rho = fv_rho
velocity = 'fv_velocity'
[]
[diff]
type = ADGenericFunctorMaterial
prop_names = 'coeff'
prop_values = '1'
[]
[]
[Executioner]
type = Steady
solve_type = NEWTON
line_search = 'none'
[]
[Outputs]
csv = true
exodus = true
[]
[Functions]
[forcing_rho]
type = ParsedFunction
value = '-1.331*sin(1.1*x)^2 + 1.331*sin(1.1*x) + 1.331*cos(1.1*x)^2'
[]
[exact_rho]
type = ParsedFunction
value = '1.1*sin(1.1*x)'
[]
[forcing_vel]
type = ParsedFunction
value = '-2.9282*sin(1.1*x)^2*cos(1.1*x) + 1.4641*cos(1.1*x)^3 + 1.331*cos(1.1*x)'
[]
[exact_vel]
type = ParsedFunction
value = '1.1*cos(1.1*x)'
[]
[]
[Postprocessors]
[./l2_rho]
type = ElementL2Error
variable = fv_rho
function = exact_rho
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./l2_vel]
type = ElementL2Error
variable = fv_vel
function = exact_vel
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_volume/ins/mms/channel-flow/2d-average-with-temp.i)
mu=1.1
rho=1.1
k=1.1
cp=1.1
advected_interp_method='average'
velocity_interp_method='average'
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 2
ymin = -1
ymax = 1
nx = 2
ny = 2
[]
[]
[GlobalParams]
rhie_chow_user_object = 'rc'
[]
[UserObjects]
[rc]
type = INSFVRhieChowInterpolator
u = u
v = v
pressure = pressure
[]
[]
[Problem]
fv_bcs_integrity_check = true
[]
[Variables]
[u]
type = INSFVVelocityVariable
initial_condition = 1
two_term_boundary_expansion = false
[]
[v]
type = INSFVVelocityVariable
initial_condition = 1
two_term_boundary_expansion = false
[]
[pressure]
type = INSFVPressureVariable
two_term_boundary_expansion = false
[]
[temperature]
type = INSFVEnergyVariable
two_term_boundary_expansion = false
[]
[]
[FVKernels]
[mass]
type = INSFVMassAdvection
variable = pressure
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
rho = ${rho}
[]
[mass_forcing]
type = FVBodyForce
variable = pressure
function = forcing_p
[]
[u_advection]
type = INSFVMomentumAdvection
variable = u
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_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
[]
[u_forcing]
type = INSFVBodyForce
variable = u
functor = forcing_u
momentum_component = 'x'
[]
[v_advection]
type = INSFVMomentumAdvection
variable = v
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_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
[]
[v_forcing]
type = INSFVBodyForce
variable = v
functor = forcing_v
momentum_component = 'y'
[]
[temp_conduction]
type = FVDiffusion
coeff = 'k'
variable = temperature
[]
[temp_advection]
type = INSFVEnergyAdvection
variable = temperature
advected_interp_method = ${advected_interp_method}
velocity_interp_method = ${velocity_interp_method}
[]
[temp_forcing]
type = FVBodyForce
variable = temperature
function = forcing_t
[]
[]
[FVBCs]
[inlet-u]
type = INSFVInletVelocityBC
boundary = 'left'
variable = u
function = 'exact_u'
[]
[inlet-v]
type = INSFVInletVelocityBC
boundary = 'left'
variable = v
function = 'exact_v'
[]
[walls-u]
type = INSFVNoSlipWallBC
boundary = 'top bottom'
variable = u
function = 'exact_u'
[]
[walls-v]
type = INSFVNoSlipWallBC
boundary = 'top bottom'
variable = v
function = 'exact_v'
[]
[inlet-and-walls-t]
type = FVFunctionDirichletBC
boundary = 'left top bottom'
variable = temperature
function = 'exact_t'
[]
[outlet_p]
type = INSFVOutletPressureBC
boundary = 'right'
variable = pressure
function = 'exact_p'
[]
[]
[Materials]
[const]
type = ADGenericFunctorMaterial
prop_names = 'k cp'
prop_values = '${k} ${cp}'
[]
[ins_fv]
type = INSFVEnthalpyMaterial
temperature = 'temperature'
rho = ${rho}
[]
[]
[Functions]
[exact_u]
type = ParsedFunction
value = 'sin((1/2)*y*pi)*cos((1/2)*x*pi)'
[]
[exact_rhou]
type = ParsedFunction
value = 'rho*sin((1/2)*y*pi)*cos((1/2)*x*pi)'
vars = 'rho'
vals = '${rho}'
[]
[forcing_u]
type = ADParsedFunction
value = '(1/2)*pi^2*mu*sin((1/2)*y*pi)*cos((1/2)*x*pi) - 1/2*pi*rho*sin((1/4)*x*pi)*sin((1/2)*y*pi)^2*cos((1/2)*x*pi) + (1/2)*pi*rho*sin((1/4)*x*pi)*cos((1/2)*x*pi)*cos((1/2)*y*pi)^2 - pi*rho*sin((1/2)*x*pi)*sin((1/2)*y*pi)^2*cos((1/2)*x*pi) - 1/4*pi*sin((1/4)*x*pi)*sin((3/2)*y*pi)'
vars = 'mu rho'
vals = '${mu} ${rho}'
[]
[exact_v]
type = ParsedFunction
value = 'sin((1/4)*x*pi)*cos((1/2)*y*pi)'
[]
[exact_rhov]
type = ParsedFunction
value = 'rho*sin((1/4)*x*pi)*cos((1/2)*y*pi)'
vars = 'rho'
vals = '${rho}'
[]
[forcing_v]
type = ADParsedFunction
value = '(5/16)*pi^2*mu*sin((1/4)*x*pi)*cos((1/2)*y*pi) - pi*rho*sin((1/4)*x*pi)^2*sin((1/2)*y*pi)*cos((1/2)*y*pi) - 1/2*pi*rho*sin((1/4)*x*pi)*sin((1/2)*x*pi)*sin((1/2)*y*pi)*cos((1/2)*y*pi) + (1/4)*pi*rho*sin((1/2)*y*pi)*cos((1/4)*x*pi)*cos((1/2)*x*pi)*cos((1/2)*y*pi) + (3/2)*pi*cos((1/4)*x*pi)*cos((3/2)*y*pi)'
vars = 'mu rho'
vals = '${mu} ${rho}'
[]
[exact_p]
type = ParsedFunction
value = 'sin((3/2)*y*pi)*cos((1/4)*x*pi)'
[]
[forcing_p]
type = ParsedFunction
value = '-1/2*pi*rho*sin((1/4)*x*pi)*sin((1/2)*y*pi) - 1/2*pi*rho*sin((1/2)*x*pi)*sin((1/2)*y*pi)'
vars = 'rho'
vals = '${rho}'
[]
[exact_t]
type = ParsedFunction
value = 'sin((1/4)*x*pi)*cos((1/2)*y*pi)'
[]
[forcing_t]
type = ParsedFunction
value = '-pi*cp*rho*sin((1/4)*x*pi)^2*sin((1/2)*y*pi)*cos((1/2)*y*pi) - 1/2*pi*cp*rho*sin((1/4)*x*pi)*sin((1/2)*x*pi)*sin((1/2)*y*pi)*cos((1/2)*y*pi) + (1/4)*pi*cp*rho*sin((1/2)*y*pi)*cos((1/4)*x*pi)*cos((1/2)*x*pi)*cos((1/2)*y*pi) + (5/16)*pi^2*k*sin((1/4)*x*pi)*cos((1/2)*y*pi)'
vars = 'k rho cp'
vals = '${k} ${rho} ${cp}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
petsc_options_value = 'asm 100 lu NONZERO'
line_search = 'none'
[]
[Outputs]
exodus = true
csv = true
[]
[Postprocessors]
[h]
type = AverageElementSize
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[./L2u]
type = ElementL2Error
variable = u
function = exact_u
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2v]
type = ElementL2Error
variable = v
function = exact_v
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2p]
variable = pressure
function = exact_p
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[./L2t]
variable = temperature
function = exact_t
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[../]
[]
(modules/navier_stokes/include/fvbcs/INSFVInletVelocityBC.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 "FVFunctionDirichletBC.h"
#include "INSFVFlowBC.h"
/**
* A class for velocity inlet boundary conditions
*/
class INSFVInletVelocityBC : public FVFunctionDirichletBC, public INSFVFlowBC
{
public:
static InputParameters validParams();
INSFVInletVelocityBC(const InputParameters & params);
};
(modules/navier_stokes/include/fvbcs/INSFVNoSlipWallBC.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 "FVFunctionDirichletBC.h"
class Function;
/**
* A class for no slip velocity boundary condtions
*/
class INSFVNoSlipWallBC : public FVFunctionDirichletBC
{
public:
static InputParameters validParams();
INSFVNoSlipWallBC(const InputParameters & params);
};