ScalarLMKernel
This Kernel
demonstrates the usage of the scalar augmentation class described in KernelScalarBase. This single kernel is an alternative to the combination of ScalarLagrangeMultiplier, AverageValueConstraint, and an Elemental Integral Postprocessor. All terms from the spatial and scalar variables are handled by this object.
This Kernel enforces the constraint of
where is a given constant, using a Lagrange multiplier approach. Since this is a single constraint, a single scalar variable is required, as shown below.
[Variables]
[lambda]
family = SCALAR
order = FIRST
[]
[]
The residual of related to the Lagrange multiplier is given as:
(1)
The residual of the Lagrange multiplier is given as:
(2)
This Kernel implements the residual and Jacobian contributions to the spatial variable equation from the in Eq. (1). Also, it implements the residual and Jacobian contributions to the scalar Lagrange multiplier constraint equation from the spatial variable in Eq. (2).
So that the input syntax matches with ScalarLagrangeMultiplier, a Postprocessor
is required that computes the total volume of the domain, assigned to the parameter pp_name
.
Currently, a NullScalarKernel is required to activate the dependency of the scalar variable within the block or subdomain of this object. See one of the example files listed below.
The detailed description of the derivation of the weak form can be found at scalar_constraint_kernel.
(test/tests/kernels/scalar_kernel_constraint/diffusion_override_scalar.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Functions]
[exact_fn]
type = ParsedFunction
value = 'x*x+y*y'
[]
[ffn]
type = ParsedFunction
value = -4
[]
[bottom_bc_fn]
type = ParsedFunction
value = -2*y
[]
[right_bc_fn]
type = ParsedFunction
value = 2*x
[]
[top_bc_fn]
type = ParsedFunction
value = 2*y
[]
[left_bc_fn]
type = ParsedFunction
value = -2*x
[]
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
# Make sure that we can derive from the scalar base class
# but actually not assign a scalar variable
[diff]
type = DiffusionNoScalar
variable = u
scalar_variable = lambda
[]
[ffnk]
type = BodyForce
variable = u
function = ffn
[]
[sk_lm]
type = ScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = FunctionNeumannBC
variable = u
boundary = 'bottom'
function = bottom_bc_fn
[]
[right]
type = FunctionNeumannBC
variable = u
boundary = 'right'
function = right_bc_fn
[]
[top]
type = FunctionNeumannBC
variable = u
boundary = 'top'
function = top_bc_fn
[]
[left]
type = FunctionNeumannBC
variable = u
boundary = 'left'
function = left_bc_fn
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-9
l_tol = 1.e-10
nl_max_its = 10
# This example builds an indefinite matrix, so "-pc_type hypre -pc_hypre_type boomeramg" cannot
# be used reliably on this problem. ILU(0) seems to do OK in both serial and parallel in my testing,
# I have not seen any zero pivot issues.
petsc_options_iname = '-pc_type -sub_pc_type'
petsc_options_value = 'bjacobi ilu'
# This is a linear problem, so we don't need to recompute the
# Jacobian. This isn't a big deal for a Steady problems, however, as
# there is only one solve.
solve_type = 'LINEAR'
[]
[Outputs]
# exodus = true
csv = true
hide = lambda
[]
(test/tests/kernels/scalar_kernel_constraint/scalar_constraint_kernel.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Functions]
[exact_fn]
type = ParsedFunction
value = 'x*x+y*y'
[]
[ffn]
type = ParsedFunction
value = -4
[]
[bottom_bc_fn]
type = ParsedFunction
value = -2*y
[]
[right_bc_fn]
type = ParsedFunction
value = 2*x
[]
[top_bc_fn]
type = ParsedFunction
value = 2*y
[]
[left_bc_fn]
type = ParsedFunction
value = -2*x
[]
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
[diff]
type = Diffusion
variable = u
[]
[ffnk]
type = BodyForce
variable = u
function = ffn
[]
[sk_lm]
type = ScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = FunctionNeumannBC
variable = u
boundary = 'bottom'
function = bottom_bc_fn
[]
[right]
type = FunctionNeumannBC
variable = u
boundary = 'right'
function = right_bc_fn
[]
[top]
type = FunctionNeumannBC
variable = u
boundary = 'top'
function = top_bc_fn
[]
[left]
type = FunctionNeumannBC
variable = u
boundary = 'left'
function = left_bc_fn
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-9
l_tol = 1.e-10
nl_max_its = 10
# This example builds an indefinite matrix, so "-pc_type hypre -pc_hypre_type boomeramg" cannot
# be used reliably on this problem
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
# This is a linear problem, so we don't need to recompute the
# Jacobian. This isn't a big deal for a Steady problems, however, as
# there is only one solve.
solve_type = 'LINEAR'
[]
[Outputs]
exodus = true
hide = lambda
[]
(test/tests/kernels/ad_scalar_kernel_constraint/diffusion_bipass_scalar.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Functions]
[exact_fn]
type = ParsedFunction
value = 'x*x+y*y'
[]
[ffn]
type = ParsedFunction
value = -4
[]
[bottom_bc_fn]
type = ParsedFunction
value = -2*y
[]
[right_bc_fn]
type = ParsedFunction
value = 2*x
[]
[top_bc_fn]
type = ParsedFunction
value = 2*y
[]
[left_bc_fn]
type = ParsedFunction
value = -2*x
[]
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
# Make sure that we can derive from the scalar base class
# but actually not assign a scalar variable
[diff]
type = ADDiffusionNoScalar
variable = u
[]
[ffnk]
type = ADBodyForce
variable = u
function = ffn
[]
[sk_lm]
type = ADScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = ADFunctionNeumannBC
variable = u
boundary = 'bottom'
function = bottom_bc_fn
[]
[right]
type = ADFunctionNeumannBC
variable = u
boundary = 'right'
function = right_bc_fn
[]
[top]
type = ADFunctionNeumannBC
variable = u
boundary = 'top'
function = top_bc_fn
[]
[left]
type = ADFunctionNeumannBC
variable = u
boundary = 'left'
function = left_bc_fn
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
residual_and_jacobian_together = true
nl_rel_tol = 1e-9
l_tol = 1.e-10
nl_max_its = 10
# This example builds an indefinite matrix, so "-pc_type hypre -pc_hypre_type boomeramg" cannot
# be used reliably on this problem
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
# This is a linear problem, so we don't need to recompute the
# Jacobian. This isn't a big deal for a Steady problems, however, as
# there is only one solve.
solve_type = 'LINEAR'
[]
[Outputs]
# exodus = true
csv = true
hide = lambda
[]
(test/tests/kernels/ad_scalar_kernel_constraint/scalar_constraint_together.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Functions]
[exact_fn]
type = ParsedFunction
value = 'x*x+y*y'
[]
[ffn]
type = ParsedFunction
value = -4
[]
[bottom_bc_fn]
type = ParsedFunction
value = -2*y
[]
[right_bc_fn]
type = ParsedFunction
value = 2*x
[]
[top_bc_fn]
type = ParsedFunction
value = 2*y
[]
[left_bc_fn]
type = ParsedFunction
value = -2*x
[]
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
[diff]
type = ADDiffusion
variable = u
[]
[ffnk]
type = ADBodyForce
variable = u
function = ffn
[]
[sk_lm]
type = ADScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = ADFunctionNeumannBC
variable = u
boundary = 'bottom'
function = bottom_bc_fn
[]
[right]
type = ADFunctionNeumannBC
variable = u
boundary = 'right'
function = right_bc_fn
[]
[top]
type = ADFunctionNeumannBC
variable = u
boundary = 'top'
function = top_bc_fn
[]
[left]
type = ADFunctionNeumannBC
variable = u
boundary = 'left'
function = left_bc_fn
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Executioner]
type = Steady
residual_and_jacobian_together = true
nl_rel_tol = 1e-9
l_tol = 1.e-10
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist'
solve_type = NEWTON
[]
[Outputs]
# exodus = true
csv = true
hide = lambda
[]
(test/tests/kernels/ad_scalar_kernel_constraint/scalar_constraint_kernel_RJ.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
[diff]
type = ADDiffusion
variable = u
[]
[sk_lm]
type = ADScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = ADDirichletBC
variable = u
boundary = 'bottom'
value = 0
[]
[right]
type = ADDirichletBC
variable = u
boundary = 'right'
value = 0
[]
[top]
type = ADDirichletBC
variable = u
boundary = 'top'
value = 0
[]
[left]
type = ADDirichletBC
variable = u
boundary = 'left'
value = 0
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[]
# Force LU decomposition, nonlinear iterations, to check Jacobian terms with single factorization
[Executioner]
type = Steady
residual_and_jacobian_together = true
nl_rel_tol = 1e-9
l_tol = 1.e-10
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist'
solve_type = NEWTON
[]
[Outputs]
exodus = true
hide = lambda
[]
(test/tests/kernels/scalar_kernel_constraint/diffusion_bipass_scalar.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Functions]
[exact_fn]
type = ParsedFunction
value = 'x*x+y*y'
[]
[ffn]
type = ParsedFunction
value = -4
[]
[bottom_bc_fn]
type = ParsedFunction
value = -2*y
[]
[right_bc_fn]
type = ParsedFunction
value = 2*x
[]
[top_bc_fn]
type = ParsedFunction
value = 2*y
[]
[left_bc_fn]
type = ParsedFunction
value = -2*x
[]
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
# Make sure that we can derive from the scalar base class
# but actually not assign a scalar variable
[diff]
type = DiffusionNoScalar
variable = u
[]
[ffnk]
type = BodyForce
variable = u
function = ffn
[]
[sk_lm]
type = ScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = FunctionNeumannBC
variable = u
boundary = 'bottom'
function = bottom_bc_fn
[]
[right]
type = FunctionNeumannBC
variable = u
boundary = 'right'
function = right_bc_fn
[]
[top]
type = FunctionNeumannBC
variable = u
boundary = 'top'
function = top_bc_fn
[]
[left]
type = FunctionNeumannBC
variable = u
boundary = 'left'
function = left_bc_fn
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-9
l_tol = 1.e-10
nl_max_its = 10
# This example builds an indefinite matrix, so "-pc_type hypre -pc_hypre_type boomeramg" cannot
# be used reliably on this problem
petsc_options_iname = '-pc_type -pc_factor_shift_type'
petsc_options_value = 'lu NONZERO'
# This is a linear problem, so we don't need to recompute the
# Jacobian. This isn't a big deal for a Steady problems, however, as
# there is only one solve.
solve_type = 'LINEAR'
[]
[Outputs]
# exodus = true
csv = true
hide = lambda
[]
(test/tests/kernels/ad_scalar_kernel_constraint/diffusion_override_scalar.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 2
ny = 2
elem_type = QUAD9
[]
[Functions]
[exact_fn]
type = ParsedFunction
value = 'x*x+y*y'
[]
[ffn]
type = ParsedFunction
value = -4
[]
[bottom_bc_fn]
type = ParsedFunction
value = -2*y
[]
[right_bc_fn]
type = ParsedFunction
value = 2*x
[]
[top_bc_fn]
type = ParsedFunction
value = 2*y
[]
[left_bc_fn]
type = ParsedFunction
value = -2*x
[]
[]
[Variables]
[u]
family = LAGRANGE
order = SECOND
[]
[lambda]
family = SCALAR
order = FIRST
[]
[]
[Kernels]
# Make sure that we can derive from the scalar base class
# but actually not assign a scalar variable
[diff]
type = ADDiffusionNoScalar
variable = u
scalar_variable = lambda
[]
[ffnk]
type = ADBodyForce
variable = u
function = ffn
[]
[sk_lm]
type = ADScalarLMKernel
variable = u
kappa = lambda
pp_name = pp
value = 2.666666666666666
[]
[]
[Problem]
kernel_coverage_check = false
error_on_jacobian_nonzero_reallocation = true
[]
[BCs]
[bottom]
type = FunctionNeumannBC
variable = u
boundary = 'bottom'
function = bottom_bc_fn
[]
[right]
type = FunctionNeumannBC
variable = u
boundary = 'right'
function = right_bc_fn
[]
[top]
type = FunctionNeumannBC
variable = u
boundary = 'top'
function = top_bc_fn
[]
[left]
type = FunctionNeumannBC
variable = u
boundary = 'left'
function = left_bc_fn
[]
[]
[Postprocessors]
# integrate the volume of domain since original objects set
# int(phi)=V0, rather than int(phi-V0)=0
[pp]
type = FunctionElementIntegral
function = 1
execute_on = initial
[]
[l2_err]
type = ElementL2Error
variable = u
function = exact_fn
execute_on = 'initial timestep_end'
[]
[]
[Preconditioning]
[pc]
type = SMP
full = true
solve_type = 'NEWTON'
[]
[]
[Executioner]
type = Steady
nl_rel_tol = 1e-9
l_tol = 1.e-10
nl_max_its = 10
# This example builds an indefinite matrix, so "-pc_type hypre -pc_hypre_type boomeramg" cannot
# be used reliably on this problem. ILU(0) seems to do OK in both serial and parallel in my testing,
# I have not seen any zero pivot issues.
petsc_options_iname = '-pc_type -sub_pc_type'
petsc_options_value = 'bjacobi ilu'
# This is a linear problem, so we don't need to recompute the
# Jacobian. This isn't a big deal for a Steady problems, however, as
# there is only one solve.
solve_type = 'LINEAR'
[]
[Outputs]
# exodus = true
csv = true
hide = lambda
[]