- ux-component
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:x-component
- vector_prop_nameThe name to give the declared vector material property
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:The name to give the declared vector material property
VectorFromComponentVariablesMaterial
Description
VectorFromComponentVariablesMaterial
computes a vector material property with name specified by the vector_prop_name
parameter from coupled variable components. The x-component is computed through the u
coupled variable. Optional coupled variables v
and w
compute the y- and z-components respectively.
An example of this object's use is shown in the listing below where in this case a velocity material property is being declared. The ability to pass constants to the coupled variables is leveraged in this example. Actual coupled variable instances would be used in, for example, a Navier-Stokes simulation in which the nonlinear system solves for velocity components and the vector velocity needs to be constructed.
[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
[vel]
type = ADVectorFromComponentVariablesMaterial<<<{"description": "Computes a vector material property from coupled variables", "href": "VectorFromComponentVariablesMaterial.html"}>>>
vector_prop_name<<<{"description": "The name to give the declared vector material property"}>>> = 'velocity'
u<<<{"description": "x-component"}>>> = 1
v<<<{"description": "y-component"}>>> = 0
[]
[]
(test/tests/dgkernels/passive-scalar-channel-flow/test.i)Input Parameters
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- boundaryThe list of boundaries (ids or names) from the mesh where this object applies
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:The list of boundaries (ids or names) from the mesh where this object applies
- computeTrueWhen false, MOOSE will not call compute methods on this material. The user must call computeProperties() after retrieving the MaterialBase via MaterialBasePropertyInterface::getMaterialBase(). Non-computed MaterialBases are not sorted for dependencies.
Default:True
C++ Type:bool
Controllable:No
Description:When false, MOOSE will not call compute methods on this material. The user must call computeProperties() after retrieving the MaterialBase via MaterialBasePropertyInterface::getMaterialBase(). Non-computed MaterialBases are not sorted for dependencies.
- constant_onNONEWhen ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped
Default:NONE
C++ Type:MooseEnum
Options:NONE, ELEMENT, SUBDOMAIN
Controllable:No
Description:When ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped
- declare_suffixAn optional suffix parameter that can be appended to any declared properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:An optional suffix parameter that can be appended to any declared properties. The suffix will be prepended with a '_' character.
- vy-component
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:y-component
- wz-component
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:z-component
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
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Controllable:No
Description:The seed for the master random number generator
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
- output_propertiesList of material properties, from this material, to output (outputs must also be defined to an output type)
C++ Type:std::vector<std::string>
Controllable:No
Description:List of material properties, from this material, to output (outputs must also be defined to an output type)
- outputsnone Vector of output names where you would like to restrict the output of variables(s) associated with this object
Default:none
C++ Type:std::vector<OutputName>
Controllable:No
Description:Vector of output names where you would like to restrict the output of variables(s) associated with this object
Outputs Parameters
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
- use_interpolated_stateFalseFor the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Default:False
C++ Type:bool
Controllable:No
Description:For the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Material Property Retrieval Parameters
Input Files
- (modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/lid-driven-fsp.i)
- (modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/mms-lid-driven.i)
- (modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/stokes-lid-driven.i)
- (modules/navier_stokes/test/tests/finite_element/ins/cg-dg-hybrid/channel-flow/channel-hybrid.i)
- (modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/lid-driven.i)
- (modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/sc-lid-driven.i)
- (modules/navier_stokes/test/tests/finite_element/ins/cg-dg-hybrid/lid-driven/hybrid-cg-dg.i)
- (modules/navier_stokes/test/tests/finite_element/ins/cg-dg-hybrid/mms/lid-driven-skewed/hybrid-skewed-vortex.i)
- (modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/mms-second-lid-driven.i)
- (modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/unstabilized-velocity-component-objects.i)
- (modules/navier_stokes/test/tests/finite_element/ins/cg-dg-hybrid/mms/lid-driven/hybrid-cg-dg-mms.i)
- (test/tests/dgkernels/passive-scalar-channel-flow/test.i)
- (modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/matrix-testing.i)
(test/tests/dgkernels/passive-scalar-channel-flow/test.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = -1
ymax = 1
nx = 20
ny = 4
[]
[]
[Variables]
[u]
family = MONOMIAL
[]
[]
[Kernels]
[convection]
type = ADConservativeAdvection
variable = u
velocity = 'velocity'
[]
[diffusion]
type = MatDiffusion
variable = u
diffusivity = 1
[]
[]
[DGKernels]
[convection]
type = ADDGAdvection
variable = u
velocity = 'velocity'
[]
[diffusion]
type = DGDiffusion
variable = u
sigma = 6
epsilon = -1
diff = 1
[]
[]
[Functions]
[v_inlet]
type = ParsedVectorFunction
expression_x = '1'
[]
[]
[BCs]
[u_walls]
type = DGFunctionDiffusionDirichletBC
boundary = 'bottom top'
variable = u
sigma = 6
epsilon = -1
function = '0'
diff = 1
[]
[u_in]
type = ADConservativeAdvectionBC
boundary = 'left'
variable = u
velocity_function = v_inlet
primal_dirichlet_value = 1
[]
[u_out]
type = ADConservativeAdvectionBC
boundary = 'right'
variable = u
velocity_mat_prop = 'velocity'
[]
[]
[Materials]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = 1
v = 0
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/lid-driven-fsp.i)
mu = 1
rho = 1
l = 1
U = 1e3
n = 200
gamma = 1e5
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = ${l}
nx = ${n}
ny = ${n}
elem_type = TRI6
[]
[]
[Problem]
type = NavierStokesProblem
extra_tag_matrices = 'mass'
mass_matrix = 'mass'
use_pressure_mass_matrix = true
[]
[Variables]
[vel_x]
family = MONOMIAL
order = FIRST
[]
[vel_y]
family = MONOMIAL
order = FIRST
[]
[pressure]
family = MONOMIAL
order = CONSTANT
[]
[vel_bar_x]
family = SIDE_HIERARCHIC
order = FIRST
[]
[vel_bar_y]
family = SIDE_HIERARCHIC
order = FIRST
[]
[pressure_bar]
family = SIDE_HIERARCHIC
order = FIRST
[]
[]
[HDGKernels]
[momentum_x_convection]
type = AdvectionIPHDGKernel
variable = vel_x
face_variable = vel_bar_x
velocity = 'velocity'
coeff = ${rho}
[]
[momentum_x_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_x
face_variable = vel_bar_x
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 0
[]
[momentum_y_convection]
type = AdvectionIPHDGKernel
variable = vel_y
face_variable = vel_bar_y
velocity = 'velocity'
coeff = ${rho}
[]
[momentum_y_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_y
face_variable = vel_bar_y
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 1
[]
[pressure_convection]
type = AdvectionIPHDGKernel
variable = pressure
face_variable = pressure_bar
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
[]
[]
[Kernels]
[mass_matrix_pressure]
type = MassMatrix
variable = pressure
matrix_tags = 'mass'
density = '${fparse -1/gamma}'
[]
[grad_div_x]
type = GradDiv
variable = vel_x
u = vel_x
v = vel_y
gamma = ${gamma}
component = 0
[]
[grad_div_y]
type = GradDiv
variable = vel_y
u = vel_x
v = vel_y
gamma = ${gamma}
component = 1
[]
[]
[DGKernels]
[pb_mass]
type = MassMatrixDGKernel
variable = pressure_bar
matrix_tags = 'mass'
density = '${fparse -1/gamma}'
[]
[u_jump]
type = MassFluxPenalty
variable = vel_x
u = vel_x
v = vel_y
component = 0
gamma = ${gamma}
[]
[v_jump]
type = MassFluxPenalty
variable = vel_y
u = vel_x
v = vel_y
component = 1
gamma = ${gamma}
[]
[]
[BCs]
[momentum_x_diffusion_walls]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right'
variable = vel_x
face_variable = vel_bar_x
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '0'
diffusivity = 'mu'
component = 0
[]
[momentum_x_diffusion_top]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'top'
variable = vel_x
face_variable = vel_bar_x
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '${U}'
diffusivity = 'mu'
component = 0
[]
[momentum_y_diffusion_all]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right top'
variable = vel_y
face_variable = vel_bar_y
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '0'
diffusivity = 'mu'
component = 1
[]
[mass_convection]
type = AdvectionIPHDGPrescribedFluxBC
face_variable = pressure_bar
variable = pressure
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
boundary = 'left bottom top right'
prescribed_normal_flux = 0
[]
[pb_mass]
type = MassMatrixIntegratedBC
variable = pressure_bar
matrix_tags = 'mass'
boundary = 'left right bottom top'
density = '${fparse -1/gamma}'
[]
[u_jump_walls]
type = MassFluxPenaltyBC
variable = vel_x
u = vel_x
v = vel_y
component = 0
boundary = 'left right bottom'
gamma = ${gamma}
dirichlet_value = 'walls'
[]
[v_jump_walls]
type = MassFluxPenaltyBC
variable = vel_y
u = vel_x
v = vel_y
component = 1
boundary = 'left right bottom'
gamma = ${gamma}
dirichlet_value = 'walls'
[]
[u_jump_top]
type = MassFluxPenaltyBC
variable = vel_x
u = vel_x
v = vel_y
component = 0
boundary = 'top'
gamma = ${gamma}
dirichlet_value = 'top'
[]
[v_jump_top]
type = MassFluxPenaltyBC
variable = vel_y
u = vel_x
v = vel_y
component = 1
boundary = 'top'
gamma = ${gamma}
dirichlet_value = 'top'
[]
[]
[Functions]
[top]
type = ParsedVectorFunction
value_x = ${U}
value_y = 0
[]
[walls]
type = ParsedVectorFunction
value_x = 0
value_y = 0
[]
[]
[Materials]
[const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${rho} ${mu}'
[]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = vel_x
v = vel_y
[]
[]
[Preconditioning]
[FSP]
type = FSP
topsplit = 'up'
[up]
splitting = 'u p'
splitting_type = schur
petsc_options_iname = '-pc_fieldsplit_schur_fact_type -pc_fieldsplit_schur_precondition -ksp_gmres_restart -ksp_type -ksp_pc_side -ksp_rtol'
petsc_options_value = 'full self 300 fgmres right 1e-4'
[]
[u]
vars = 'vel_x vel_y vel_bar_x vel_bar_y'
petsc_options = '-ksp_converged_reason'
petsc_options_iname = '-pc_type -ksp_type -ksp_rtol -ksp_gmres_restart -ksp_pc_side -pc_factor_mat_solver_type'
petsc_options_value = 'ilu gmres 1e-2 300 right strumpack'
[]
[p]
vars = 'pressure pressure_bar'
petsc_options = '-ksp_converged_reason'
petsc_options_iname = '-ksp_type -ksp_gmres_restart -ksp_rtol -pc_type -ksp_pc_side'
petsc_options_value = 'gmres 300 1e-2 ilu right'
[]
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
active = ''
[out]
type = Exodus
hide = 'pressure_average vel_bar_x vel_bar_y pressure_bar'
[]
[]
[Postprocessors]
[Re]
type = ParsedPostprocessor
pp_names = ''
expression = '${rho} * ${U} * ${l} / ${mu}'
[]
[pressure_average]
type = ElementAverageValue
variable = pressure
[]
[]
[Correctors]
[set_pressure]
type = NSPressurePin
pin_type = 'average'
variable = pressure
pressure_average = 'pressure_average'
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/mms-lid-driven.i)
mu = 2
rho = 2
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = TRI6
[]
[]
[Variables]
[vel_x]
family = MONOMIAL
order = FIRST
[]
[vel_y]
family = MONOMIAL
order = FIRST
[]
[pressure]
family = MONOMIAL
order = CONSTANT
[]
[vel_bar_x]
family = SIDE_HIERARCHIC
order = FIRST
[]
[vel_bar_y]
family = SIDE_HIERARCHIC
order = FIRST
[]
[pressure_bar]
family = SIDE_HIERARCHIC
order = FIRST
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[HDGKernels]
[momentum_x_convection]
type = AdvectionIPHDGKernel
variable = vel_x
face_variable = vel_bar_x
velocity = 'velocity'
coeff = ${rho}
[]
[momentum_x_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_x
face_variable = vel_bar_x
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 0
[]
[momentum_y_convection]
type = AdvectionIPHDGKernel
variable = vel_y
face_variable = vel_bar_y
velocity = 'velocity'
coeff = ${rho}
[]
[momentum_y_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_y
face_variable = vel_bar_y
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 1
[]
[pressure_convection]
type = AdvectionIPHDGKernel
variable = pressure
face_variable = pressure_bar
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
[]
[]
[Kernels]
[momentum_x_ffn]
type = BodyForce
variable = vel_x
function = forcing_u
[]
[momentum_y_ffn]
type = BodyForce
variable = vel_y
function = forcing_v
[]
[mass_ffn]
type = BodyForce
variable = pressure
function = forcing_p
[]
[mean_zero_pressure]
type = ScalarLagrangeMultiplier
variable = pressure
lambda = lambda
[]
[]
[ScalarKernels]
[mean_zero_pressure_lm]
type = AverageValueConstraint
variable = lambda
pp_name = pressure_integral
value = 0
[]
[]
[BCs]
[momentum_x_diffusion_all]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right top'
variable = vel_x
face_variable = vel_bar_x
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = exact_u
diffusivity = 'mu'
component = 0
[]
[momentum_y_diffusion_all]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right top'
variable = vel_y
face_variable = vel_bar_y
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = exact_v
diffusivity = 'mu'
component = 1
[]
[mass_convection]
type = AdvectionIPHDGPrescribedFluxBC
face_variable = pressure_bar
variable = pressure
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
boundary = 'left bottom top right'
prescribed_normal_flux = 0
[]
[]
[Functions]
[exact_u]
type = ParsedFunction
expression = 'sin(y)*cos((1/2)*x*pi)'
[]
[forcing_u]
type = ParsedFunction
expression = 'mu*sin(y)*cos((1/2)*x*pi) + (1/4)*pi^2*mu*sin(y)*cos((1/2)*x*pi) - 1/2*pi*rho*sin(x)*sin(y)*sin((1/2)*y*pi)*cos((1/2)*x*pi) + rho*sin(x)*cos(y)*cos((1/2)*x*pi)*cos((1/2)*y*pi) - pi*rho*sin(y)^2*sin((1/2)*x*pi)*cos((1/2)*x*pi) + sin(y)*cos(x)'
symbol_names = 'mu rho'
symbol_values = '${mu} ${rho}'
[]
[exact_v]
type = ParsedFunction
expression = 'sin(x)*cos((1/2)*y*pi)'
[]
[forcing_v]
type = ParsedFunction
expression = 'mu*sin(x)*cos((1/2)*y*pi) + (1/4)*pi^2*mu*sin(x)*cos((1/2)*y*pi) - pi*rho*sin(x)^2*sin((1/2)*y*pi)*cos((1/2)*y*pi) - 1/2*pi*rho*sin(x)*sin(y)*sin((1/2)*x*pi)*cos((1/2)*y*pi) + rho*sin(y)*cos(x)*cos((1/2)*x*pi)*cos((1/2)*y*pi) + sin(x)*cos(y)'
symbol_names = 'mu rho'
symbol_values = '${mu} ${rho}'
[]
[exact_p]
type = ParsedFunction
expression = 'sin(x)*sin(y)'
[]
[forcing_p]
type = ParsedFunction
expression = '(1/2)*pi*rho*sin(x)*sin((1/2)*y*pi) + (1/2)*pi*rho*sin(y)*sin((1/2)*x*pi)'
symbol_names = 'rho'
symbol_values = '${rho}'
[]
[]
[Materials]
[const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${rho} ${mu}'
[]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = vel_x
v = vel_y
[]
[]
[AuxVariables]
[vel_exact_x]
[]
[vel_exact_y]
[]
[pressure_exact]
[]
[]
[AuxKernels]
[vel_exact_x]
type = FunctionAux
variable = vel_exact_x
function = exact_u
[]
[vel_exact_y]
type = FunctionAux
variable = vel_exact_y
function = exact_v
[]
[pressure_exact]
type = FunctionAux
variable = pressure_exact
function = exact_p
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_rel_tol = 1e-12
[]
[Outputs]
[out]
type = Exodus
hide = 'lambda pressure_integral'
[]
csv = true
[]
[Postprocessors]
[pressure_integral]
type = ElementIntegralVariablePostprocessor
variable = pressure
execute_on = linear
[]
[h]
type = AverageElementSize
execute_on = 'timestep_end'
[]
[L2u]
type = ElementL2Error
variable = vel_x
function = exact_u
execute_on = 'timestep_end'
[]
[L2v]
type = ElementL2Error
variable = vel_y
function = exact_v
execute_on = 'timestep_end'
[]
[L2p]
type = ElementL2Error
variable = pressure
function = exact_p
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/stokes-lid-driven.i)
mu = 1
rho = 1
l = 1
U = 1
n = 4
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = ${l}
nx = ${n}
ny = ${n}
elem_type = QUAD9
[]
[]
[Variables]
[vel_x]
family = L2_HIERARCHIC
order = SECOND
[]
[vel_y]
family = L2_HIERARCHIC
order = SECOND
[]
[pressure]
family = L2_HIERARCHIC
order = FIRST
[]
[vel_bar_x]
family = SIDE_HIERARCHIC
order = SECOND
[]
[vel_bar_y]
family = SIDE_HIERARCHIC
order = SECOND
[]
[pressure_bar]
family = SIDE_HIERARCHIC
order = FIRST
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[HDGKernels]
[momentum_x_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_x
face_variable = vel_bar_x
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 0
[]
[momentum_y_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_y
face_variable = vel_bar_y
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 1
[]
[pressure_convection]
type = AdvectionIPHDGKernel
variable = pressure
face_variable = pressure_bar
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
[]
[]
[Kernels]
[mean_zero_pressure]
type = ScalarLagrangeMultiplier
variable = pressure
lambda = lambda
[]
[]
[ScalarKernels]
[mean_zero_pressure_lm]
type = AverageValueConstraint
variable = lambda
pp_name = pressure_integral
value = 0
[]
[]
[BCs]
[momentum_x_diffusion_walls]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right'
variable = vel_x
face_variable = vel_bar_x
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '0'
diffusivity = 'mu'
component = 0
[]
[momentum_x_diffusion_top]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'top'
variable = vel_x
face_variable = vel_bar_x
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '${U}'
diffusivity = 'mu'
component = 0
[]
[momentum_y_diffusion_all]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right top'
variable = vel_y
face_variable = vel_bar_y
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '0'
diffusivity = 'mu'
component = 1
[]
[mass_convection]
type = AdvectionIPHDGPrescribedFluxBC
face_variable = pressure_bar
variable = pressure
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
boundary = 'left bottom top right'
prescribed_normal_flux = 0
[]
[]
[Materials]
[const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${rho} ${mu}'
[]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = vel_x
v = vel_y
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_rel_tol = 1e-12
[]
[Outputs]
[out]
type = CSV
show = 'symmetric'
execute_on = 'timestep_end'
[]
[]
[Postprocessors]
[symmetric]
type = MatrixSymmetryCheck
execute_on = 'timestep_end'
[]
[pressure_integral]
type = ElementIntegralVariablePostprocessor
variable = pressure
execute_on = linear
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/cg-dg-hybrid/channel-flow/channel-hybrid.i)
mu = 1.1
rho = 1.1
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = -1
ymax = 1
nx = 100
ny = 20
[]
[]
[Variables]
[u]
family = MONOMIAL
[]
[v]
family = MONOMIAL
[]
[pressure][]
[]
[Kernels]
[momentum_x_convection]
type = ADConservativeAdvection
variable = u
velocity = 'velocity'
advected_quantity = 'rhou'
[]
[momentum_x_diffusion]
type = MatDiffusion
variable = u
diffusivity = 'mu'
[]
[momentum_x_pressure]
type = PressureGradient
integrate_p_by_parts = false
variable = u
pressure = pressure
component = 0
[]
[momentum_y_convection]
type = ADConservativeAdvection
variable = v
velocity = 'velocity'
advected_quantity = 'rhov'
[]
[momentum_y_diffusion]
type = MatDiffusion
variable = v
diffusivity = 'mu'
[]
[momentum_y_pressure]
type = PressureGradient
integrate_p_by_parts = false
variable = v
pressure = pressure
component = 1
[]
[mass]
type = ADConservativeAdvection
variable = pressure
velocity = velocity
advected_quantity = -1
[]
[]
[DGKernels]
[momentum_x_convection]
type = ADDGAdvection
variable = u
velocity = 'velocity'
advected_quantity = 'rhou'
[]
[momentum_x_diffusion]
type = DGDiffusion
variable = u
sigma = 6
epsilon = -1
diff = 'mu'
[]
[momentum_y_convection]
type = ADDGAdvection
variable = v
velocity = 'velocity'
advected_quantity = 'rhov'
[]
[momentum_y_diffusion]
type = DGDiffusion
variable = v
sigma = 6
epsilon = -1
diff = 'mu'
[]
[]
[Functions]
[v_inlet]
type = ParsedVectorFunction
expression_x = '1'
[]
[]
[BCs]
[u_walls]
type = DGFunctionDiffusionDirichletBC
boundary = 'bottom top'
variable = u
sigma = 6
epsilon = -1
function = '0'
diff = 'mu'
[]
[v_walls]
type = DGFunctionDiffusionDirichletBC
boundary = 'bottom top'
variable = v
sigma = 6
epsilon = -1
function = '0'
diff = 'mu'
[]
[u_in]
type = ADConservativeAdvectionBC
boundary = 'left'
variable = u
velocity_function = v_inlet
primal_dirichlet_value = 1
primal_coefficient = 'rho'
[]
[v_in]
type = ADConservativeAdvectionBC
boundary = 'left'
variable = v
velocity_function = v_inlet
primal_dirichlet_value = 0
primal_coefficient = 'rho'
[]
[p_in]
type = ADConservativeAdvectionBC
boundary = 'left'
variable = pressure
velocity_function = v_inlet
advected_quantity = -1
[]
[u_out]
type = ADConservativeAdvectionBC
boundary = 'right'
variable = u
velocity_mat_prop = 'velocity'
advected_quantity = 'rhou'
[]
[v_out]
type = ADConservativeAdvectionBC
boundary = 'right'
variable = v
velocity_mat_prop = 'velocity'
advected_quantity = 'rhov'
[]
[p_out]
type = DirichletBC
variable = pressure
boundary = 'right'
value = 0
[]
[]
[Materials]
[const]
type = ADGenericConstantMaterial
prop_names = 'rho'
prop_values = '${rho}'
[]
[const_reg]
type = GenericConstantMaterial
prop_names = 'mu'
prop_values = '${mu}'
[]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = u
v = v
[]
[rhou]
type = ADParsedMaterial
property_name = 'rhou'
coupled_variables = 'u'
material_property_names = 'rho'
expression = 'rho*u'
[]
[rhov]
type = ADParsedMaterial
property_name = 'rhov'
coupled_variables = 'v'
material_property_names = 'rho'
expression = 'rho*v'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/lid-driven.i)
mu = 1
rho = 1
l = 1
U = 100
n = 8
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = ${l}
nx = ${n}
ny = ${n}
elem_type = TRI6
[]
[]
[Variables]
[vel_x]
family = L2_HIERARCHIC
order = SECOND
[]
[vel_y]
family = L2_HIERARCHIC
order = SECOND
[]
[pressure]
family = L2_HIERARCHIC
order = FIRST
[]
[vel_bar_x]
family = SIDE_HIERARCHIC
order = SECOND
[]
[vel_bar_y]
family = SIDE_HIERARCHIC
order = SECOND
[]
[pressure_bar]
family = SIDE_HIERARCHIC
order = FIRST
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[HDGKernels]
[momentum_x_convection]
type = AdvectionIPHDGKernel
variable = vel_x
face_variable = vel_bar_x
velocity = 'velocity'
coeff = ${rho}
[]
[momentum_x_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_x
face_variable = vel_bar_x
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 0
[]
[momentum_y_convection]
type = AdvectionIPHDGKernel
variable = vel_y
face_variable = vel_bar_y
velocity = 'velocity'
coeff = ${rho}
[]
[momentum_y_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_y
face_variable = vel_bar_y
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 1
[]
[pressure_convection]
type = AdvectionIPHDGKernel
variable = pressure
face_variable = pressure_bar
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
[]
[]
[Kernels]
[mean_zero_pressure]
type = ScalarLagrangeMultiplier
variable = pressure
lambda = lambda
[]
[]
[ScalarKernels]
[mean_zero_pressure_lm]
type = AverageValueConstraint
variable = lambda
pp_name = pressure_integral
value = 0
[]
[]
[BCs]
[momentum_x_diffusion_walls]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right'
variable = vel_x
face_variable = vel_bar_x
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '0'
diffusivity = 'mu'
component = 0
[]
[momentum_x_diffusion_top]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'top'
variable = vel_x
face_variable = vel_bar_x
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '${U}'
diffusivity = 'mu'
component = 0
[]
[momentum_y_diffusion_all]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right top'
variable = vel_y
face_variable = vel_bar_y
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '0'
diffusivity = 'mu'
component = 1
[]
[mass_convection]
type = AdvectionIPHDGPrescribedFluxBC
face_variable = pressure_bar
variable = pressure
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
boundary = 'left bottom top right'
prescribed_normal_flux = 0
[]
[]
[Materials]
[const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${rho} ${mu}'
[]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = vel_x
v = vel_y
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_rel_tol = 1e-12
[]
[Outputs]
[out]
type = Exodus
hide = 'lambda pressure_integral vel_bar_x vel_bar_y pressure_bar'
[]
[]
[Postprocessors]
[Re]
type = ParsedPostprocessor
pp_names = ''
expression = '${rho} * ${U} * ${l} / ${mu}'
[]
[pressure_integral]
type = ElementIntegralVariablePostprocessor
variable = pressure
execute_on = linear
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/sc-lid-driven.i)
mu = 1
rho = 1
l = 1
U = 100
n = 8
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = ${l}
nx = ${n}
ny = ${n}
elem_type = TRI6
[]
[]
[Variables]
[vel_x]
family = L2_HIERARCHIC
order = SECOND
[]
[vel_y]
family = L2_HIERARCHIC
order = SECOND
[]
[pressure]
family = L2_HIERARCHIC
order = FIRST
[]
[vel_bar_x]
family = SIDE_HIERARCHIC
order = SECOND
[]
[vel_bar_y]
family = SIDE_HIERARCHIC
order = SECOND
[]
[pressure_bar]
family = SIDE_HIERARCHIC
order = FIRST
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[HDGKernels]
[momentum_x_advection]
type = AdvectionIPHDGKernel
variable = vel_x
face_variable = vel_bar_x
velocity = 'velocity'
coeff = ${rho}
[]
[momentum_x_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_x
face_variable = vel_bar_x
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 0
[]
[momentum_y_advection]
type = AdvectionIPHDGKernel
variable = vel_y
face_variable = vel_bar_y
velocity = 'velocity'
coeff = ${rho}
[]
[momentum_y_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_y
face_variable = vel_bar_y
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 1
[]
[mass_advection]
type = AdvectionIPHDGKernel
variable = pressure
face_variable = pressure_bar
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
[]
[]
[Kernels]
[mean_zero_pressure]
type = ScalarLagrangeMultiplier
variable = pressure
lambda = lambda
[]
[]
[ScalarKernels]
[mean_zero_pressure_lm]
type = AverageValueConstraint
variable = lambda
pp_name = pressure_integral
value = 0
[]
[]
[BCs]
[momentum_x_diffusion_walls]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right'
variable = vel_x
face_variable = vel_bar_x
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '0'
diffusivity = 'mu'
component = 0
[]
[momentum_x_diffusion_top]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'top'
variable = vel_x
face_variable = vel_bar_x
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '${U}'
diffusivity = 'mu'
component = 0
[]
[momentum_y_diffusion_all]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right top'
variable = vel_y
face_variable = vel_bar_y
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '0'
diffusivity = 'mu'
component = 1
[]
[mass_advection]
type = AdvectionIPHDGPrescribedFluxBC
face_variable = pressure_bar
variable = pressure
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
boundary = 'left bottom top right'
prescribed_normal_flux = 0
[]
[]
[Materials]
[const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${rho} ${mu}'
[]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = vel_x
v = vel_y
[]
[]
[Preconditioning]
[sc]
type = StaticCondensation
petsc_options_iname = '-pc_type -pc_factor_shift_type -ksp_view_pmat'
petsc_options_value = 'lu NONZERO binary'
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-12
petsc_options_iname = '-ksp_type'
petsc_options_value = 'preonly'
[]
[Outputs]
[out]
type = Exodus
hide = 'lambda pressure_integral symmetric vel_bar_x vel_bar_y pressure_bar'
[]
[csv]
type = CSV
hide = 'lambda pressure_integral'
[]
[]
[Postprocessors]
[Re]
type = ParsedPostprocessor
pp_names = ''
expression = '${rho} * ${U} * ${l} / ${mu}'
[]
[symmetric]
type = MatrixSymmetryCheck
mat = binaryoutput
mat_number_to_load = 2
[]
[pressure_integral]
type = ElementIntegralVariablePostprocessor
variable = pressure
execute_on = linear
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/cg-dg-hybrid/lid-driven/hybrid-cg-dg.i)
mu = 1
rho = 1
l = 200
U = 1
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = ${l}
nx = 20
ny = 20
[]
[corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[]
[]
[Variables]
[u]
family = MONOMIAL
[]
[v]
family = MONOMIAL
[]
[pressure][]
[]
[Kernels]
[momentum_x_convection]
type = ADConservativeAdvection
variable = u
velocity = 'velocity'
advected_quantity = 'rhou'
[]
[momentum_x_diffusion]
type = MatDiffusion
variable = u
diffusivity = 'mu'
[]
[momentum_x_pressure]
type = PressureGradient
integrate_p_by_parts = false
variable = u
pressure = pressure
component = 0
[]
[momentum_y_convection]
type = ADConservativeAdvection
variable = v
velocity = 'velocity'
advected_quantity = 'rhov'
[]
[momentum_y_diffusion]
type = MatDiffusion
variable = v
diffusivity = 'mu'
[]
[momentum_y_pressure]
type = PressureGradient
integrate_p_by_parts = false
variable = v
pressure = pressure
component = 1
[]
[mass]
type = ADConservativeAdvection
variable = pressure
velocity = velocity
advected_quantity = -1
[]
[]
[DGKernels]
[momentum_x_convection]
type = ADDGAdvection
variable = u
velocity = 'velocity'
advected_quantity = 'rhou'
[]
[momentum_x_diffusion]
type = DGDiffusion
variable = u
sigma = 6
epsilon = -1
diff = 'mu'
[]
[momentum_y_convection]
type = ADDGAdvection
variable = v
velocity = 'velocity'
advected_quantity = 'rhov'
[]
[momentum_y_diffusion]
type = DGDiffusion
variable = v
sigma = 6
epsilon = -1
diff = 'mu'
[]
[]
[BCs]
[u_walls]
type = DGFunctionDiffusionDirichletBC
boundary = 'left bottom right'
variable = u
sigma = 6
epsilon = -1
function = '0'
diff = 'mu'
[]
[v_walls]
type = DGFunctionDiffusionDirichletBC
boundary = 'left bottom right top'
variable = v
sigma = 6
epsilon = -1
function = '0'
diff = 'mu'
[]
[u_top]
type = DGFunctionDiffusionDirichletBC
boundary = 'top'
variable = u
sigma = 6
epsilon = -1
function = '${U}'
diff = 'mu'
[]
[pressure_pin]
type = DirichletBC
variable = pressure
boundary = 'pinned_node'
value = 0
[]
[]
[Materials]
[const]
type = ADGenericConstantMaterial
prop_names = 'rho'
prop_values = '${rho}'
[]
[const_reg]
type = GenericConstantMaterial
prop_names = 'mu'
prop_values = '${mu}'
[]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = u
v = v
[]
[rhou]
type = ADParsedMaterial
property_name = 'rhou'
coupled_variables = 'u'
material_property_names = 'rho'
expression = 'rho*u'
[]
[rhov]
type = ADParsedMaterial
property_name = 'rhov'
coupled_variables = 'v'
material_property_names = 'rho'
expression = 'rho*v'
[]
[]
[AuxVariables]
[vel_x]
family = MONOMIAL
order = CONSTANT
[]
[vel_y]
family = MONOMIAL
order = CONSTANT
[]
[p]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[vel_x]
type = ProjectionAux
variable = vel_x
v = u
[]
[vel_y]
type = ProjectionAux
variable = vel_y
v = v
[]
[p]
type = ProjectionAux
variable = p
v = pressure
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
[Postprocessors]
[Re]
type = ParsedPostprocessor
expression = '${rho} * ${U} * ${l} / ${mu}'
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/cg-dg-hybrid/mms/lid-driven-skewed/hybrid-skewed-vortex.i)
rho=1
mu=1
[Mesh]
[gen_mesh]
type = FileMeshGenerator
file = skewed.msh
[]
[corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen_mesh
[]
[]
[Variables]
[u]
family = MONOMIAL
order = SECOND
[]
[v]
family = MONOMIAL
order = SECOND
[]
[pressure][]
[]
[Kernels]
[momentum_x_convection]
type = ADConservativeAdvection
variable = u
velocity = 'velocity'
[]
[momentum_x_diffusion]
type = Diffusion
variable = u
[]
[momentum_x_pressure]
type = PressureGradient
integrate_p_by_parts = false
variable = u
pressure = pressure
component = 0
[]
[u_forcing]
type = BodyForce
variable = u
function = forcing_u
[]
[momentum_y_convection]
type = ADConservativeAdvection
variable = v
velocity = 'velocity'
[]
[momentum_y_diffusion]
type = Diffusion
variable = v
[]
[momentum_y_pressure]
type = PressureGradient
integrate_p_by_parts = false
variable = v
pressure = pressure
component = 1
[]
[v_forcing]
type = BodyForce
variable = v
function = forcing_v
[]
[mass]
type = ADConservativeAdvection
variable = pressure
velocity = velocity
advected_quantity = -1
[]
[]
[DGKernels]
[momentum_x_convection]
type = ADDGAdvection
variable = u
velocity = 'velocity'
[]
[momentum_x_diffusion]
type = DGDiffusion
variable = u
sigma = 6
epsilon = -1
[]
[momentum_y_convection]
type = ADDGAdvection
variable = v
velocity = 'velocity'
[]
[momentum_y_diffusion]
type = DGDiffusion
variable = v
sigma = 6
epsilon = -1
[]
[]
[BCs]
[u_walls]
type = DGFunctionDiffusionDirichletBC
boundary = 'left bottom right top'
variable = u
sigma = 6
epsilon = -1
function = exact_u
[]
[v_walls]
type = DGFunctionDiffusionDirichletBC
boundary = 'left bottom right top'
variable = v
sigma = 6
epsilon = -1
function = exact_v
[]
[pressure_pin]
type = FunctionDirichletBC
variable = pressure
boundary = 'pinned_node'
function = 'exact_p'
[]
[]
[Materials]
[rho]
type = ADGenericConstantMaterial
prop_names = 'rho'
prop_values = '${rho}'
[]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = u
v = v
[]
[rhou]
type = ADParsedMaterial
property_name = 'rhou'
coupled_variables = 'u'
material_property_names = 'rho'
expression = 'rho*u'
[]
[rhov]
type = ADParsedMaterial
property_name = 'rhov'
coupled_variables = 'v'
material_property_names = 'rho'
expression = 'rho*v'
[]
[]
[Functions]
[exact_u]
type = ParsedFunction
expression = 'x^2*(1-x)^2*(2*y-6*y^2+4*y^3)'
[]
[exact_v]
type = ParsedFunction
expression = '-y^2*(1-y)^2*(2*x-6*x^2+4*x^3)'
[]
[exact_p]
type = ParsedFunction
expression = 'x*(1-x)-2/12'
[]
[forcing_u]
type = ParsedFunction
expression = '-4*mu/rho*(-1+2*y)*(y^2-6*x*y^2+6*x^2*y^2-y+6*x*y-6*x^2*y+3*x^2-6*x^3+3*x^4)+1-2*x+4*x^3'
'*y^2*(2*y^2-2*y+1)*(y-1)^2*(-1+2*x)*(x-1)^3'
symbol_names = 'mu rho'
symbol_values = '${mu} ${rho}'
[]
[forcing_v]
type = ParsedFunction
expression = '4*mu/rho*(-1+2*x)*(x^2-6*y*x^2+6*x^2*y^2-x+6*x*y-6*x*y^2+3*y^2-6*y^3+3*y^4)+4*y^3*x^2*(2'
'*x^2-2*x+1)*(x-1)^2*(-1+2*y)*(y-1)^3'
symbol_names = 'mu rho'
symbol_values = '${mu} ${rho}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu NONZERO mumps'
nl_rel_tol = 1e-12
[]
[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]
variable = v
function = exact_v
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2p]
variable = pressure
function = exact_p
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/mms-second-lid-driven.i)
mu = 2
rho = 2
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = TRI6
[]
[]
[Variables]
[vel_x]
family = MONOMIAL
order = SECOND
[]
[vel_y]
family = MONOMIAL
order = SECOND
[]
[pressure]
family = MONOMIAL
order = FIRST
[]
[vel_bar_x]
family = SIDE_HIERARCHIC
order = SECOND
[]
[vel_bar_y]
family = SIDE_HIERARCHIC
order = SECOND
[]
[pressure_bar]
family = SIDE_HIERARCHIC
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[HDGKernels]
[momentum_x_convection]
type = AdvectionIPHDGKernel
variable = vel_x
face_variable = vel_bar_x
velocity = 'velocity'
coeff = ${rho}
[]
[momentum_x_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_x
face_variable = vel_bar_x
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 0
[]
[momentum_y_convection]
type = AdvectionIPHDGKernel
variable = vel_y
face_variable = vel_bar_y
velocity = 'velocity'
coeff = ${rho}
[]
[momentum_y_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_y
face_variable = vel_bar_y
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 1
[]
[pressure_convection]
type = AdvectionIPHDGKernel
variable = pressure
face_variable = pressure_bar
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
[]
[]
[Kernels]
[momentum_x_ffn]
type = BodyForce
variable = vel_x
function = forcing_u
[]
[momentum_y_ffn]
type = BodyForce
variable = vel_y
function = forcing_v
[]
[mass_ffn]
type = BodyForce
variable = pressure
function = forcing_p
[]
[mean_zero_pressure]
type = ScalarLagrangeMultiplier
variable = pressure
lambda = lambda
[]
[]
[ScalarKernels]
[mean_zero_pressure_lm]
type = AverageValueConstraint
variable = lambda
pp_name = pressure_integral
value = 0
[]
[]
[BCs]
[momentum_x_diffusion_all]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right top'
variable = vel_x
face_variable = vel_bar_x
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = exact_u
diffusivity = 'mu'
component = 0
[]
[momentum_y_diffusion_all]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right top'
variable = vel_y
face_variable = vel_bar_y
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = exact_v
diffusivity = 'mu'
component = 1
[]
[mass_convection]
type = AdvectionIPHDGPrescribedFluxBC
face_variable = pressure_bar
variable = pressure
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
boundary = 'left bottom top right'
prescribed_normal_flux = 0
[]
[]
[Functions]
[exact_u]
type = ParsedFunction
expression = 'sin(y)*cos((1/2)*x*pi)'
[]
[forcing_u]
type = ParsedFunction
expression = 'mu*sin(y)*cos((1/2)*x*pi) + (1/4)*pi^2*mu*sin(y)*cos((1/2)*x*pi) - 1/2*pi*rho*sin(x)*sin(y)*sin((1/2)*y*pi)*cos((1/2)*x*pi) + rho*sin(x)*cos(y)*cos((1/2)*x*pi)*cos((1/2)*y*pi) - pi*rho*sin(y)^2*sin((1/2)*x*pi)*cos((1/2)*x*pi) + sin(y)*cos(x)'
symbol_names = 'mu rho'
symbol_values = '${mu} ${rho}'
[]
[exact_v]
type = ParsedFunction
expression = 'sin(x)*cos((1/2)*y*pi)'
[]
[forcing_v]
type = ParsedFunction
expression = 'mu*sin(x)*cos((1/2)*y*pi) + (1/4)*pi^2*mu*sin(x)*cos((1/2)*y*pi) - pi*rho*sin(x)^2*sin((1/2)*y*pi)*cos((1/2)*y*pi) - 1/2*pi*rho*sin(x)*sin(y)*sin((1/2)*x*pi)*cos((1/2)*y*pi) + rho*sin(y)*cos(x)*cos((1/2)*x*pi)*cos((1/2)*y*pi) + sin(x)*cos(y)'
symbol_names = 'mu rho'
symbol_values = '${mu} ${rho}'
[]
[exact_p]
type = ParsedFunction
expression = 'sin(x)*sin(y)'
[]
[forcing_p]
type = ParsedFunction
expression = '(1/2)*pi*rho*sin(x)*sin((1/2)*y*pi) + (1/2)*pi*rho*sin(y)*sin((1/2)*x*pi)'
symbol_names = 'rho'
symbol_values = '${rho}'
[]
[]
[Materials]
[const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${rho} ${mu}'
[]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = vel_x
v = vel_y
[]
[]
[AuxVariables]
[vel_exact_x]
[]
[vel_exact_y]
[]
[pressure_exact]
[]
[]
[AuxKernels]
[vel_exact_x]
type = FunctionAux
variable = vel_exact_x
function = exact_u
[]
[vel_exact_y]
type = FunctionAux
variable = vel_exact_y
function = exact_v
[]
[pressure_exact]
type = FunctionAux
variable = pressure_exact
function = exact_p
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_rel_tol = 1e-12
[]
[Outputs]
[out]
type = Exodus
hide = 'lambda pressure_integral'
[]
csv = true
[]
[Postprocessors]
[pressure_integral]
type = ElementIntegralVariablePostprocessor
variable = pressure
execute_on = linear
[]
[h]
type = AverageElementSize
execute_on = 'timestep_end'
[]
[L2u]
type = ElementL2Error
variable = vel_x
function = exact_u
execute_on = 'timestep_end'
[]
[L2v]
type = ElementL2Error
variable = vel_y
function = exact_v
execute_on = 'timestep_end'
[]
[L2p]
type = ElementL2Error
variable = pressure
function = exact_p
execute_on = 'timestep_end'
[]
[]
(modules/navier_stokes/test/tests/finite_element/ins/RZ_cone/unstabilized-velocity-component-objects.i)
[Mesh]
file = '2d_cone.msh'
coord_type = RZ
[]
[Variables]
[vel_x]
order = SECOND
[]
[vel_y]
order = SECOND
[]
[p][]
[]
[Kernels]
[momentum_x_time]
type = TimeDerivative
variable = vel_x
[]
[momentum_x_convection]
type = ADAdvection
variable = vel_x
velocity = 'velocity'
[]
[momentum_x_diffusion]
type = MatDiffusion
variable = vel_x
diffusivity = 1
[]
[momentum_x_diffusion_rz]
type = ADMomentumViscousRZ
variable = vel_x
mu_name = 1
component = 0
[]
[momentum_x_pressure]
type = PressureGradient
integrate_p_by_parts = true
variable = vel_x
pressure = p
component = 0
[]
[momentum_y_time]
type = TimeDerivative
variable = vel_y
[]
[momentum_y_convection]
type = ADAdvection
variable = vel_y
velocity = 'velocity'
[]
[momentum_y_diffusion]
type = MatDiffusion
variable = vel_y
diffusivity = 1
[]
[momentum_y_diffusion_rz]
type = ADMomentumViscousRZ
variable = vel_y
mu_name = 1
component = 1
[]
[momentum_y_pressure]
type = PressureGradient
integrate_p_by_parts = true
variable = vel_y
pressure = p
component = 1
[]
[mass]
type = ADMassAdvection
variable = p
vel_x = vel_x
vel_y = vel_y
[]
[]
[BCs]
[u_in]
type = DirichletBC
boundary = bottom
variable = vel_x
value = 0
[]
[v_in]
type = FunctionDirichletBC
boundary = bottom
variable = vel_y
function = 'inlet_func'
[]
[u_axis_and_walls]
type = DirichletBC
boundary = 'left right'
variable = vel_x
value = 0
[]
[v_no_slip]
type = DirichletBC
boundary = 'right'
variable = vel_y
value = 0
[]
[]
[Materials]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = vel_x
v = vel_y
[]
[]
[Functions]
[inlet_func]
type = ParsedFunction
expression = '-4 * x^2 + 1'
[]
[]
[Executioner]
type = Transient
dt = 0.005
dtmin = 0.005
num_steps = 5
l_max_its = 100
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_levels'
petsc_options_value = 'bjacobi ilu 4'
nl_rel_tol = 1e-12
nl_max_its = 6
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/cg-dg-hybrid/mms/lid-driven/hybrid-cg-dg-mms.i)
rho=1.1
mu=1.1
cp=1.1
k=1.1
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = -1
xmax = 1.0
ymin = -1
ymax = 1.0
nx = 2
ny = 2
[]
[corner_node]
type = ExtraNodesetGenerator
new_boundary = 'pinned_node'
nodes = '0'
input = gen
[]
[]
[Variables]
[u]
family = MONOMIAL
[]
[v]
family = MONOMIAL
[]
[pressure][]
[T]
family = MONOMIAL
[]
[]
[Kernels]
[momentum_x_convection]
type = ADConservativeAdvection
variable = u
velocity = 'velocity'
advected_quantity = 'rhou'
[]
[momentum_x_diffusion]
type = MatDiffusion
variable = u
diffusivity = 'mu'
[]
[momentum_x_pressure]
type = PressureGradient
integrate_p_by_parts = false
variable = u
pressure = pressure
component = 0
[]
[u_forcing]
type = BodyForce
variable = u
function = forcing_u
[]
[momentum_y_convection]
type = ADConservativeAdvection
variable = v
velocity = 'velocity'
advected_quantity = 'rhov'
[]
[momentum_y_diffusion]
type = MatDiffusion
variable = v
diffusivity = 'mu'
[]
[momentum_y_pressure]
type = PressureGradient
integrate_p_by_parts = false
variable = v
pressure = pressure
component = 1
[]
[v_forcing]
type = BodyForce
variable = v
function = forcing_v
[]
[mass]
type = ADConservativeAdvection
variable = pressure
velocity = velocity
advected_quantity = -1
[]
[p_forcing]
type = BodyForce
variable = pressure
function = forcing_p
[]
[T_convection]
type = ADConservativeAdvection
variable = T
velocity = 'velocity'
advected_quantity = 'rho_cp_temp'
[]
[T_diffusion]
type = MatDiffusion
variable = T
diffusivity = 'k'
[]
[T_forcing]
type = BodyForce
variable = T
function = forcing_T
[]
[]
[DGKernels]
[momentum_x_convection]
type = ADDGAdvection
variable = u
velocity = 'velocity'
advected_quantity = 'rhou'
[]
[momentum_x_diffusion]
type = DGDiffusion
variable = u
sigma = 6
epsilon = -1
diff = 'mu'
[]
[momentum_y_convection]
type = ADDGAdvection
variable = v
velocity = 'velocity'
advected_quantity = 'rhov'
[]
[momentum_y_diffusion]
type = DGDiffusion
variable = v
sigma = 6
epsilon = -1
diff = 'mu'
[]
[T_convection]
type = ADDGAdvection
variable = T
velocity = 'velocity'
advected_quantity = 'rho_cp_temp'
[]
[T_diffusion]
type = DGDiffusion
variable = T
sigma = 6
epsilon = -1
diff = 'k'
[]
[]
[BCs]
[u_walls]
type = DGFunctionDiffusionDirichletBC
boundary = 'left bottom right top'
variable = u
sigma = 6
epsilon = -1
function = exact_u
diff = 'mu'
[]
[v_walls]
type = DGFunctionDiffusionDirichletBC
boundary = 'left bottom right top'
variable = v
sigma = 6
epsilon = -1
function = exact_v
diff = 'mu'
[]
[pressure_pin]
type = FunctionDirichletBC
variable = pressure
boundary = 'pinned_node'
function = 'exact_p'
[]
[T_walls]
type = DGFunctionDiffusionDirichletBC
boundary = 'left bottom right top'
variable = T
sigma = 6
epsilon = -1
function = exact_T
diff = 'k'
[]
[]
[Materials]
[const]
type = ADGenericConstantMaterial
prop_names = 'rho cp'
prop_values = '${rho} ${cp}'
[]
[const_reg]
type = GenericConstantMaterial
prop_names = 'mu k'
prop_values = '${mu} ${k}'
[]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = u
v = v
[]
[rhou]
type = ADParsedMaterial
property_name = 'rhou'
coupled_variables = 'u'
material_property_names = 'rho'
expression = 'rho*u'
[]
[rhov]
type = ADParsedMaterial
property_name = 'rhov'
coupled_variables = 'v'
material_property_names = 'rho'
expression = 'rho*v'
[]
[rho_cp]
type = ADParsedMaterial
property_name = 'rho_cp'
material_property_names = 'rho cp'
expression = 'rho*cp'
[]
[rho_cp_temp]
type = ADParsedMaterial
property_name = 'rho_cp_temp'
material_property_names = 'rho cp'
coupled_variables = 'T'
expression = 'rho*cp*T'
[]
[]
[Functions]
[exact_u]
type = ParsedFunction
expression = 'sin(y)*cos((1/2)*x*pi)'
[]
[forcing_u]
type = ParsedFunction
expression = 'mu*sin(y)*cos((1/2)*x*pi) + (1/4)*pi^2*mu*sin(y)*cos((1/2)*x*pi) - 1/2*pi*rho*sin(x)*sin(y)*sin((1/2)*y*pi)*cos((1/2)*x*pi) + rho*sin(x)*cos(y)*cos((1/2)*x*pi)*cos((1/2)*y*pi) - pi*rho*sin(y)^2*sin((1/2)*x*pi)*cos((1/2)*x*pi) + sin(y)*cos(x)'
symbol_names = 'mu rho'
symbol_values = '${mu} ${rho}'
[]
[exact_v]
type = ParsedFunction
expression = 'sin(x)*cos((1/2)*y*pi)'
[]
[forcing_v]
type = ParsedFunction
expression = 'mu*sin(x)*cos((1/2)*y*pi) + (1/4)*pi^2*mu*sin(x)*cos((1/2)*y*pi) - pi*rho*sin(x)^2*sin((1/2)*y*pi)*cos((1/2)*y*pi) - 1/2*pi*rho*sin(x)*sin(y)*sin((1/2)*x*pi)*cos((1/2)*y*pi) + rho*sin(y)*cos(x)*cos((1/2)*x*pi)*cos((1/2)*y*pi) + sin(x)*cos(y)'
symbol_names = 'mu rho'
symbol_values = '${mu} ${rho}'
[]
[exact_p]
type = ParsedFunction
expression = 'sin(x)*sin(y)'
[]
[forcing_p]
type = ParsedFunction
expression = '(1/2)*pi*sin(x)*sin((1/2)*y*pi) + (1/2)*pi*sin(y)*sin((1/2)*x*pi)'
[]
[exact_T]
type = ParsedFunction
expression = 'cos(x)*cos(y)'
[]
[forcing_T]
type = ParsedFunction
expression = '-cp*rho*sin(x)*sin(y)*cos(x)*cos((1/2)*y*pi) - cp*rho*sin(x)*sin(y)*cos(y)*cos((1/2)*x*pi) - 1/2*pi*cp*rho*sin(x)*sin((1/2)*y*pi)*cos(x)*cos(y) - 1/2*pi*cp*rho*sin(y)*sin((1/2)*x*pi)*cos(x)*cos(y) + 2*k*cos(x)*cos(y)'
symbol_names = 'rho cp k'
symbol_values = '${rho} ${cp} ${k}'
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type -pc_factor_mat_solver_type -mat_mumps_icntl_14'
petsc_options_value = 'lu NONZERO mumps 300'
nl_rel_tol = 1e-12
[]
[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]
variable = v
function = exact_v
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2T]
variable = T
function = exact_T
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[L2p]
variable = pressure
function = exact_p
type = ElementL2Error
outputs = 'console csv'
execute_on = 'timestep_end'
[]
[]
(test/tests/dgkernels/passive-scalar-channel-flow/test.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = -1
ymax = 1
nx = 20
ny = 4
[]
[]
[Variables]
[u]
family = MONOMIAL
[]
[]
[Kernels]
[convection]
type = ADConservativeAdvection
variable = u
velocity = 'velocity'
[]
[diffusion]
type = MatDiffusion
variable = u
diffusivity = 1
[]
[]
[DGKernels]
[convection]
type = ADDGAdvection
variable = u
velocity = 'velocity'
[]
[diffusion]
type = DGDiffusion
variable = u
sigma = 6
epsilon = -1
diff = 1
[]
[]
[Functions]
[v_inlet]
type = ParsedVectorFunction
expression_x = '1'
[]
[]
[BCs]
[u_walls]
type = DGFunctionDiffusionDirichletBC
boundary = 'bottom top'
variable = u
sigma = 6
epsilon = -1
function = '0'
diff = 1
[]
[u_in]
type = ADConservativeAdvectionBC
boundary = 'left'
variable = u
velocity_function = v_inlet
primal_dirichlet_value = 1
[]
[u_out]
type = ADConservativeAdvectionBC
boundary = 'right'
variable = u
velocity_mat_prop = 'velocity'
[]
[]
[Materials]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = 1
v = 0
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(modules/navier_stokes/test/tests/finite_element/ins/hdg/ip/lid-driven/matrix-testing.i)
mu = 1
rho = 1
U = 1
n = 5
l = 1
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = ${l}
ymin = 0
ymax = ${l}
nx = ${n}
ny = ${n}
elem_type = TRI6
[]
[]
[Problem]
type = PrintMatricesNSProblem
extra_tag_matrices = 'mass jump combined grad_div'
pressure_mass_matrix = 'mass'
velocity_mass_matrix = 'mass'
augmented_lagrange_matrices = 'jump grad_div combined'
u = vel_x
v = vel_y
pressure = pressure
pressure_bar = pressure_bar
print = false
[]
[Variables]
[vel_x]
family = MONOMIAL
order = FIRST
[]
[vel_y]
family = MONOMIAL
order = FIRST
[]
[pressure]
family = MONOMIAL
order = CONSTANT
[]
[vel_bar_x]
family = SIDE_HIERARCHIC
order = FIRST
[]
[vel_bar_y]
family = SIDE_HIERARCHIC
order = FIRST
[]
[pressure_bar]
family = SIDE_HIERARCHIC
order = FIRST
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[HDGKernels]
[momentum_x_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_x
face_variable = vel_bar_x
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 0
[]
[momentum_y_diffusion]
type = NavierStokesStressIPHDGKernel
variable = vel_y
face_variable = vel_bar_y
diffusivity = 'mu'
alpha = 6
pressure_variable = pressure
pressure_face_variable = pressure_bar
component = 1
[]
[pressure_convection]
type = AdvectionIPHDGKernel
variable = pressure
face_variable = pressure_bar
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
[]
[]
[Kernels]
[mean_zero_pressure]
type = ScalarLagrangeMultiplier
variable = pressure
lambda = lambda
[]
[mass_matrix_vel_x]
type = MassMatrix
variable = vel_x
matrix_tags = 'mass'
[]
[mass_matrix_vel_y]
type = MassMatrix
variable = vel_y
matrix_tags = 'mass'
[]
[mass_matrix_pressure]
type = MassMatrix
variable = pressure
matrix_tags = 'mass'
[]
[u_jump]
type = GradDiv
matrix_only = true
variable = vel_x
u = vel_x
v = vel_y
component = 0
vector_tags = ''
matrix_tags = 'grad_div combined'
[]
[v_jump]
type = GradDiv
matrix_only = true
variable = vel_y
u = vel_x
v = vel_y
component = 1
vector_tags = ''
matrix_tags = 'grad_div combined'
[]
[]
[ScalarKernels]
[mean_zero_pressure_lm]
type = AverageValueConstraint
variable = lambda
pp_name = pressure_integral
value = 0
[]
[]
[DGKernels]
[pb_mass]
type = MassMatrixDGKernel
variable = pressure_bar
matrix_tags = 'mass'
[]
[u_jump]
type = MassFluxPenalty
matrix_only = true
variable = vel_x
u = vel_x
v = vel_y
component = 0
vector_tags = ''
matrix_tags = 'jump combined'
[]
[v_jump]
type = MassFluxPenalty
matrix_only = true
variable = vel_y
u = vel_x
v = vel_y
component = 1
vector_tags = ''
matrix_tags = 'jump combined'
[]
[]
[BCs]
[momentum_x_diffusion_walls]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right'
variable = vel_x
face_variable = vel_bar_x
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '0'
diffusivity = 'mu'
component = 0
[]
[momentum_x_diffusion_top]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'top'
variable = vel_x
face_variable = vel_bar_x
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '${U}'
diffusivity = 'mu'
component = 0
[]
[momentum_y_diffusion_all]
type = NavierStokesStressIPHDGDirichletBC
boundary = 'left bottom right top'
variable = vel_y
face_variable = vel_bar_y
pressure_variable = pressure
pressure_face_variable = pressure_bar
alpha = 6
functor = '0'
diffusivity = 'mu'
component = 1
[]
[mass_convection]
type = AdvectionIPHDGPrescribedFluxBC
face_variable = pressure_bar
variable = pressure
velocity = 'velocity'
coeff = '${fparse -rho}'
self_advection = false
boundary = 'left bottom top right'
prescribed_normal_flux = 0
[]
[pb_mass]
type = MassMatrixIntegratedBC
variable = pressure_bar
matrix_tags = 'mass'
boundary = 'left right bottom top'
[]
[u_jump_walls]
type = MassFluxPenaltyBC
matrix_only = true
variable = vel_x
u = vel_x
v = vel_y
component = 0
vector_tags = ''
matrix_tags = 'jump combined'
boundary = 'left right bottom'
dirichlet_value = 'walls'
[]
[v_jump_walls]
type = MassFluxPenaltyBC
matrix_only = true
variable = vel_y
u = vel_x
v = vel_y
component = 1
vector_tags = ''
matrix_tags = 'jump combined'
boundary = 'left right bottom'
dirichlet_value = 'walls'
[]
[u_jump_top]
type = MassFluxPenaltyBC
matrix_only = true
variable = vel_x
u = vel_x
v = vel_y
component = 0
vector_tags = ''
matrix_tags = 'jump combined'
boundary = 'top'
dirichlet_value = 'top'
[]
[v_jump_top]
type = MassFluxPenaltyBC
matrix_only = true
variable = vel_y
u = vel_x
v = vel_y
component = 1
vector_tags = ''
matrix_tags = 'jump combined'
boundary = 'top'
dirichlet_value = 'top'
[]
[]
[Functions]
[top]
type = ParsedVectorFunction
value_x = ${U}
value_y = 0
[]
[walls]
type = ParsedVectorFunction
value_x = 0
value_y = 0
[]
[]
[Materials]
[const]
type = ADGenericConstantMaterial
prop_names = 'rho mu'
prop_values = '${rho} ${mu}'
[]
[vel]
type = ADVectorFromComponentVariablesMaterial
vector_prop_name = 'velocity'
u = vel_x
v = vel_y
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
nl_rel_tol = 1e-12
[]
[Outputs]
[out]
type = CSV
hide = 'pressure_integral lambda'
execute_on = 'timestep_end'
[]
[]
[Postprocessors]
[symmetric]
type = MatrixSymmetryCheck
execute_on = 'timestep_end'
[]
[pressure_integral]
type = ElementIntegralVariablePostprocessor
variable = pressure
execute_on = linear
[]
[fe_jump_and_upb_equiv]
type = MatrixEqualityCheck
mat1 = 'vel_pb_grad_div.mat'
mat2 = 'jump.mat'
[]
[fe_grad_div_and_up_equiv]
type = MatrixEqualityCheck
mat1 = 'vel_p_grad_div.mat'
mat2 = 'grad_div.mat'
[]
[fe_combined_and_upall_equiv]
type = MatrixEqualityCheck
mat1 = 'vel_all_p_grad_div.mat'
mat2 = 'combined.mat'
[]
[upb_grad_div_num_zero_eig]
type = MatrixEigenvalueCheck
mat = vel_pb_grad_div.mat
[]
[upb_div_grad_num_zero_eig]
type = MatrixEigenvalueCheck
mat = vel_pb_div_grad.mat
[]
[up_grad_div_num_zero_eig]
type = MatrixEigenvalueCheck
mat = vel_p_grad_div.mat
[]
[up_div_grad_num_zero_eig]
type = MatrixEigenvalueCheck
mat = vel_p_div_grad.mat
[]
[upall_grad_div_num_zero_eig]
type = MatrixEigenvalueCheck
mat = vel_all_p_grad_div.mat
[]
[upall_div_grad_num_zero_eig]
type = MatrixEigenvalueCheck
mat = vel_all_p_div_grad.mat
[]
[]