- 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
- valueValue of the BC
C++ Type:double
Unit:(no unit assumed)
Controllable:Yes
Description:Value of the BC
- 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
KokkosDirichletBC
This is the Kokkos version of DirichletBC. See the original document for details.
Example Input Syntax
[left_u]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 1
[](test/tests/kokkos/bcs/coupled_dirichlet_bc/kokkos_coupled_dirichlet_bc.i)Input Parameters
- 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)
- presetTrueWhether or not to preset the BC (apply the value before the solve begins).
Default:True
C++ Type:bool
Controllable:No
Description:Whether or not to preset the BC (apply the value before the solve begins).
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.
- 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
Input Files
- (test/tests/kokkos/nodalkernels/reaction/kokkos_reaction.i)
- (test/tests/kokkos/bcs/misc_bcs/kokkos_vacuum_bc_test.i)
- (test/tests/kokkos/kernels/2d_diffusion/kokkos_2d_diffusion_test.i)
- (test/tests/kokkos/nodalkernels/constraint_enforcement/kokkos_lower_bound.i)
- (test/tests/kokkos/nodalkernels/jac_test/kokkos_block_jacobian_test.i)
- (test/tests/kokkos/tag/kokkos_2d_diffusion_vector_tag_test.i)
- (test/tests/kokkos/tag/kokkos_tag_neumann.i)
- (test/tests/kokkos/auxkernels/pp_depend/kokkos_pp_depend.i)
- (test/tests/kokkos/bcs/matched_value_bc/kokkos_matched_value_bc_test.i)
- (test/tests/kokkos/functions/piecewise_constant/kokkos_piecewise_constant.i)
- (test/tests/kokkos/materials/error/kokkos_multid.i)
- (test/tests/kokkos/bcs/coupled_var_neumann/kokkos_coupled_var_neumann.i)
- (test/tests/kokkos/kernels/2d_diffusion/kokkos_2d_diffusion_bodyforce_test.i)
- (test/tests/kokkos/materials/stateful_prop/kokkos_stateful_prop_on_bnd_only.i)
- (test/tests/kokkos/materials/coupling/kokkos_var_coupling.i)
- (modules/heat_transfer/test/tests/kokkos/kokkos_conduction.i)
- (test/tests/kokkos/nodalkernels/jac_test/kokkos_bc_jacobian_test.i)
- (test/tests/kokkos/nodalkernels/constraint_enforcement/kokkos_upper_and_lower_bound.i)
- (test/tests/kokkos/functions/concrete_type/kokkos_concrete_function.i)
- (test/tests/kokkos/materials/stateful_prop/kokkos_stateful_prop_test.i)
- (test/tests/kokkos/kernels/2d_diffusion/kokkos_2d_diffusion_neumannbc_test.i)
- (test/tests/kokkos/restart/stateful/kokkos_stateful_prop_spatial_test_restart.i)
- (test/tests/kokkos/coord_type/kokkos_diffusion_rspherical.i)
- (test/tests/kokkos/restart/stateful/kokkos_stateful_prop_spatial_test.i)
- (test/tests/kokkos/nodalkernels/constraint_enforcement/kokkos_upper_bound.i)
- (test/tests/kokkos/kernels/block_kernel/kokkos_block_vars.i)
- (test/tests/kokkos/materials/stateful_prop/kokkos_stateful_prop_spatial_test.i)
- (test/tests/kokkos/tag/kokkos_2d_diffusion_tag_vector.i)
- (test/tests/kokkos/bcs/misc_bcs/kokkos_convective_flux_bc.i)
- (test/tests/kokkos/bcs/pp_neumann/kokkos_pp_neumann.i)
- (test/tests/kokkos/tag/kokkos_2d_diffusion_tag_matrix.i)
- (test/tests/kokkos/petsc_gpu/kokkos_2d_diffusion_tag_vector.i)
- (test/tests/kokkos/bcs/1d_neumann/kokkos_1d_neumann.i)
- (test/tests/kokkos/materials/coupling/kokkos_prop_stateful_coupling.i)
- (modules/heat_transfer/test/tests/kokkos/kokkos_convective_heat_flux.i)
- (test/tests/kokkos/restart/kernel_restartable/kokkos_kernel_restartable.i)
- (test/tests/kokkos/functions/constant_function/kokkos_constant_function.i)
- (test/tests/kokkos/coord_type/kokkos_diffusion_rz.i)
- (test/tests/kokkos/bcs/coupled_var_neumann/kokkos_on_off.i)
- (test/tests/kokkos/tag/kokkos_tag_nodal_kernels.i)
- (test/tests/kokkos/restart/kernel_restartable/kokkos_kernel_restartable_second.i)
- (test/tests/kokkos/bcs/coupled_dirichlet_bc/kokkos_coupled_dirichlet_bc.i)
- (test/tests/kokkos/kernels/kernel_precompute/kokkos_kernel_precompute_test.i)
- (test/tests/kokkos/nodalkernels/constant_rate/kokkos_constant_rate.i)
- (test/tests/kokkos/materials/error/kokkos_overlap.i)
- (test/tests/kokkos/kernels/block_kernel/kokkos_block_kernel_test.i)
- (test/tests/kokkos/functions/default_function/kokkos_default_function.i)
- (test/tests/kokkos/kernels/coupled_time_derivative/kokkos_coupled_time_derivative_test.i)
- (test/tests/kokkos/bcs/coupled_var_neumann/kokkos_coupled_var_neumann_nl.i)
- (test/tests/kokkos/tag/kokkos_2d_diffusion_matrix_tag_test.i)
- (test/tests/kokkos/materials/error/kokkos_dependency.i)
Child Objects
(test/tests/kokkos/bcs/coupled_dirichlet_bc/kokkos_coupled_dirichlet_bc.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[v]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff_u]
type = KokkosDiffusion
variable = u
[]
[coupled_force_u]
type = KokkosCoupledForce
variable = u
v = v
[]
[diff_v]
type = KokkosDiffusion
variable = v
[]
[]
[BCs]
# BCs on left
# u: u=1
# v: v=2
[left_u]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 1
[]
[left_v]
type = KokkosDirichletBC
variable = v
boundary = 3
value = 2
[]
# BCs on right
# u: c*u + u^2 + v^2 = 9
# v: no flux
[right_u]
type = KokkosCoupledDirichletBC
variable = u
boundary = 1
value = 9
v=v
[]
[]
[Preconditioning]
[precond]
type = SMP
# 'full = true' is required for computeOffDiagJacobian() to get
# called. If you comment this out, you should see that this test
# requires a different number of linear and nonlinear iterations.
full = true
[]
[]
[Executioner]
type = Steady
# solve_type = 'PJFNK'
solve_type = 'NEWTON'
# Uncomment next line to disable line search. With line search enabled, you must use full=true with Newton or else it will fail.
# line_search = 'none'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
nl_rel_tol = 1e-10
l_tol = 1e-12
nl_max_its = 10
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/nodalkernels/reaction/kokkos_reaction.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[nodal_ode]
[]
[]
[Kernels]
[diff]
type = KokkosCoefDiffusion
variable = u
coef = 0.1
[]
[time]
type = KokkosTimeDerivative
variable = u
[]
[]
[NodalKernels]
[td]
type = KokkosTimeDerivativeNodalKernel
variable = nodal_ode
[]
[reaction]
type = KokkosReactionNodalKernel
variable = nodal_ode
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = right
value = 1
[]
[]
[Executioner]
type = Transient
num_steps = 20
dt = 0.1
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/bcs/misc_bcs/kokkos_vacuum_bc_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
nz = 0
zmin = 0
zmax = 0
elem_type = QUAD4
[]
[Variables]
active = 'u'
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
active = 'diff'
[diff]
type = KokkosDiffusion
variable = u
[]
[]
[BCs]
active = 'left right top'
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0.0
[]
[right]
type = KokkosNeumannBC
variable = u
boundary = 1
value = 2.0
[]
[top]
type = KokkosVacuumBC
variable = u
boundary = 2
alpha = 5.0
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/kernels/2d_diffusion/kokkos_2d_diffusion_test.i)
###########################################################
# This is a simple test of the KokkosKernel System.
# It solves the Laplacian equation on a small 2x2 grid.
# The "KokkosDiffusion" kernel is used to calculate the
# residuals of the weak form of this operator.
#
# @Requirement F3.30
###########################################################
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[]
[]
[Variables]
active = 'u'
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
active = 'diff'
[diff]
type = KokkosDiffusion
variable = u
[]
[]
[BCs]
# BCs cannot be preset due to Jacobian test
active = 'left right'
[left]
type = KokkosDirichletBC
variable = u
preset = false
boundary = 3
value = 0
[]
[right]
type = KokkosDirichletBC
variable = u
preset = false
boundary = 1
value = 1
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/nodalkernels/constraint_enforcement/kokkos_lower_bound.i)
l=10
nx=100
num_steps=10
[Mesh]
type = GeneratedMesh
dim = 1
xmax = ${l}
nx = ${nx}
[]
[Variables]
[u]
[]
[lm]
[]
[]
[ICs]
[u]
type = FunctionIC
variable = u
function = '${l} - x'
[]
[]
[Kernels]
[time]
type = KokkosTimeDerivative
variable = u
[]
[diff]
type = KokkosDiffusion
variable = u
[]
[ffn]
type = KokkosBodyForce
variable = u
value = -1
[]
[]
[NodalKernels]
[positive_constraint]
type = KokkosLowerBoundNodalKernel
variable = lm
v = u
exclude_boundaries = 'left right'
[]
[forces]
type = KokkosCoupledForceNodalKernel
variable = u
v = lm
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
boundary = left
value = ${l}
variable = u
[]
[right]
type = KokkosDirichletBC
boundary = right
value = 0
variable = u
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
num_steps = ${num_steps}
solve_type = NEWTON
dtmin = 1
petsc_options_iname = '-snes_max_linear_solve_fail -ksp_max_it -pc_type -sub_pc_factor_levels -snes_linesearch_type'
petsc_options_value = '0 30 asm 16 basic'
[]
[Outputs]
exodus = true
[csv]
type = CSV
execute_on = 'nonlinear timestep_end'
[]
[dof]
type = DOFMap
execute_on = 'initial'
[]
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[active_lm]
type = GreaterThanLessThanPostprocessor
variable = lm
execute_on = 'nonlinear timestep_end'
value = 1e-8
[]
[violations]
type = GreaterThanLessThanPostprocessor
variable = u
execute_on = 'nonlinear timestep_end'
value = -1e-8
comparator = 'less'
[]
[]
(test/tests/kokkos/nodalkernels/jac_test/kokkos_block_jacobian_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u_x]
[]
[u_y]
[]
[]
[Kernels]
[diff_x]
type = KokkosCoefDiffusion
variable = u_x
coef = 0.1
[]
[diff_y]
type = KokkosCoefDiffusion
variable = u_y
coef = 0.1
[]
[]
[NodalKernels]
[test_y]
type = KokkosJacobianCheck
variable = u_y
[]
[test_x]
type = KokkosJacobianCheck
variable = u_x
[]
[]
[BCs]
[left_x]
type = KokkosDirichletBC
variable = u_x
preset = false
boundary = left
value = 0
[]
[right_x]
type = KokkosDirichletBC
variable = u_x
preset = false
boundary = right
value = 1
[]
[left_y]
type = KokkosDirichletBC
variable = u_y
preset = false
boundary = left
value = 0
[]
[]
[Executioner]
type = Transient
num_steps = 1
dt = 0.1
solve_type = NEWTON
# petsc_options = '-snes_check_jacobian -snes_check_jacobian_view'
nl_max_its = 1
nl_abs_tol = 1e0
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/tag/kokkos_2d_diffusion_vector_tag_test.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[]
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
initial_condition = 1.0
[]
[]
[AuxVariables]
[tag_variable1]
order = FIRST
family = LAGRANGE
[]
[tag_variable2]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[reaction1]
type = KokkosReaction
variable = u
extra_vector_tags = 'vec_tag1 vec_tag2'
[]
[reaction2]
type = KokkosReaction
variable = u
extra_vector_tags = 'vec_tag1 vec_tag2'
[]
[reaction3]
type = KokkosReaction
variable = u
[]
[reaction4]
type = KokkosReaction
variable = u
[]
[]
[AuxKernels]
[TagVectorAux1]
type = KokkosTagVectorAux
variable = tag_variable1
v = u
vector_tag = vec_tag1
[]
[TagVectorAux2]
type = KokkosTagVectorAux
variable = tag_variable2
v = u
vector_tag = vec_tag2
[]
[]
[BCs]
active = 'left right'
[left]
type = KokkosDirichletBC
variable = u
preset = false
boundary = 3
value = 10
extra_vector_tags = vec_tag1
[]
[right]
type = KokkosDirichletBC
variable = u
preset = false
boundary = 1
value = 100
extra_vector_tags = vec_tag2
[]
[right1]
type = KokkosDirichletBC
variable = u
preset = false
boundary = 1
value = 100
[]
[right2]
type = KokkosDirichletBC
variable = u
preset = false
boundary = 1
value = 100
[]
[]
[Problem]
type = TagTestProblem
extra_tag_vectors = 'vec_tag1 vec_tag2'
test_tag_vectors = 'vec_tag1 vec_tag2'
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/tag/kokkos_tag_neumann.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
construct_side_list_from_node_list = true
[]
[Variables]
[u]
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[]
[]
[AuxVariables]
[tag_variable1]
order = FIRST
family = LAGRANGE
[]
[tag_variable2]
order = FIRST
family = LAGRANGE
[]
[]
[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
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[]
[right]
type = KokkosNeumannBC
variable = u
boundary = right
value = 2
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 = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/auxkernels/pp_depend/kokkos_pp_depend.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[AuxVariables]
[pp_aux]
[]
[]
[Functions]
[t_func]
type = ParsedFunction
expression = t
[]
[]
[Kernels]
[diff]
type = KokkosCoefDiffusion
variable = u
coef = 0.01
[]
[]
[AuxKernels]
[pp_aux]
type = KokkosPostprocessorAux
variable = pp_aux
execute_on = timestep_end
pp = t_pp
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = right
value = 1
[]
[]
[Postprocessors]
[t_pp]
type = FunctionValuePostprocessor
function = t_func
[]
[]
[Problem]
type = FEProblem
solve = false
[]
[Executioner]
type = Transient
solve_type = PJFNK
dt = 1
num_steps = 5
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/bcs/matched_value_bc/kokkos_matched_value_bc_test.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[]
[]
# Solves a pair of coupled diffusion equations where u=v on the boundary
[Variables]
active = 'u v'
[u]
order = FIRST
family = LAGRANGE
initial_condition = 3
[]
[v]
order = FIRST
family = LAGRANGE
initial_condition = 2
[]
[]
[Kernels]
active = 'diff_u diff_v'
[diff_u]
type = KokkosDiffusion
variable = u
[]
[diff_v]
type = KokkosDiffusion
variable = v
[]
[]
[BCs]
active = 'right_v left_u'
[right_v]
type = KokkosDirichletBC
variable = v
boundary = 1
value = 3
[]
[left_u]
type = KokkosMatchedValueBC
variable = u
boundary = 3
v = v
[]
[]
[Preconditioning]
[precond]
type = SMP
full = true
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
nl_rel_tol = 1e-10
l_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/functions/piecewise_constant/kokkos_piecewise_constant.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[AuxVariables]
[coef]
[]
[]
[UserObjects]
[json]
type = JSONFileReader
filename = 'xy.json'
[]
[]
[Functions]
[func_x_y]
type = KokkosPiecewiseConstant
x = '0.0 0.5'
y = '2.0 3.0'
[]
[func_xy_data]
type = KokkosPiecewiseConstant
xy_data = '0.0 2.0
0.5 3.0'
[]
[func_csv]
type = KokkosPiecewiseConstant
data_file = xy.csv
[]
[func_json]
type = KokkosPiecewiseConstant
json_uo = json
x_keys = 'data x value'
y_keys = 'data y value'
[]
[]
[Kernels]
[diff]
type = KokkosFuncCoefDiffusion
variable = u
coef = func_x_y
[]
[time]
type = KokkosTimeDerivative
variable = u
[]
[]
[AuxKernels]
[coef_aux]
type = KokkosFunctionAux
variable = coef
function = func_x_y
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = KokkosNeumannBC
variable = u
boundary = right
value = 1
[]
[]
[Postprocessors]
[coef]
type = ElementIntegralVariablePostprocessor
variable = coef
[]
[]
[Executioner]
type = Transient
num_steps = 10
dt = 0.1
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/materials/error/kokkos_multid.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 1
[]
[left_domain]
input = gen
type = SubdomainBoundingBoxGenerator
bottom_left = '0 0 0'
top_right = '0.5 1 0'
block_id = 10
[]
[]
[Variables]
[u]
initial_condition = 2
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 2
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = right
value = 3
[]
[]
[Materials]
[left_real_1D]
type = Kokkos1DRealProperty
block = 0
name = 'prop_a'
dims = '1'
[]
[right_real_1D]
type = Kokkos1DRealProperty
block = 10
name = 'prop_a'
dims = '2'
[]
[right_real_1D_correct]
type = Kokkos1DRealProperty
block = 10
name = 'prop_a'
dims = '1'
[]
[right_real_2D]
type = Kokkos2DRealProperty
block = 10
name = 'prop_a'
dims = '1 1'
[]
[right_int_1D]
type = Kokkos1DIntProperty
block = 10
name = 'prop_a'
dims = '1'
[]
[]
[Executioner]
type = Steady
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
print_linear_residuals = true
perf_graph = true
[]
(test/tests/kokkos/bcs/coupled_var_neumann/kokkos_coupled_var_neumann.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 10
ny = 10
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
[]
[]
[AuxVariables]
[coupled_bc_var]
[]
[]
[ICs]
[coupled_bc_var]
type = FunctionIC
variable = coupled_bc_var
function = set_coupled_bc_var
[]
[]
[Functions]
[set_coupled_bc_var]
type = ParsedFunction
expression = 'y - 0.5'
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0
[]
[right]
type = KokkosCoupledVarNeumannBC
variable = u
boundary = 1
v = coupled_bc_var
[]
[]
[Executioner]
type = Steady
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/kernels/2d_diffusion/kokkos_2d_diffusion_bodyforce_test.i)
###########################################################
# This is a simple test of the KokkosKernel System.
# It solves the Laplacian equation on a small 2x2 grid.
# The "KokkosDiffusion" kernel is used to calculate the
# residuals of the weak form of this operator. The
# "KokkosBodyForce" kernel is used to apply a time-dependent
# volumetric source.
###########################################################
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[]
[]
[Variables]
active = 'u'
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
[]
[bf]
type = KokkosBodyForce
variable = u
postprocessor = ramp
[]
[]
[Functions]
[ramp]
type = ParsedFunction
expression = 't'
[]
[]
[Postprocessors]
[ramp]
type = FunctionValuePostprocessor
function = ramp
execute_on = linear
[]
[]
[BCs]
active = 'left right'
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = 1
value = 0
[]
[]
[Executioner]
type = Transient
dt = 1.0
end_time = 1.0
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/materials/stateful_prop/kokkos_stateful_prop_on_bnd_only.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = 10
nx = 10
ny = 10
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[heat]
type = KokkosMatDiffusionTest
variable = u
prop_name = thermal_conductivity
[]
[ie]
type = KokkosTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0.0
[]
[right]
type = KokkosMTBC
variable = u
boundary = right
grad = 1.0
prop_name = thermal_conductivity
[]
[]
[Materials]
[volatile]
type = KokkosGenericConstantMaterial
prop_names = 'thermal_conductivity'
prop_values = 10
block = 0
[]
[stateful_on_boundary]
type = KokkosStatefulSpatialTest
boundary = right
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
start_time = 0.0
num_steps = 5
dt = .1
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/materials/coupling/kokkos_var_coupling.i)
# The diffusion coefficient is given as the variable itself, which makes it a non-linear problem.
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
initial_condition = 1
[]
[]
[AuxVariables]
[v]
initial_condition = 1
[]
[]
[Kernels]
[diff]
type = KokkosMatDiffusionTest
variable = u
prop_name = diffusion
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = right
value = 1
[]
[]
[Materials]
[coupling_u]
type = KokkosVarCouplingMaterial
block = 0
var = u
[]
[]
[Executioner]
type = Transient
num_steps = 1
dt = 1
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(modules/heat_transfer/test/tests/kokkos/kokkos_conduction.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[temp]
initial_condition = 100.0
[]
[]
[Kernels]
[heat_conduction]
type = KokkosHeatConduction
variable = temp
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = temp
boundary = left
value = 100.0
[]
[right]
type = KokkosDirichletBC
variable = temp
boundary = right
value = 200.0
[]
[]
[Materials]
[conduction]
type = KokkosHeatConductionMaterial
thermal_conductivity = 10
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/nodalkernels/jac_test/kokkos_bc_jacobian_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 4
ny = 4
[]
[Variables]
[u_x]
[]
[u_y]
[]
[]
[Kernels]
[diff_x]
type = KokkosCoefDiffusion
variable = u_x
coef = 0.1
[]
[diff_y]
type = KokkosCoefDiffusion
variable = u_y
coef = 0.1
[]
[]
[NodalKernels]
[test_y]
type = KokkosJacobianCheck
variable = u_y
boundary = top
[]
[test_x]
type = KokkosJacobianCheck
variable = u_x
boundary = top
[]
[]
[BCs]
[left_x]
type = KokkosDirichletBC
variable = u_x
preset = false
boundary = left
value = 0
[]
[right_x]
type = KokkosDirichletBC
variable = u_x
preset = false
boundary = right
value = 1
[]
[left_y]
type = KokkosDirichletBC
variable = u_y
preset = false
boundary = left
value = 0
[]
[]
[Executioner]
type = Transient
num_steps = 1
dt = 0.1
solve_type = NEWTON
# petsc_options = '-snes_check_jacobian -snes_check_jacobian_view'
nl_max_its = 1
nl_abs_tol = 1e0
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/nodalkernels/constraint_enforcement/kokkos_upper_and_lower_bound.i)
l=10
nx=100
num_steps=10
[Mesh]
type = GeneratedMesh
dim = 1
xmax = ${l}
nx = ${nx}
[]
[Variables]
[u]
[]
[lm_upper]
[]
[lm_lower]
[]
[]
[AuxVariables]
[force]
family = MONOMIAL
order = CONSTANT
[]
[]
[ICs]
[u]
type = FunctionIC
variable = u
function = 'x'
[]
[force]
type = FunctionIC
variable = force
function = 'if(x<5,-1,1)'
[]
[]
[Kernels]
[time]
type = KokkosTimeDerivative
variable = u
[]
[diff]
type = KokkosDiffusion
variable = u
[]
[ffn]
type = KokkosCoupledForce
variable = u
v = force
[]
[]
[NodalKernels]
[upper_bound]
type = KokkosUpperBoundNodalKernel
variable = lm_upper
v = u
exclude_boundaries = 'left right'
upper_bound = 10
[]
[forces_from_upper]
type = KokkosCoupledForceNodalKernel
variable = u
v = lm_upper
coef = -1
[]
[lower_bound]
type = KokkosLowerBoundNodalKernel
variable = lm_lower
v = u
exclude_boundaries = 'left right'
lower_bound = 0
[]
[forces_from_lower]
type = KokkosCoupledForceNodalKernel
variable = u
v = lm_lower
coef = 1
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
boundary = left
value = 0
variable = u
[]
[right]
type = KokkosDirichletBC
boundary = right
value = ${l}
variable = u
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
num_steps = ${num_steps}
solve_type = NEWTON
dtmin = 1
petsc_options_iname = '-snes_max_linear_solve_fail -ksp_max_it -pc_type -sub_pc_factor_levels -snes_linesearch_type'
petsc_options_value = '0 30 asm 16 basic'
[]
[Outputs]
exodus = true
hide = force
[csv]
type = CSV
execute_on = 'nonlinear timestep_end'
[]
[dof]
type = DOFMap
execute_on = 'initial'
[]
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[active_upper_lm]
type = GreaterThanLessThanPostprocessor
variable = lm_upper
execute_on = 'nonlinear timestep_end'
value = 1e-8
comparator = 'greater'
[]
[upper_violations]
type = GreaterThanLessThanPostprocessor
variable = u
execute_on = 'nonlinear timestep_end'
value = ${fparse 10+1e-8}
comparator = 'greater'
[]
[active_lower_lm]
type = GreaterThanLessThanPostprocessor
variable = lm_lower
execute_on = 'nonlinear timestep_end'
value = 1e-8
comparator = 'greater'
[]
[lower_violations]
type = GreaterThanLessThanPostprocessor
variable = u
execute_on = 'nonlinear timestep_end'
value = -1e-8
comparator = 'less'
[]
[nls]
type = NumNonlinearIterations
[]
[cum_nls]
type = CumulativeValuePostprocessor
postprocessor = nls
[]
[]
(test/tests/kokkos/functions/concrete_type/kokkos_concrete_function.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[Functions]
[constant]
type = KokkosConstantFunction
value = 2
[]
[piecewise_constant]
type = KokkosPiecewiseConstant
x = '0 0.5 1'
y = '1 2 3'
[]
[]
[Kernels]
[diff]
type = KokkosConstantFuncCoefDiffusion
variable = u
coef = constant
[]
[time]
type = KokkosTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = KokkosNeumannBC
variable = u
boundary = right
value = 1
[]
[]
[Executioner]
type = Transient
num_steps = 10
dt = 0.1
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/materials/stateful_prop/kokkos_stateful_prop_test.i)
[Mesh]
dim = 3
file = cube.e
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[]
[AuxVariables]
[prop1]
order = CONSTANT
family = MONOMIAL
[]
[]
[Kernels]
[heat]
type = KokkosMatDiffusionTest
variable = u
prop_name = thermal_conductivity
prop_state = 'older' # Use the "Older" value to compute conductivity
[]
[ie]
type = KokkosTimeDerivative
variable = u
[]
[]
[AuxKernels]
[prop1_output_init]
type = KokkosMaterialRealAux
variable = prop1
property = thermal_conductivity
execute_on = initial
[]
[prop1_output]
type = KokkosMaterialRealAux
variable = prop1
property = thermal_conductivity
[]
[]
[BCs]
[bottom]
type = KokkosDirichletBC
variable = u
boundary = 1
value = 0.0
[]
[top]
type = KokkosDirichletBC
variable = u
boundary = 2
value = 1.0
[]
[]
[Materials]
[stateful]
type = KokkosStatefulTest
prop_names = thermal_conductivity
prop_values = 1.0
[]
[]
[Postprocessors]
[integral]
type = ElementAverageValue
variable = prop1
execute_on = 'initial timestep_end'
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
start_time = 0.0
num_steps = 5
dt = .1
[]
[Outputs]
exodus = true
csv = true
[]
(test/tests/kokkos/kernels/2d_diffusion/kokkos_2d_diffusion_neumannbc_test.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[]
[]
[Variables]
active = 'u'
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
active = 'diff'
[diff]
type = KokkosDiffusion
variable = u
[]
[]
[BCs]
active = 'left right'
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0
[]
[right]
type = KokkosNeumannBC
variable = u
boundary = 1
value = 1
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/restart/stateful/kokkos_stateful_prop_spatial_test_restart.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = 10
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[AuxVariables]
[prop1]
order = SECOND
family = MONOMIAL
[]
[]
[AuxKernels]
[prop1_output]
type = KokkosMaterialRealAux
variable = prop1
property = thermal_conductivity
[]
[]
[Kernels]
[heat]
type = KokkosMatDiffusionTest
variable = u
prop_name = thermal_conductivity
[]
[ie]
type = KokkosTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0.0
[]
[right]
type = KokkosMTBC
variable = u
boundary = 1
grad = 1.0
prop_name = thermal_conductivity
[]
[]
[Materials]
active = stateful
[stateful]
type = KokkosStatefulSpatialTest
block = 0
[]
[another_stateful]
type = KokkosStatefulSpatialTest
block = 0
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
start_time = 0.3
num_steps = 2
dt = .1
[]
[Outputs]
[out]
type = Exodus
elemental_as_nodal = true
execute_elemental_on = none
[]
[]
[Problem]
restart_file_base = kokkos_stateful_prop_spatial_test_out_cp/LATEST
[]
(test/tests/kokkos/coord_type/kokkos_diffusion_rspherical.i)
[Mesh]
[sphere]
type = GeneratedMeshGenerator
nx = 10
dim = 1
[]
coord_type = RSPHERICAL
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 0
value = 1
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = 1
value = 0
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/restart/stateful/kokkos_stateful_prop_spatial_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = 10
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[AuxVariables]
[prop1]
order = SECOND
family = MONOMIAL
[]
[]
[AuxKernels]
[prop1_output]
type = KokkosMaterialRealAux
variable = prop1
property = thermal_conductivity
[]
[]
[Kernels]
[heat]
type = KokkosMatDiffusionTest
variable = u
prop_name = thermal_conductivity
[]
[ie]
type = KokkosTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0.0
[]
[right]
type = KokkosMTBC
variable = u
boundary = 1
grad = 1.0
prop_name = thermal_conductivity
[]
[]
[Materials]
[stateful]
type = KokkosStatefulSpatialTest
block = 0
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
start_time = 0.0
num_steps = 3
dt = .1
[]
[Outputs]
checkpoint = true
[]
(test/tests/kokkos/nodalkernels/constraint_enforcement/kokkos_upper_bound.i)
l=10
nx=100
num_steps=10
[Mesh]
type = GeneratedMesh
dim = 1
xmax = ${l}
nx = ${nx}
[]
[Variables]
[u]
[]
[lm]
[]
[]
[ICs]
[u]
type = FunctionIC
variable = u
function = '${l} - x'
[]
[]
[Kernels]
[time]
type = KokkosTimeDerivative
variable = u
[]
[diff]
type = KokkosDiffusion
variable = u
[]
[ffn]
type = KokkosBodyForce
variable = u
[]
[]
[NodalKernels]
[positive_constraint]
type = KokkosUpperBoundNodalKernel
variable = lm
v = u
exclude_boundaries = 'left right'
upper_bound = 10
[]
[forces]
type = KokkosCoupledForceNodalKernel
variable = u
v = lm
coef = -1
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
boundary = left
value = ${l}
variable = u
[]
[right]
type = KokkosDirichletBC
boundary = right
value = 0
variable = u
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
num_steps = ${num_steps}
solve_type = NEWTON
dtmin = 1
petsc_options_iname = '-snes_max_linear_solve_fail -ksp_max_it -pc_type -sub_pc_factor_levels -snes_linesearch_type'
petsc_options_value = '0 30 asm 16 basic'
[]
[Outputs]
exodus = true
[csv]
type = CSV
execute_on = 'nonlinear timestep_end'
[]
[dof]
type = DOFMap
execute_on = 'initial'
[]
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[active_lm]
type = GreaterThanLessThanPostprocessor
variable = lm
execute_on = 'nonlinear timestep_end'
value = 1e-8
[]
[violations]
type = GreaterThanLessThanPostprocessor
variable = u
execute_on = 'nonlinear timestep_end'
value = ${fparse 10+1e-8}
comparator = 'greater'
[]
[]
(test/tests/kokkos/kernels/block_kernel/kokkos_block_vars.i)
[Mesh]
file = rect-2blk.e
[]
[Variables]
active = 'u v'
[u]
order = FIRST
family = LAGRANGE
block = 1
[]
[v]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
active = 'diff_u diff_v'
[diff_u]
type = KokkosDiffusion
variable = u
[]
[diff_v]
type = KokkosDiffusion
variable = v
[]
[]
[BCs]
active = 'left_u right_u left_v right_v'
[left_u]
type = KokkosDirichletBC
variable = u
boundary = 6
value = 0
[]
[right_u]
type = KokkosNeumannBC
variable = u
boundary = 8
value = 4
[]
[left_v]
type = KokkosDirichletBC
variable = v
boundary = 6
value = 1
[]
[right_v]
type = KokkosDirichletBC
variable = v
boundary = 3
value = 6
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Outputs]
execute_on = 'timestep_end'
exodus = true
[]
(test/tests/kokkos/materials/stateful_prop/kokkos_stateful_prop_spatial_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = 10
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[AuxVariables]
[prop1]
order = SECOND
family = MONOMIAL
[]
[]
[AuxKernels]
[prop1_output]
type = KokkosMaterialRealAux
variable = prop1
property = thermal_conductivity
[]
[]
[Kernels]
[heat]
type = KokkosMatDiffusionTest
variable = u
prop_name = thermal_conductivity
[]
[ie]
type = KokkosTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0.0
[]
[right]
type = KokkosMTBC
variable = u
boundary = 1
grad = 1.0
prop_name = thermal_conductivity
[]
[]
[Materials]
[stateful]
type = KokkosStatefulSpatialTest
block = 0
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
start_time = 0.0
num_steps = 5
dt = .1
[]
[Outputs]
[out]
type = Exodus
elemental_as_nodal = true
execute_elemental_on = none
[]
[]
(test/tests/kokkos/tag/kokkos_2d_diffusion_tag_vector.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[]
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[]
[AuxVariables]
[tag_variable1]
order = FIRST
family = LAGRANGE
[]
[tag_variable2]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
extra_vector_tags = 'vec_tag1 vec_tag2'
[]
[]
[AuxKernels]
[TagVectorAux1]
type = KokkosTagVectorAux
variable = tag_variable1
v = u
vector_tag = vec_tag1
[]
[TagVectorAux2]
type = KokkosTagVectorAux
variable = tag_variable2
v = u
vector_tag = vec_tag2
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0
preset = false
extra_vector_tags = vec_tag1
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = 1
value = 1
preset = false
extra_vector_tags = vec_tag2
[]
[]
[Problem]
type = FEProblem
extra_tag_vectors = 'vec_tag1 vec_tag2'
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/bcs/misc_bcs/kokkos_convective_flux_bc.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
nz = 0
zmin = 0
zmax = 0
elem_type = QUAD4
[]
[Variables]
active = 'u'
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
active = 'diff'
[diff]
type = KokkosDiffusion
variable = u
[]
[]
[BCs]
active = 'left right'
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0.0
[]
[right]
type = KokkosConvectiveFluxBC
variable = u
boundary = 1
rate = 100
initial = 10
final = 20
duration = 10
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
start_time = 0.0
num_steps = 10
dt = 1.0
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/bcs/pp_neumann/kokkos_pp_neumann.i)
# NOTE: This file is used within the documentation, so please do not change names within the file
# without checking that associated documentation is not affected, see syntax/Postprocessors/index.md.
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[AuxVariables]
[aux]
initial_condition = 5
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = KokkosPostprocessorNeumannBC
variable = u
boundary = right
postprocessor = right_pp
[]
[]
[Postprocessors]
[right_pp]
type = PointValue
point = '0.5 0.5 0'
variable = aux
execute_on = 'initial'
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
execute_on = 'timestep_end'
exodus = true
[]
(test/tests/kokkos/tag/kokkos_2d_diffusion_tag_matrix.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[]
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[]
[AuxVariables]
[tag_variable1]
order = FIRST
family = LAGRANGE
[]
[tag_variable2]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
[]
[]
[AuxKernels]
[TagMatrixAux1]
type = TagMatrixAux
variable = tag_variable1
v = u
matrix_tag = mat_tag1
[]
[TagMatrixAux2]
type = TagMatrixAux
variable = tag_variable2
v = u
matrix_tag = mat_tag2
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0
extra_matrix_tags = mat_tag1
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = 1
value = 1
extra_matrix_tags = mat_tag2
[]
[]
[Problem]
type = FEProblem
extra_tag_matrices = 'mat_tag1 mat_tag2'
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/petsc_gpu/kokkos_2d_diffusion_tag_vector.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[]
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[]
[AuxVariables]
[tag_variable1]
order = FIRST
family = LAGRANGE
[]
[tag_variable2]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
extra_vector_tags = 'vec_tag1 vec_tag2'
[]
[]
[AuxKernels]
[TagVectorAux1]
type = KokkosTagVectorAux
variable = tag_variable1
v = u
vector_tag = vec_tag1
[]
[TagVectorAux2]
type = KokkosTagVectorAux
variable = tag_variable2
v = u
vector_tag = vec_tag2
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0
extra_vector_tags = vec_tag1
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = 1
value = 1
extra_vector_tags = vec_tag2
[]
[]
[Problem]
type = FEProblem
extra_tag_vectors = 'vec_tag1 vec_tag2'
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
petsc_options_iname = '-vec_type -nl0_mat_type'
petsc_options_value = 'kokkos hypre'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/bcs/1d_neumann/kokkos_1d_neumann.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
construct_side_list_from_node_list = true
[]
[Variables]
[u]
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = KokkosNeumannBC
variable = u
boundary = right
value = 2
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/materials/coupling/kokkos_prop_stateful_coupling.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
active = 'u'
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = KokkosMatDiffusionTest
variable = u
prop_name = 'some_prop'
[]
[time]
type = KokkosTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 'left'
value = 1
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = 'right'
value = 2
[]
[]
[Materials]
# This material couples in a stateful property from StatefulTest
[coupled_mat]
type = KokkosCoupledMaterial
mat_prop = 'some_prop'
coupled_mat_prop = 'thermal_conductivity'
use_old_prop = true
[]
[stateful_mat]
type = KokkosStatefulTest
prop_names = thermal_conductivity
prop_values = 1.0
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
num_steps = 4
[]
[Outputs]
exodus = true
[]
[Debug]
show_material_props = true
[]
(modules/heat_transfer/test/tests/kokkos/kokkos_convective_heat_flux.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[temp]
initial_condition = 100.0
[]
[]
[Kernels]
[heat_conduction]
type = KokkosHeatConduction
variable = temp
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = temp
boundary = left
value = 100.0
[]
[right]
type = KokkosConvectiveHeatFluxBC
variable = temp
boundary = right
T_infinity = 200.0
heat_transfer_coefficient = 10
[]
[]
[Materials]
[conduction]
type = KokkosHeatConductionMaterial
thermal_conductivity = 10
[]
[bc]
type = KokkosGenericConstantMaterial
prop_names = 'Tinf htc'
prop_values = '200 10'
[]
[]
[Postprocessors]
[right_flux]
type = SideDiffusiveFluxAverage
variable = temp
boundary = right
diffusivity = 10
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/restart/kernel_restartable/kokkos_kernel_restartable.i)
###########################################################
# This test exercises the restart system and verifies
# correctness with parallel computation, but distributed
# and with threading.
#
# See kernel_restartable_second.i
#
# @Requirement F1.60
# @Requirement P1.10
# @Requirement P1.20
###########################################################
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[Kernels]
[diff]
type = KokkosRestartDiffusion
variable = u
[]
[td]
type = KokkosTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = right
value = 1
[]
[]
[Executioner]
type = Transient
num_steps = 5
dt = 1e-2
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[restart]
type = Checkpoint
num_files = 100
[]
[]
(test/tests/kokkos/functions/constant_function/kokkos_constant_function.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[Functions]
[constant]
type = KokkosConstantFunction
value = 2
[]
[]
[Kernels]
[diff]
type = KokkosFuncCoefDiffusion
variable = u
coef = constant
[]
[time]
type = KokkosTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = KokkosNeumannBC
variable = u
boundary = right
value = 1
[]
[]
[Executioner]
type = Transient
num_steps = 10
dt = 0.1
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/coord_type/kokkos_diffusion_rz.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 10
ny = 10
dim = 2
[]
coord_type = RZ
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 1
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = 1
value = 0
[]
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/bcs/coupled_var_neumann/kokkos_on_off.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 10
ny = 10
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
[]
[]
[AuxVariables]
[coupled_bc_var]
[]
[active]
initial_condition = 1
[]
[]
[AuxKernels]
[active_right]
type = KokkosConstantAux
variable = active
value = 0.5
boundary = 1
[]
[]
[ICs]
[coupled_bc_var]
type = FunctionIC
variable = coupled_bc_var
function = set_coupled_bc_var
[]
[]
[Functions]
[set_coupled_bc_var]
type = ParsedFunction
expression = 'y - 0.5'
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0
[]
[right]
type = KokkosCoupledVarNeumannBC
variable = u
boundary = 1
v = coupled_bc_var
scale_factor = active
[]
[]
[Executioner]
type = Steady
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/tag/kokkos_tag_nodal_kernels.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[nodal_ode]
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[]
[time]
type = KokkosTimeDerivative
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[]
[]
[NodalKernels]
[td]
type = KokkosTimeDerivativeNodalKernel
variable = nodal_ode
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[]
[constant_rate]
type = KokkosConstantRate
variable = nodal_ode
rate = 1.0
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1 vec_tag2'
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = right
value = 10
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[]
[]
[Problem]
type = TagTestProblem
test_tag_vectors = 'time 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'
[]
[AuxVariables]
[tag_variable1]
order = FIRST
family = LAGRANGE
[]
[tag_variable2]
order = FIRST
family = LAGRANGE
[]
[]
[AuxKernels]
[TagVectorAux1]
type = TagVectorAux
variable = tag_variable1
v = nodal_ode
vector_tag = vec_tag2
[]
[TagVectorAux2]
type = TagMatrixAux
variable = tag_variable2
v = u
matrix_tag = mat_tag2
[]
[]
[Executioner]
type = Transient
num_steps = 10
nl_rel_tol = 1e-08
dt = 0.01
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/restart/kernel_restartable/kokkos_kernel_restartable_second.i)
###########################################################
# This test exercises the restart system and verifies
# correctness with parallel computation, but distributed
# and with threading.
#
# See kernel_restartable.i
#
# @Requirement F1.60
# @Requirement P1.10
# @Requirement P1.20
###########################################################
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[Kernels]
[diff]
type = KokkosRestartDiffusion
variable = u
[]
[td]
type = KokkosTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = right
value = 1
[]
[]
[Executioner]
type = Transient
num_steps = 5
dt = 1e-2
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
[Problem]
restart_file_base = kokkos_kernel_restartable_restart_cp/LATEST
[]
(test/tests/kokkos/bcs/coupled_dirichlet_bc/kokkos_coupled_dirichlet_bc.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[v]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff_u]
type = KokkosDiffusion
variable = u
[]
[coupled_force_u]
type = KokkosCoupledForce
variable = u
v = v
[]
[diff_v]
type = KokkosDiffusion
variable = v
[]
[]
[BCs]
# BCs on left
# u: u=1
# v: v=2
[left_u]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 1
[]
[left_v]
type = KokkosDirichletBC
variable = v
boundary = 3
value = 2
[]
# BCs on right
# u: c*u + u^2 + v^2 = 9
# v: no flux
[right_u]
type = KokkosCoupledDirichletBC
variable = u
boundary = 1
value = 9
v=v
[]
[]
[Preconditioning]
[precond]
type = SMP
# 'full = true' is required for computeOffDiagJacobian() to get
# called. If you comment this out, you should see that this test
# requires a different number of linear and nonlinear iterations.
full = true
[]
[]
[Executioner]
type = Steady
# solve_type = 'PJFNK'
solve_type = 'NEWTON'
# Uncomment next line to disable line search. With line search enabled, you must use full=true with Newton or else it will fail.
# line_search = 'none'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
nl_rel_tol = 1e-10
l_tol = 1e-12
nl_max_its = 10
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/kernels/kernel_precompute/kokkos_kernel_precompute_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
nz = 0
zmin = 0
zmax = 0
elem_type = QUAD4
[]
[Variables]
active = 'convected'
[convected]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
active = 'diff conv'
[diff]
type = KokkosDiffusionPrecompute
variable = convected
[]
[conv]
type = KokkosConvectionPrecompute
variable = convected
velocity = '1.0 0.0 0.0'
[]
[]
[BCs]
active = 'bottom top'
[bottom]
type = KokkosDirichletBC
variable = convected
boundary = 'left'
value = 0
[]
[top]
type = KokkosDirichletBC
variable = convected
boundary = 'right'
value = 1
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Outputs]
execute_on = 'timestep_end'
exodus = true
[]
(test/tests/kokkos/nodalkernels/constant_rate/kokkos_constant_rate.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[nodal_ode]
[]
[]
[Kernels]
[diff]
type = KokkosCoefDiffusion
variable = u
coef = 0.1
[]
[time]
type = KokkosTimeDerivative
variable = u
[]
[]
[NodalKernels]
[td]
type = KokkosTimeDerivativeNodalKernel
variable = nodal_ode
[]
[constant_rate]
type = KokkosConstantRate
variable = nodal_ode
rate = 1.0
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = right
value = 1
[]
[]
[Executioner]
type = Transient
num_steps = 20
dt = 0.1
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/materials/error/kokkos_overlap.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 10
ny = 1
[]
[left_domain]
input = gen
type = SubdomainBoundingBoxGenerator
bottom_left = '0 0 0'
top_right = '0.5 1 0'
block_id = 10
[]
[]
[Variables]
[u]
initial_condition = 2
[]
[]
[Kernels]
[diff]
type = KokkosMatDiffusionTest
variable = u
prop_name = 'p'
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 2
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = right
value = 3
[]
[]
[Materials]
[all]
type = KokkosGenericConstantMaterial
prop_names = 'f f_prime p'
prop_values = '2 2.5 2.468'
[]
[left]
type = KokkosGenericConstantMaterial
prop_names = 'f f_prime p'
prop_values = '1 0.5 1.2345'
[]
[]
[Executioner]
type = Steady
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
print_linear_residuals = true
perf_graph = true
[]
(test/tests/kokkos/kernels/block_kernel/kokkos_block_kernel_test.i)
[Mesh]
file = rectangle.e
# uniform_refine = 1
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
initial_condition = 1.0
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
[]
[body_force]
# Corresponds to BodyForce with function = 'x + y'
type = KokkosXYBodyForce
variable = u
block = 1
value = 10
[]
[time]
type = KokkosTimeDerivative
variable = u
[]
[]
[BCs]
[right]
type = KokkosDirichletBC
variable = u
boundary = 2
value = 1
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
dt = 0.1
num_steps = 10
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/functions/default_function/kokkos_default_function.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[Kernels]
[diff]
type = KokkosFuncCoefDiffusion
variable = u
# No default function supplied
[]
[time]
type = KokkosTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = KokkosNeumannBC
variable = u
boundary = right
value = 1
[]
[]
[Executioner]
type = Transient
num_steps = 10
dt = 0.1
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/kernels/coupled_time_derivative/kokkos_coupled_time_derivative_test.i)
###########################################################
# This is a simple test of the CoupledTimeDerivative kernel.
# The expected solution for the variable v is
# v(x) = 1/2 * (x^2 + x)
###########################################################
[Mesh]
type = GeneratedMesh
nx = 5
ny = 5
dim = 2
[]
[Variables]
[u]
[]
[v]
[]
[]
[Kernels]
[time_u]
type = KokkosTimeDerivative
variable = u
[]
[fn_u]
type = KokkosBodyForce
variable = u
[]
[time_v]
type = KokkosCoupledTimeDerivative
variable = v
v = u
[]
[diff_v]
type = KokkosDiffusion
variable = v
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = v
boundary = 'left'
value = 0
[]
[right]
type = KokkosDirichletBC
variable = v
boundary = 'right'
value = 1
[]
[]
[Executioner]
type = Transient
num_steps = 1
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/bcs/coupled_var_neumann/kokkos_coupled_var_neumann_nl.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 20
[]
[Variables]
[u][]
[v][]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
[]
[diff_v]
type = KokkosDiffusion
variable = v
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 'left'
value = 0
[]
[right]
type = KokkosCoupledVarNeumannBC
variable = u
boundary = 'right'
v = v
[]
[v_left]
type = KokkosDirichletBC
variable = v
boundary = 'left'
value = 0
[]
[v_right]
type = KokkosDirichletBC
variable = v
boundary = 'right'
value = 1
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Steady
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/tag/kokkos_2d_diffusion_matrix_tag_test.i)
[Mesh]
[square]
type = GeneratedMeshGenerator
nx = 2
ny = 2
dim = 2
[]
[]
[Variables]
[u]
order = FIRST
family = LAGRANGE
[]
[]
[AuxVariables]
[tag_variable]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = KokkosDiffusion
variable = u
extra_matrix_tags = 'mat_tag1 mat_tag2'
[]
[diff1]
type = KokkosDiffusion
variable = u
extra_matrix_tags = 'mat_tag2'
vector_tags = vec_tag1
[]
[diff2]
type = KokkosDiffusion
variable = u
vector_tags = vec_tag1
[]
[diff3]
type = KokkosDiffusion
variable = u
vector_tags = vec_tag1
[]
[]
[AuxKernels]
[TagMatrixAux]
type = TagMatrixAux
variable = tag_variable
v = u
matrix_tag = mat_tag2
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 3
value = 0
extra_matrix_tags = mat_tag1
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = 1
value = 1
extra_matrix_tags = mat_tag1
[]
[]
[Problem]
type = TagTestProblem
test_tag_vectors = 'nontime residual'
test_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_vectors = 'vec_tag1'
[]
[Executioner]
type = Steady
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(test/tests/kokkos/materials/error/kokkos_dependency.i)
# This test checks that the usage of an old/older (stateful) material property
# does not create a dependency on that property for the purposes of
# dependency resolution for material property evaluation.
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
active = 'u'
[u]
order = FIRST
family = LAGRANGE
[]
[]
[Kernels]
[diff]
type = KokkosCoefDiffusion
variable = u
coef = 0.1
[]
[]
[BCs]
[left]
type = KokkosDirichletBC
variable = u
boundary = 'left'
value = 1
[]
[right]
type = KokkosDirichletBC
variable = u
boundary = 'right'
value = 2
[]
[]
[Materials]
[mat1]
type = KokkosCoupledMaterial
mat_prop = 'prop-a'
coupled_mat_prop = 'prop-b'
use_old_prop = true
block = 0
[]
[mat2]
type = KokkosCoupledMaterial
mat_prop = 'prop-b'
coupled_mat_prop = 'prop-a'
use_old_prop = false
block = 0
[]
[]
[Executioner]
type = Steady
solve_type = 'PJFNK'
[]
[Debug]
show_material_props = true
[]
(framework/include/kokkos/bcs/KokkosDirichletBC.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "KokkosDirichletBCBase.h"
template <typename DirichletBC>
class KokkosDirichletBC : public Moose::Kokkos::DirichletBCBase<DirichletBC>
{
usingKokkosDirichletBCBaseMembers(DirichletBC);
public:
static InputParameters validParams();
KokkosDirichletBC(const InputParameters & parameters);
KOKKOS_FUNCTION Real computeValue(const unsigned int /* qp */, AssemblyDatum & /* datum */) const
{
return _value;
}
protected:
const Moose::Kokkos::Scalar<const Real> _value;
};
template <typename DirichletBC>
InputParameters
KokkosDirichletBC<DirichletBC>::validParams()
{
InputParameters params = Moose::Kokkos::DirichletBCBase<DirichletBC>::validParams();
params.addRequiredParam<Real>("value", "Value of the BC");
params.declareControllable("value");
return params;
}
template <typename DirichletBC>
KokkosDirichletBC<DirichletBC>::KokkosDirichletBC(const InputParameters & parameters)
: Moose::Kokkos::DirichletBCBase<DirichletBC>(parameters),
_value(this->template getParam<Real>("value"))
{
}
class KokkosDirichletBCWrapper final : public KokkosDirichletBC<KokkosDirichletBCWrapper>
{
public:
static InputParameters validParams();
KokkosDirichletBCWrapper(const InputParameters & parameters);
};
(test/include/kokkos/bcs/KokkosCoupledDirichletBC.h)
// This file is part of the MOOSE framework
// https://mooseframework.inl.gov
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "KokkosDirichletBC.h"
/**
* Implements the Dirichlet boundary condition
* c*u + u^2 + v^2 = _value
* where "u" is the current variable, and "v" is a coupled variable.
* Note: without the constant term, a zero initial guess gives you a
* zero row in the Jacobian, which is a bad thing.
*/
class KokkosCoupledDirichletBC final : public KokkosDirichletBC<KokkosCoupledDirichletBC>
{
public:
static InputParameters validParams();
KokkosCoupledDirichletBC(const InputParameters & parameters);
KOKKOS_FUNCTION Real computeQpResidual(const unsigned int qp, AssemblyDatum & datum) const;
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int qp, AssemblyDatum & datum) const;
KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int jvar,
const unsigned int qp,
AssemblyDatum & datum) const;
protected:
// The coupled variable
const Moose::Kokkos::VariableValue _v;
/// The id of the coupled variable
const unsigned int _v_num;
// The constant (not user-selectable for now)
const Real _c;
};
KOKKOS_FUNCTION inline Real
KokkosCoupledDirichletBC::computeQpResidual(const unsigned int qp, AssemblyDatum & datum) const
{
return _c * _u(datum, qp) + _u(datum, qp) * _u(datum, qp) + _v(datum, qp) * _v(datum, qp) -
_value;
}
KOKKOS_FUNCTION inline Real
KokkosCoupledDirichletBC::computeQpJacobian(const unsigned int qp, AssemblyDatum & datum) const
{
return _c + 2 * _u(datum, qp);
}
KOKKOS_FUNCTION inline Real
KokkosCoupledDirichletBC::computeQpOffDiagJacobian(const unsigned int jvar,
const unsigned int qp,
AssemblyDatum & datum) const
{
if (jvar == _v_num)
return 2 * _v(datum, qp);
else
return 0;
}