- boundaryThe list of boundary IDs from the mesh where this object applies
C++ Type:std::vector<BoundaryName>
Controllable:No
Description:The list of boundary IDs from the mesh where this object applies
- epsilonEpsilon
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:Epsilon
- functionThe forcing function.
C++ Type:FunctionName
Unit:(no unit assumed)
Controllable:No
Description:The forcing function.
- sigmaSigma
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:Sigma
- variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable that this residual object operates on
DGFunctionDiffusionDirichletBC
Diffusion Dirichlet boundary condition for discontinuous Galerkin method.
Note that these boundary conditions are specific to DG and to a diffusion problem. Using a Functions System for the Dirichlet boundary conditions means that the spatial and time dependence is either known or imposed.
More information about Dirichlet boundary conditions and their mathematical meaning may be found in the DirichletBC documentation, and more information may be found about the discontinuous Galerkin discretization in the DGKernels documentation
Example input syntax
In this example, a 2D diffusion problem is solved with DG. A Dirichlet boundary condition is imposed on all boundaries using the exact_fn
function. The epsilon
and sigma
parameters are DG parameters.
[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
active<<<{"description": "If specified only the blocks named will be visited and made active"}>>> = 'all'
[all]
type = DGFunctionDiffusionDirichletBC<<<{"description": "Diffusion Dirichlet boundary condition for discontinuous Galerkin method.", "href": "DGFunctionDiffusionDirichletBC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = u
boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = '0 1 2 3'
function<<<{"description": "The forcing function."}>>> = exact_fn
epsilon<<<{"description": "Epsilon"}>>> = -1
sigma<<<{"description": "Sigma"}>>> = 6
[]
[]
(test/tests/dgkernels/2d_diffusion_dg/2d_diffusion_dg_test.i)Input Parameters
- diff1The diffusion (or thermal conductivity or viscosity) coefficient.
Default:1
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:The diffusion (or thermal conductivity or viscosity) coefficient.
- displacementsThe displacements
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The displacements
- matrix_onlyFalseWhether this object is only doing assembly to matrices (no vectors)
Default:False
C++ Type:bool
Controllable:No
Description:Whether this object is only doing assembly to matrices (no vectors)
- value0The value the variable should have on the boundary
Default:0
C++ Type:double
Unit:(no unit assumed)
Controllable:No
Description:The value the variable should have on the boundary
Optional Parameters
- absolute_value_vector_tagsThe tags for the vectors this residual object should fill with the absolute value of the residual contribution
C++ Type:std::vector<TagName>
Controllable:No
Description:The tags for the vectors this residual object should fill with the absolute value of the residual contribution
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
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
Contribution To Tagged Field Data 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.
- diag_save_inThe name of auxiliary variables to save this BC's diagonal jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The name of auxiliary variables to save this BC's diagonal jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- 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
- save_inThe name of auxiliary variables to save this BC's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The name of auxiliary variables to save this BC's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
- 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
- skip_execution_outside_variable_domainFalseWhether to skip execution of this boundary condition when the variable it applies to is not defined on the boundary. This can facilitate setups with moving variable domains and fixed boundaries. Note that the FEProblem boundary-restricted integrity checks will also need to be turned off if using this option
Default:False
C++ Type:bool
Controllable:No
Description:Whether to skip execution of this boundary condition when the variable it applies to is not defined on the boundary. This can facilitate setups with moving variable domains and fixed boundaries. Note that the FEProblem boundary-restricted integrity checks will also need to be turned off if using this option
- 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
- 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
- (test/tests/controls/time_periods/dgkernels/dgkernels.i)
- (test/tests/time_integrators/implicit-euler/ie-monomials.i)
- (test/tests/dgkernels/dg_displacement/dg_displacement.i)
- (test/tests/dgkernels/2d_diffusion_dg/dg_stateful.i)
- (test/tests/adaptivity/dont-p-refine/test.i)
- (test/tests/dgkernels/advection_diffusion_mixed_bcs_test_resid_jac/dg_advection_diffusion_test.i)
- (test/tests/dgkernels/2d_diffusion_dg/2d_diffusion_dg_test.i)
- (test/tests/tag/2d_diffusion_dg_tag.i)
- (modules/heat_transfer/test/tests/sideset_heat_transfer/gap_thermal_ktemp_1D.i)
- (modules/navier_stokes/test/tests/finite_element/ins/cg-dg-hybrid/lid-driven/hybrid-cg-dg.i)
- (test/tests/interfacekernels/1d_interface/mixed_shapes.i)
- (modules/navier_stokes/test/tests/finite_element/ins/cg-dg-hybrid/mms/lid-driven-skewed/hybrid-skewed-vortex.i)
- (test/tests/preconditioners/pbp/pbp_dg_test.i)
- (test/tests/hdgkernels/block-restricted/test.i)
- (test/tests/dgkernels/2d_diffusion_dg/no_mallocs_with_action.i)
- (test/tests/restart/p_refinement_restart/steady.i)
- (test/tests/outputs/exodus/exodus_discontinuous.i)
- (modules/heat_transfer/test/tests/sideset_heat_transfer/gap_thermal_1D.i)
- (test/tests/dgkernels/2d_diffusion_dg/no_functor_additions.i)
- (test/tests/coord_type/coord_type_rz_integrated.i)
- (test/tests/restart/p_refinement_restart/restarted_steady.i)
- (test/tests/userobjects/domain-user-object/measure-conservation.i)
- (test/tests/dgkernels/passive-scalar-channel-flow/test.i)
- (test/tests/dgkernels/3d_diffusion_dg/3d_diffusion_dg_test.i)
- (test/tests/dgkernels/ad_dg_diffusion/2d_diffusion_ad_dg_test.i)
- (test/tests/dgkernels/3d_diffusion_dg/3d_diffusion_p_refinement.i)
- (test/tests/bounds/constant_bounds_elem.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/cg-dg-hybrid/mms/lid-driven/hybrid-cg-dg-mms.i)
- (test/tests/dgkernels/stateful-coupled-var/test.i)
(test/tests/dgkernels/2d_diffusion_dg/2d_diffusion_dg_test.i)
###########################################################
# This is a test of the Discontinuous Galerkin System.
# Discontinous basis functions are used (Monomials) and
# a Laplacian DGKernel contributes to the
# internal edges around each element. Jumps are allowed
# but penalized by this method.
#
# @Requirement F3.60
###########################################################
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
# xmin = -1
# xmax = 1
# ymin = -1
# ymax = 1
xmin = 0
xmax = 1
ymin = 0
ymax = 1
elem_type = QUAD4
[]
[Variables]
active = 'u'
[u]
order = FIRST
family = MONOMIAL
[InitialCondition]
type = ConstantIC
value = 1
[]
[]
[]
[Functions]
active = 'forcing_fn exact_fn'
[forcing_fn]
type = ParsedFunction
# function = -4.0+(x*x)+(y*y)
# function = x
# function = (x*x)-2.0
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
# function = (x*x*x)-6.0*x
[]
[exact_fn]
type = ParsedGradFunction
# function = x
# grad_x = 1
# grad_y = 0
# function = (x*x)+(y*y)
# grad_x = 2*x
# grad_y = 2*y
# function = (x*x)
# grad_x = 2*x
# grad_y = 0
expression = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
# function = (x*x*x)
# grad_x = 3*x*x
# grad_y = 0
[]
[]
[Kernels]
active = 'diff abs forcing'
[diff]
type = Diffusion
variable = u
[]
[abs] # u * v
type = Reaction
variable = u
[]
[forcing]
type = BodyForce
variable = u
function = forcing_fn
[]
[]
[DGKernels]
active = 'dg_diff'
[dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[]
[]
[BCs]
active = 'all'
[all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
# petsc_options = '-snes_mf'
# petsc_options_iname = '-pc_type -pc_hypre_type'
# petsc_options_value = 'hypre boomeramg'
# petsc_options = '-snes_mf'
# max_r_steps = 2
[Adaptivity]
steps = 2
refine_fraction = 1.0
coarsen_fraction = 0
max_h_level = 8
[]
nl_rel_tol = 1e-10
# nl_rel_tol = 1e-12
[]
[Postprocessors]
active = 'h dofs l2_err'
[h]
type = AverageElementSize
[]
[dofs]
type = NumDOFs
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
[]
[]
[Outputs]
file_base = out
exodus = true
[]
(test/tests/controls/time_periods/dgkernels/dgkernels.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
[]
[Adaptivity]
marker = uniform_marker
[Markers]
[uniform_marker]
type = UniformMarker
mark = REFINE
[]
[]
[]
[Variables]
[u]
order = FIRST
family = MONOMIAL
initial_condition = 1
[]
[]
[Functions]
[forcing_fn]
type = ParsedFunction
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
[]
[exact_fn]
type = ParsedGradFunction
expression = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[]
[]
[Kernels]
[diff]
type = Diffusion
variable = u
[]
[abs] # u * v
type = Reaction
variable = u
[]
[forcing]
type = BodyForce
variable = u
function = forcing_fn
[]
[]
[DGKernels]
[dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[]
[dg_diff2]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 4
[]
[]
[BCs]
[all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
num_steps = 4
dt = 1
nl_rel_tol = 1e-10
[]
[Outputs]
exodus = true
[]
[Controls]
[dg_problem]
type = TimePeriod
enable_objects = 'DGKernels/dg_diff2'
disable_objects = 'DGKernel::dg_diff'
start_time = '2'
execute_on = 'initial timestep_begin'
[]
[]
(test/tests/time_integrators/implicit-euler/ie-monomials.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
elem_type = QUAD4
[]
[Variables]
[./u]
order = FIRST
family = MONOMIAL
[../]
[]
[ICs]
[./u_ic]
type = ConstantIC
variable = u
value = 1
[../]
[]
[Functions]
active = 'forcing_fn exact_fn'
[./forcing_fn]
type = ParsedFunction
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
[../]
[./exact_fn]
type = ParsedGradFunction
value = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[../]
[]
[Kernels]
[./time]
type = TimeDerivative
variable = u
[../]
[./diff]
type = Diffusion
variable = u
[../]
[./abs] # u * v
type = Reaction
variable = u
[../]
[./forcing]
type = BodyForce
variable = u
function = forcing_fn
[../]
[]
[DGKernels]
[./dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[../]
[]
[BCs]
[./all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
solve_type = 'NEWTON'
[../]
[]
[Executioner]
type = Transient
nl_rel_tol = 1e-10
num_steps = 1
[]
[Outputs]
execute_on = 'timestep_end'
console = true
[]
(test/tests/dgkernels/dg_displacement/dg_displacement.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
displacements = 'disp_x disp_y'
[]
[Variables]
[./u]
order = FIRST
family = MONOMIAL
[../]
[]
[AuxVariables]
[./disp_x]
[../]
[./disp_y]
[../]
[]
[Functions]
[./forcing_fn]
type = ParsedFunction
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
[../]
[./exact_fn]
type = ParsedGradFunction
expression = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[../]
[./disp_func]
type = ParsedFunction
expression = x
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[./abs]
type = Reaction
variable = u
[../]
[./forcing]
type = BodyForce
variable = u
function = forcing_fn
[../]
[]
[DGKernels]
[./dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
use_displaced_mesh = true
[../]
[]
[BCs]
[./all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
[../]
[]
[Executioner]
type = Steady
solve_type = PJFNK
nl_rel_tol = 1e-10
[]
[Outputs]
execute_on = 'timestep_end'
file_base = out
exodus = true
[]
[ICs]
[./disp_x_ic]
function = disp_func
variable = disp_x
type = FunctionIC
[../]
[]
(test/tests/dgkernels/2d_diffusion_dg/dg_stateful.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
elem_type = QUAD4
[]
[Variables]
[./u]
order = FIRST
family = MONOMIAL
[./InitialCondition]
type = ConstantIC
value = 1
[../]
[../]
[]
[Functions]
[./forcing_fn]
type = ParsedFunction
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
[../]
[./exact_fn]
type = ParsedGradFunction
expression = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[./abs]
type = Reaction
variable = u
[../]
[./forcing]
type = BodyForce
variable = u
function = forcing_fn
[../]
[]
[DGKernels]
[./dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[../]
[]
[BCs]
[./all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
[../]
[]
[Materials]
[./stateful]
type = StatefulMaterial
initial_diffusivity = 1
boundary = 'left'
[../]
[./general]
type = GenericConstantMaterial
block = '0'
prop_names = 'dummy'
prop_values = '1'
[../]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
nl_rel_tol = 1e-10
[]
(test/tests/adaptivity/dont-p-refine/test.i)
[Mesh]
[cmg]
type = CartesianMeshGenerator
dim = 2
dx = '1 1'
dy = '1 1'
ix = '2 2'
iy = '2 2'
subdomain_id = '0 0
0 1'
[]
[]
[Adaptivity]
switch_h_to_p_refinement = true
initial_marker = uniform
initial_steps = 1
[Markers/uniform]
type = UniformMarker
mark = REFINE
block = 1
[]
[]
[Variables]
[u]
family = MONOMIAL
order = FIRST
[]
[]
[AuxVariables]
[test][]
[]
[Kernels]
[diff]
type = Diffusion
variable = u
[]
[src]
type = BodyForce
variable = u
value = 1
[]
[]
[DGKernels]
[dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[]
[]
[BCs]
[left_u]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = 0
epsilon = -1
sigma = 6
[]
[]
[Postprocessors]
[avg]
type = ElementAverageValue
variable = u
[]
[]
[Executioner]
type = Steady
[]
[Outputs]
csv = true
[]
(test/tests/dgkernels/advection_diffusion_mixed_bcs_test_resid_jac/dg_advection_diffusion_test.i)
[Mesh]
type = GeneratedMesh
nx = 2
dim = 1
[]
[Kernels]
[./source]
type = BodyForce
variable = u
function = 'forcing_func'
[../]
[./convection]
type = ConservativeAdvection
variable = u
velocity = '1 0 0'
[../]
[./diffusion]
type = MatDiffusionTest
variable = u
prop_name = 'k'
[../]
[]
[DGKernels]
[./convection]
type = DGConvection
variable = u
velocity = '1 0 0'
[../]
[./diffusion]
type = DGDiffusion
variable = u
diff = 'k'
sigma = 6
epsilon = -1
[../]
[]
[BCs]
[./advection]
type = OutflowBC
boundary = 'right'
variable = u
velocity = '1 0 0'
[../]
[./diffusion_left]
type = DGFunctionDiffusionDirichletBC
boundary = 'left'
variable = u
sigma = 6
epsilon = -1
function = 'boundary_left_func'
diff = 'k'
[../]
[]
[Variables]
[./u]
family = MONOMIAL
order = THIRD
[../]
[]
[Materials]
[./test]
block = 0
type = GenericFunctionMaterial
prop_names = 'k'
prop_values = 'k_func'
[../]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
[../]
[]
[Functions]
[./forcing_func]
type = ParsedFunction
expression = '1'
[../]
[./boundary_left_func]
type = ParsedFunction
expression = '0'
[../]
[./k_func]
type = ParsedFunction
expression = '1 + x'
[../]
[]
[Outputs]
exodus = true
execute_on = 'timestep_end'
[]
(test/tests/dgkernels/2d_diffusion_dg/2d_diffusion_dg_test.i)
###########################################################
# This is a test of the Discontinuous Galerkin System.
# Discontinous basis functions are used (Monomials) and
# a Laplacian DGKernel contributes to the
# internal edges around each element. Jumps are allowed
# but penalized by this method.
#
# @Requirement F3.60
###########################################################
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
# xmin = -1
# xmax = 1
# ymin = -1
# ymax = 1
xmin = 0
xmax = 1
ymin = 0
ymax = 1
elem_type = QUAD4
[]
[Variables]
active = 'u'
[u]
order = FIRST
family = MONOMIAL
[InitialCondition]
type = ConstantIC
value = 1
[]
[]
[]
[Functions]
active = 'forcing_fn exact_fn'
[forcing_fn]
type = ParsedFunction
# function = -4.0+(x*x)+(y*y)
# function = x
# function = (x*x)-2.0
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
# function = (x*x*x)-6.0*x
[]
[exact_fn]
type = ParsedGradFunction
# function = x
# grad_x = 1
# grad_y = 0
# function = (x*x)+(y*y)
# grad_x = 2*x
# grad_y = 2*y
# function = (x*x)
# grad_x = 2*x
# grad_y = 0
expression = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
# function = (x*x*x)
# grad_x = 3*x*x
# grad_y = 0
[]
[]
[Kernels]
active = 'diff abs forcing'
[diff]
type = Diffusion
variable = u
[]
[abs] # u * v
type = Reaction
variable = u
[]
[forcing]
type = BodyForce
variable = u
function = forcing_fn
[]
[]
[DGKernels]
active = 'dg_diff'
[dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[]
[]
[BCs]
active = 'all'
[all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
# petsc_options = '-snes_mf'
# petsc_options_iname = '-pc_type -pc_hypre_type'
# petsc_options_value = 'hypre boomeramg'
# petsc_options = '-snes_mf'
# max_r_steps = 2
[Adaptivity]
steps = 2
refine_fraction = 1.0
coarsen_fraction = 0
max_h_level = 8
[]
nl_rel_tol = 1e-10
# nl_rel_tol = 1e-12
[]
[Postprocessors]
active = 'h dofs l2_err'
[h]
type = AverageElementSize
[]
[dofs]
type = NumDOFs
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
[]
[]
[Outputs]
file_base = out
exodus = true
[]
(test/tests/tag/2d_diffusion_dg_tag.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
xmin = 0
xmax = 1
ymin = 0
ymax = 1
elem_type = QUAD4
[]
[Variables]
[./u]
order = FIRST
family = MONOMIAL
[./InitialCondition]
type = ConstantIC
value = 1
[../]
[../]
[]
[AuxVariables]
[./tag_variable1]
order = FIRST
family = MONOMIAL
[../]
[./tag_variable2]
order = FIRST
family = MONOMIAL
[../]
[]
[AuxKernels]
[./TagVectorAux1]
type = TagVectorAux
variable = tag_variable1
v = u
vector_tag = vec_tag2
[../]
[./TagVectorAux2]
type = TagMatrixAux
variable = tag_variable2
v = u
matrix_tag = mat_tag2
[../]
[]
[Functions]
[./forcing_fn]
type = ParsedFunction
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
[../]
[./exact_fn]
type = ParsedGradFunction
value = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[./abs]
type = Reaction
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[./forcing]
type = BodyForce
variable = u
function = forcing_fn
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[]
[DGKernels]
[./dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[]
[BCs]
[./all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[../]
[]
[Problem]
type = TagTestProblem
test_tag_vectors = 'nontime residual vec_tag1 vec_tag2'
test_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_vectors = 'vec_tag1 vec_tag2'
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
nl_rel_tol = 1e-10
[]
[Postprocessors]
[./h]
type = AverageElementSize
[../]
[./dofs]
type = NumDOFs
[../]
[./l2_err]
type = ElementL2Error
variable = u
function = exact_fn
[../]
[]
[Outputs]
exodus = true
[]
(modules/heat_transfer/test/tests/sideset_heat_transfer/gap_thermal_ktemp_1D.i)
[Mesh]
[mesh]
type = GeneratedMeshGenerator
dim = 1
nx = 2
xmax = 2
[]
[split]
type = SubdomainBoundingBoxGenerator
input = mesh
block_id = 1
bottom_left = '1 0 0'
top_right = '2 0 0'
[]
[interface]
type = SideSetsBetweenSubdomainsGenerator
input = split
primary_block = 1
paired_block = 0
new_boundary = 'interface0'
[]
uniform_refine = 4
[]
[Variables]
[T]
order = FIRST
family = MONOMIAL
[]
[]
[AuxVariables]
[Tbulk]
order = FIRST
family = LAGRANGE
initial_condition = 300 # K
[]
[]
[Kernels]
[diff]
type = MatDiffusion
variable = T
diffusivity = conductivity
[]
[source]
type = BodyForce
variable = T
value = 1.0
[]
[]
[DGKernels]
[dg_diff]
type = DGDiffusion
variable = T
epsilon = -1
sigma = 6
diff = conductivity
exclude_boundary = 'interface0'
[]
[]
[InterfaceKernels]
[gap_var]
type = SideSetHeatTransferKernel
variable = T
neighbor_var = T
boundary = 'interface0'
Tbulk_var = Tbulk
[]
[]
[Functions]
# Defining temperature dependent fucntion for conductivity across side set
[kgap]
type = ParsedFunction
expression = 't / 200'
[]
[bc_func]
type = ConstantFunction
value = 300
[]
[exact]
type = ParsedFunction
expression = '
A := if(x < 1, -0.5, -0.25);
B := if(x < 1, -0.293209850655001, 0.0545267662299068);
C := if(x < 1, 300.206790149345, 300.19547323377);
d := -1;
A * (x+d) * (x+d) + B * (x+d) + C'
[]
[]
[BCs]
[bc_left]
type = DGFunctionDiffusionDirichletBC
boundary = 'left'
variable = T
diff = 'conductivity'
epsilon = -1
sigma = 6
function = bc_func
[]
[bc_right]
type = DGFunctionDiffusionDirichletBC
boundary = 'right'
variable = T
diff = 'conductivity'
epsilon = -1
sigma = 6
function = bc_func
[]
[]
[Materials]
[k0]
type = GenericConstantMaterial
prop_names = 'conductivity'
prop_values = 1.0
block = 0
[]
[k1]
type = GenericConstantMaterial
prop_names = 'conductivity'
prop_values = 2.0
block = 1
[]
[gap_mat]
type = SideSetHeatTransferMaterial
boundary = 'interface0'
# Using temperature dependent function for gap conductivity
conductivity_temperature_function = kgap
# Variable to evaluate conductivity with
gap_temperature = Tbulk
gap_length = 1.0
h_primary = 1
h_neighbor = 1
emissivity_primary = 1
emissivity_neighbor = 1
[]
[]
[Postprocessors]
[error]
type = ElementL2Error
variable = T
function = exact
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(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}'
[]
[]
(test/tests/interfacekernels/1d_interface/mixed_shapes.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 1
nx = 10
xmax = 2
[]
[./subdomain1]
input = gen
type = SubdomainBoundingBoxGenerator
bottom_left = '1.0 0 0'
block_id = 1
top_right = '2.0 1.0 0'
[../]
[./interface]
input = subdomain1
type = SideSetsBetweenSubdomainsGenerator
primary_block = '0'
paired_block = '1'
new_boundary = 'primary0_interface'
[../]
[./interface_again]
input = interface
type = SideSetsBetweenSubdomainsGenerator
primary_block = '1'
paired_block = '0'
new_boundary = 'primary1_interface'
[../]
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
block = '0'
[../]
[./v]
order = FIRST
family = MONOMIAL
block = '1'
[../]
[]
[Kernels]
[./diff_u]
type = CoeffParamDiffusion
variable = u
D = 4
block = 0
[../]
[./diff_v]
type = CoeffParamDiffusion
variable = v
D = 2
block = 1
[../]
[./body_u]
type = BodyForce
variable = u
block = 0
function = 'x^3+x^2+x+1'
[../]
[./body_v]
type = BodyForce
variable = v
block = 1
function = 'x^3+x^2+x+1'
[../]
[]
[DGKernels]
[./dg_diff_v]
type = DGDiffusion
variable = v
block = 1
diff = 2
sigma = 6
epsilon = -1
[../]
[]
[InterfaceKernels]
[./interface]
type = OneSideDiffusion
variable = u
neighbor_var = v
boundary = primary0_interface
D = 4
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = 'left'
value = 1
[../]
# [./right]
# type = DirichletBC
# variable = v
# boundary = 'right'
# value = 0
# [../]
[./right]
type = DGFunctionDiffusionDirichletBC
variable = v
boundary = 'right'
function = 0
epsilon = -1
sigma = 6
[../]
[./middle]
type = NeumannBC
variable = u
boundary = 'primary0_interface'
value = '.5'
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Steady
solve_type = NEWTON
[]
[Outputs]
exodus = true
print_linear_residuals = true
[]
[Debug]
show_var_residual_norms = true
[]
(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'
[]
[]
(test/tests/preconditioners/pbp/pbp_dg_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 20
ny = 20
xmin = 0
xmax = 1
ymin = 0
ymax = 1
elem_type = QUAD4
[]
[Variables]
[./u]
order = FIRST
family = MONOMIAL
[../]
[./v]
order = FIRST
family = MONOMIAL
[../]
[]
[Preconditioning]
[./PBP]
type = PBP
solve_order = 'u v'
preconditioner = 'AMG AMG'
[../]
[]
[Kernels]
[./diff_u]
type = Diffusion
variable = u
[../]
[./abs_u]
type = Reaction
variable = u
[../]
[./forcing_u]
type = BodyForce
variable = u
function = forcing_fn
[../]
[./diff_v]
type = Diffusion
variable = v
[../]
[./abs_v]
type = Reaction
variable = v
[../]
[./forcing_v]
type = BodyForce
variable = v
function = forcing_fn
[../]
[./conv_v]
type = CoupledForce
variable = v
v = u
[../]
[]
[DGKernels]
[./dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[../]
[./dg_diff_2]
type = DGDiffusion
variable = v
epsilon = -1
sigma = 6
[../]
[]
[Functions]
[./forcing_fn]
type = ParsedFunction
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
[../]
[./exact_fn]
type = ParsedGradFunction
value = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[../]
[]
[BCs]
[./all_u]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
[../]
[./all_v]
type = DGFunctionDiffusionDirichletBC
variable = v
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
[../]
[]
[Problem]
type = FEProblem
error_on_jacobian_nonzero_reallocation = true
[]
[Executioner]
type = Steady
l_max_its = 10
nl_max_its = 10
solve_type = JFNK
[]
[Outputs]
exodus = true
[]
(test/tests/hdgkernels/block-restricted/test.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 4
ny = 1
elem_type = QUAD9
xmax = 4
[]
[different_sub]
type = SubdomainBoundingBoxGenerator
block_id = 1
bottom_left = '2 0 0'
input = gen
top_right = '4 1 0'
[]
[left_sideset_for_right_sub]
type = SideSetsBetweenSubdomainsGenerator
input = different_sub
new_boundary = left_1
paired_block = 0
primary_block = 1
[]
[break]
type = BreakBoundaryOnSubdomainGenerator
input = left_sideset_for_right_sub
[]
[right_sideset_for_left_sub]
type = SideSetsBetweenSubdomainsGenerator
input = break
new_boundary = right_0
paired_block = 1
primary_block = 0
[]
[]
[Problem]
nl_sys_names = 'nl0 nl1'
kernel_coverage_check = false
[]
[Variables]
[u]
order = FIRST
family = MONOMIAL
block = 0
solver_sys = 'nl0'
[]
[side_u]
order = FIRST
family = SIDE_HIERARCHIC
block = 0
solver_sys = 'nl0'
[]
[v]
order = FIRST
family = MONOMIAL
block = 1
solver_sys = 'nl1'
[]
[]
[Kernels]
[diff_v]
type = Diffusion
variable = v
[]
[]
[HDGKernels]
[diff]
type = DiffusionIPHDGKernel
variable = u
face_variable = side_u
diffusivity = 1
alpha = 6
[]
[]
[DGKernels]
[diff]
type = DGDiffusion
epsilon = -1
sigma = 6
variable = v
[]
[]
[BCs]
[left_u]
type = DiffusionIPHDGDirichletBC
functor = 0
boundary = 'left'
alpha = 6
variable = u
face_variable = side_u
diffusivity = 1
[]
[zero_flux_u]
type = DiffusionIPHDGPrescribedFluxBC
boundary = 'top_to_0 bottom_to_0'
prescribed_normal_flux = 0
variable = u
face_variable = side_u
alpha = 6
diffusivity = 1
[]
[right_u]
type = DiffusionIPHDGDirichletBC
functor = 1
boundary = 'right_0'
alpha = 6
variable = u
face_variable = side_u
diffusivity = 1
[]
[left_v]
type = DGFunctionDiffusionDirichletBC
boundary = left_1
epsilon = -1
function = 0
sigma = 6
variable = v
[]
[right_v]
type = DGFunctionDiffusionDirichletBC
boundary = right
epsilon = -1
function = 1
sigma = 6
variable = v
[]
[]
[Executioner]
type = Steady
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_type'
petsc_options_value = 'lu mumps'
[]
[Outputs]
csv = true
[]
[Postprocessors]
[avg_u]
type = ElementAverageValue
variable = u
block = 0
[]
[avg_side_u]
type = ElementAverageValue
variable = side_u
block = 0
[]
[avg_v]
type = ElementAverageValue
variable = v
block = 1
[]
[]
(test/tests/dgkernels/2d_diffusion_dg/no_mallocs_with_action.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
elem_type = QUAD4
[]
[Variables]
[./u]
order = FIRST
family = MONOMIAL
[./InitialCondition]
type = ConstantIC
value = 1
[../]
[../]
[]
[AuxVariables]
[v]
order = FIRST
family = MONOMIAL
[]
[]
[Functions]
[./forcing_fn]
type = ParsedFunction
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
[../]
[./exact_fn]
type = ParsedGradFunction
expression = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[./abs] # u * v
type = Reaction
variable = u
[../]
[./forcing]
type = BodyForce
variable = u
function = forcing_fn
[../]
[]
[DGDiffusionAction]
variable = u
epsilon = -1
sigma = 6
# We couple in an auxiliary variable in order to ensure that we've properly
# ghosted both non-linear and auxiliary solution vectors
coupled_var = v
[]
[BCs]
[./all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
[../]
[]
[Executioner]
type = Steady
[]
[Outputs]
exodus = true
[console]
type = Console
system_info = 'framework mesh aux nonlinear relationship execution'
[]
[]
[Problem]
error_on_jacobian_nonzero_reallocation = true
[]
[Postprocessors]
active = 'num_rm'
[num_rm]
type = NumRelationshipManagers
[]
[num_internal_sides]
type = NumInternalSides
[]
[]
(test/tests/restart/p_refinement_restart/steady.i)
[Mesh]
[cmg]
type = CartesianMeshGenerator
dim = 2
dx = '1 1'
dy = '1 1'
ix = '2 2'
iy = '2 2'
subdomain_id = '0 0
0 1'
[]
[]
[Adaptivity]
switch_h_to_p_refinement = true
initial_marker = uniform
initial_steps = 1
[Markers/uniform]
type = UniformMarker
mark = REFINE
block = 1
[]
[]
[Variables]
[u]
family = MONOMIAL
order = FIRST
[]
[]
[AuxVariables]
[test][]
[]
[Kernels]
[diff]
type = Diffusion
variable = u
[]
[src]
type = BodyForce
variable = u
value = 1
[]
[]
[DGKernels]
[dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[]
[]
[BCs]
[left_u]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = 0
epsilon = -1
sigma = 6
[]
[]
[Postprocessors]
[avg]
type = ElementAverageValue
variable = u
[]
[]
[Executioner]
type = Steady
nl_abs_tol = 1e-10
[]
[Outputs]
checkpoint = true
exodus = true
print_mesh_changed_info = true
[]
(test/tests/outputs/exodus/exodus_discontinuous.i)
##
# \file exodus/exodus_discontinuous.i
# \example exodus/exodus_discontinuous.i
# Input file for testing discontinuous data output
#
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[./disc_u]
family = monomial
order = first
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = disc_u
[../]
[./forcing]
type = BodyForce
variable = disc_u
value = 7
[../]
[]
[DGKernels]
[./diff_dg]
type = DGDiffusion
variable = disc_u
sigma = 1
epsilon = 1
[../]
[]
[Functions]
[./zero_fn]
type = ParsedFunction
expression = 0.0
[../]
[]
[BCs]
[./all]
type = DGFunctionDiffusionDirichletBC
variable = disc_u
boundary = 'left right top bottom'
function = zero_fn
sigma = 1
epsilon = 1
[../]
[]
[Executioner]
type = Steady
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
execute_on = 'timestep_end'
[./exo_out]
type = Exodus
discontinuous = true
file_base = 'exodus_discontinuous_out'
[../]
[]
(modules/heat_transfer/test/tests/sideset_heat_transfer/gap_thermal_1D.i)
[Mesh]
[mesh]
type = GeneratedMeshGenerator
dim = 1
nx = 2
xmax = 2
[]
[split]
type = SubdomainBoundingBoxGenerator
input = mesh
block_id = 1
bottom_left = '1 0 0'
top_right = '2 0 0'
[]
[interface]
type = SideSetsBetweenSubdomainsGenerator
input = split
primary_block = 1
paired_block = 0
new_boundary = 'interface0'
[]
uniform_refine = 4
[]
[Variables]
# Defining a DFEM variable to handle gap discontinuity
[T]
order = FIRST
family = MONOMIAL
[]
[]
[AuxVariables]
# Auxvariable containing bulk temperature of gap
[Tbulk]
order = FIRST
family = LAGRANGE
initial_condition = 300 # K
[]
[]
[Kernels]
[diff]
type = MatDiffusion
variable = T
diffusivity = conductivity
[]
[source]
type = BodyForce
variable = T
value = 1.0
[]
[]
[DGKernels]
# DG kernel to represent diffusion accross element faces
[./dg_diff]
type = DGDiffusion
variable = T
epsilon = -1
sigma = 6
diff = conductivity
# Ignoring gap side set because no diffusion accross there
exclude_boundary = 'interface0'
[../]
[]
[InterfaceKernels]
active = 'gap'
# Heat transfer kernel using Tbulk as material
[gap]
type = SideSetHeatTransferKernel
variable = T
neighbor_var = T
boundary = 'interface0'
[]
# Heat transfer kernel using Tbulk as auxvariable
[gap_var]
type = SideSetHeatTransferKernel
variable = T
neighbor_var = T
boundary = 'interface0'
Tbulk_var = Tbulk
[]
[]
[Functions]
[bc_func]
type = ConstantFunction
value = 300
[]
[exact]
type = ParsedFunction
expression = '
A := if(x < 1, -0.5, -0.25);
B := if(x < 1, -0.293209850655001, 0.0545267662299068);
C := if(x < 1, 300.206790149345, 300.19547323377);
d := -1;
A * (x+d) * (x+d) + B * (x+d) + C'
[]
[]
[BCs]
[bc_left]
type = DGFunctionDiffusionDirichletBC
boundary = 'left'
variable = T
diff = 'conductivity'
epsilon = -1
sigma = 6
function = bc_func
[]
[bc_right]
type = DGFunctionDiffusionDirichletBC
boundary = 'right'
variable = T
diff = 'conductivity'
epsilon = -1
sigma = 6
function = bc_func
[]
[]
[Materials]
[k0]
type = GenericConstantMaterial
prop_names = 'conductivity'
prop_values = 1.0
block = 0
[]
[k1]
type = GenericConstantMaterial
prop_names = 'conductivity'
prop_values = 2.0
block = 1
[]
[gap_mat]
type = SideSetHeatTransferMaterial
boundary = 'interface0'
conductivity = 1.5
gap_length = 1.0
h_primary = 1
h_neighbor = 1
Tbulk = 300
emissivity_primary = 1
emissivity_neighbor = 1
[]
[]
[Postprocessors]
[error]
type = ElementL2Error
variable = T
function = exact
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(test/tests/dgkernels/2d_diffusion_dg/no_functor_additions.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
elem_type = QUAD4
[]
[Variables]
[./u]
order = FIRST
family = MONOMIAL
[./InitialCondition]
type = ConstantIC
value = 1
[../]
[../]
[]
[AuxVariables]
[v]
order = FIRST
family = MONOMIAL
[]
[]
[Functions]
[./forcing_fn]
type = ParsedFunction
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
[../]
[./exact_fn]
type = ParsedGradFunction
expression = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[./abs] # u * v
type = Reaction
variable = u
[../]
[./forcing]
type = BodyForce
variable = u
function = forcing_fn
[../]
[]
[DGKernels]
[regular_dg_diffusion]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[]
[]
[DGDiffusionAction]
variable = u
kernels_to_add = 'COUPLED'
coupled_var = v
[]
[BCs]
[./all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
[../]
[]
[Executioner]
type = Steady
[]
[Outputs]
exodus = true
[console]
type = Console
system_info = 'framework mesh aux nonlinear relationship execution'
[]
[]
[Problem]
error_on_jacobian_nonzero_reallocation = true
[]
[Postprocessors]
[num_rm]
type = NumRelationshipManagers
[]
[]
(test/tests/coord_type/coord_type_rz_integrated.i)
[Mesh]
type = GeneratedMesh
nx = 10
xmax = 1
ny = 10
ymax = 1
dim = 2
allow_renumbering = false
[]
[Problem]
type = FEProblem
coord_type = RZ
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Steady
solve_type = NEWTON
[]
[Outputs]
[./out]
type = Exodus
[../]
[]
[Kernels]
[./diff_u]
type = Diffusion
variable = u
[../]
[]
[DGKernels]
[./dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[../]
[]
[Variables]
[./u]
order = FIRST
family = MONOMIAL
[../]
[]
[BCs]
[./source]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = 'right'
function = exact_fn
epsilon = -1
sigma = 6
[../]
[./vacuum]
boundary = 'top'
type = VacuumBC
variable = u
[../]
[]
[Functions]
[./exact_fn]
type = ConstantFunction
value = 1
[../]
[]
[ICs]
[./u]
type = ConstantIC
value = 1
variable = u
[../]
[]
(test/tests/restart/p_refinement_restart/restarted_steady.i)
[Mesh]
[cmg]
type = FileMeshGenerator
file = steady_out_cp/LATEST
skip_partitioning = true
[]
[]
[Problem]
restart_file_base = steady_out_cp/LATEST
[]
[Variables]
[u]
family = MONOMIAL
order = FIRST
[]
[]
[AuxVariables]
[test][]
[]
[Kernels]
[diff]
type = Diffusion
variable = u
[]
[src]
type = BodyForce
variable = u
value = 1
[]
[]
[DGKernels]
[dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[]
[]
[BCs]
[left_u]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = 0
epsilon = -1
sigma = 6
[]
[]
[Postprocessors]
[avg]
type = ElementAverageValue
variable = u
[]
[]
[Executioner]
type = Steady
nl_abs_tol = 1e-10
[]
[Outputs]
exodus = true
[]
(test/tests/userobjects/domain-user-object/measure-conservation.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
[]
[Variables]
[u]
order = FIRST
family = MONOMIAL
[]
[]
[UserObjects]
[test]
type = DGDiffusionDomainUserObject
function = 'x'
epsilon = -1
sigma = 6
u = u
diff = 'diff'
ad_diff = 'ad_diff'
[]
[]
[Kernels]
[diff]
type = MatDiffusion
variable = u
diffusivity = 'diff'
[]
[]
[DGKernels]
[dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
diff = 'diff'
[]
[]
[BCs]
[all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = 'left right'
function = 'x'
epsilon = -1
sigma = 6
diff = 'diff'
[]
[]
[Materials]
[constant]
type = GenericConstantMaterial
prop_names = 'diff'
prop_values = '2'
[]
[ad_constant]
type = ADGenericConstantMaterial
prop_names = 'ad_diff'
prop_values = '2'
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-14
[]
[Outputs]
exodus = true
[]
(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
[]
(test/tests/dgkernels/3d_diffusion_dg/3d_diffusion_dg_test.i)
[Mesh]
type = GeneratedMesh
dim = 3
nx = 5
ny = 5
nz = 5
xmin = 0
xmax = 1
ymin = 0
ymax = 1
zmin = 0
zmax = 1
elem_type = HEX8
[]
[Variables]
active = 'u'
[u]
order = FIRST
family = MONOMIAL
[InitialCondition]
type = ConstantIC
value = 0.5
[]
[]
[]
[Functions]
active = 'forcing_fn exact_fn'
[forcing_fn]
type = ParsedFunction
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
[]
[exact_fn]
type = ParsedGradFunction
expression = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[]
[]
[Kernels]
active = 'diff abs forcing'
[diff]
type = Diffusion
variable = u
[]
[abs] # u * v
type = Reaction
variable = u
[]
[forcing]
type = BodyForce
variable = u
function = forcing_fn
[]
[]
[DGKernels]
active = 'dg_diff'
[dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[]
[]
[BCs]
active = 'all'
[all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3 4 5'
function = exact_fn
epsilon = -1
sigma = 6
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Postprocessors]
active = 'h dofs l2_err'
[h]
type = AverageElementSize
execute_on = 'initial timestep_end'
[]
[dofs]
type = NumDOFs
execute_on = 'initial timestep_end'
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Outputs]
file_base = out
exodus = true
[]
(test/tests/dgkernels/ad_dg_diffusion/2d_diffusion_ad_dg_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
elem_type = QUAD4
[]
[Variables]
[u]
order = FIRST
family = MONOMIAL
[InitialCondition]
type = ConstantIC
value = 1
[]
[]
[]
[Functions]
[forcing_fn]
type = ParsedFunction
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
[]
[exact_fn]
type = ParsedGradFunction
expression = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[]
[]
[Kernels]
[diff]
type = Diffusion
variable = u
[]
[abs] # u * v
type = Reaction
variable = u
[]
[forcing]
type = BodyForce
variable = u
function = forcing_fn
[]
[]
[DGKernels]
[dg_diff]
type = ADDGDiffusion
variable = u
epsilon = -1
sigma = 6
diff = diff
[]
[]
[Materials]
[ad_coupled_mat]
type = ADCoupledMaterial
coupled_var = u
ad_mat_prop = diff
regular_mat_prop = diff_regular
[]
[]
[BCs]
[all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[Adaptivity]
steps = 2
refine_fraction = 1.0
coarsen_fraction = 0
max_h_level = 8
[]
nl_rel_tol = 1e-10
[]
[Postprocessors]
[h]
type = AverageElementSize
[]
[dofs]
type = NumDOFs
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
[]
[]
[Outputs]
exodus = true
[]
(test/tests/dgkernels/3d_diffusion_dg/3d_diffusion_p_refinement.i)
[Mesh]
type = GeneratedMesh
dim = 3
nx = 1
ny = 1
nz = 1
xmin = 0
xmax = 1
ymin = 0
ymax = 1
zmin = 0
zmax = 1
elem_type = HEX8
[]
[Variables]
[u]
order = FIRST
family = MONOMIAL
[InitialCondition]
type = ConstantIC
value = 0.5
[]
[]
[]
[Functions]
[forcing_fn]
type = ParsedFunction
expression = 2*pow(e,-x-(y*y))*(1-2*y*y)
[]
[exact_fn]
type = ParsedGradFunction
expression = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[]
[]
[Kernels]
[diff]
type = Diffusion
variable = u
[]
[abs]
type = Reaction
variable = u
[]
[forcing]
type = BodyForce
variable = u
function = forcing_fn
[]
[]
[DGKernels]
[dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
[]
[]
[BCs]
[all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3 4 5'
function = exact_fn
epsilon = -1
sigma = 6
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[Adaptivity]
switch_h_to_p_refinement = true
steps = 2
refine_fraction = 1.0
coarsen_fraction = 0
max_h_level = 8
[]
[]
[Postprocessors]
[h]
type = AverageElementSize
execute_on = 'initial timestep_end'
[]
[dofs]
type = NumDOFs
execute_on = 'initial timestep_end'
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Outputs]
exodus = true
[]
(test/tests/bounds/constant_bounds_elem.i)
[Mesh]
type = GeneratedMesh
dim = 1
xmin = 0
xmax = 1
nx = 10
[]
[Variables]
[u]
order = CONSTANT
family = MONOMIAL
[]
[v]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxVariables]
[bounds_dummy]
order = CONSTANT
family = MONOMIAL
[]
[]
[Kernels]
[diff_u]
type = Diffusion
variable = u
[]
[reaction_u]
type = Reaction
variable = u
[]
[diff_v]
type = Diffusion
variable = v
[]
[reaction_v]
type = Reaction
variable = v
[]
[]
[DGKernels]
[dg_diff_u]
type = ADDGDiffusion
variable = u
epsilon = -1
sigma = 6
diff = 3
[]
[dg_diff_v]
type = ADDGDiffusion
variable = v
epsilon = -1
sigma = 6
diff = 4
[]
[]
[BCs]
[left_u]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0'
function = -0.5
epsilon = -1
sigma = 6
[]
[right_u]
type = NeumannBC
variable = u
boundary = 1
value = 30
[]
[left_v]
type = DGFunctionDiffusionDirichletBC
variable = v
boundary = '0'
function = 4
epsilon = -1
sigma = 6
[]
[right_v]
type = NeumannBC
variable = v
boundary = 1
value = -40
[]
[]
[Bounds]
[u_upper_bound]
type = ConstantBounds
variable = bounds_dummy
bounded_variable = u
bound_type = upper
bound_value = 1
[]
[u_lower_bound]
type = ConstantBounds
variable = bounds_dummy
bounded_variable = u
bound_type = lower
bound_value = 0
[]
[v_upper_bound]
type = ConstantBounds
variable = bounds_dummy
bounded_variable = v
bound_type = upper
bound_value = 3
[]
[v_lower_bound]
type = ConstantBounds
variable = bounds_dummy
bounded_variable = v
bound_type = lower
bound_value = -1
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
petsc_options_iname = '-snes_type'
petsc_options_value = 'vinewtonrsls'
[]
[Outputs]
exodus = true
[]
(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/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/stateful-coupled-var/test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
elem_type = QUAD4
[]
[Variables]
[u]
order = FIRST
family = MONOMIAL
[]
[]
[Functions]
[exact_fn]
type = ParsedGradFunction
expression = pow(e,-x-(y*y))
grad_x = -pow(e,-x-(y*y))
grad_y = -2*y*pow(e,-x-(y*y))
[]
[]
[Kernels]
[diff]
type = Diffusion
variable = u
[]
[]
[DGKernels]
[dg_diff]
type = DGDiffusion
variable = u
epsilon = -1
sigma = 6
diff = diffusion
[]
[]
[Materials]
[coupled_mat]
type = VarCouplingMaterial
var = u
declare_old = true
use_tag = false
[]
[]
[BCs]
[all]
type = DGFunctionDiffusionDirichletBC
variable = u
boundary = '0 1 2 3'
function = exact_fn
epsilon = -1
sigma = 6
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]