Grad-div Problem
Summary
Solves a diffusion problem for a vector field on a cuboid domain, discretized using conforming Raviart-Thomas elements. This example is based on MFEM Example 4 and is relevant for solving Maxwell's equations using potentials without the Coulomb gauge.
Description
This problem solves a grad-div equation with strong form:
where
In this example, we solve this using the weak form
where
Example File
# Grad-div problem using method of manufactured solutions,
# based on MFEM Example 4.
[Mesh<<<{"href": "../Mesh/index.html"}>>>]
type = MFEMMesh
file = ../mesh/beam-tet.mesh
dim = 3
uniform_refine<<<{"description": "Specify the level of uniform refinement applied to the initial mesh"}>>> = 1
[]
[Problem<<<{"href": "../Problem/index.html"}>>>]
type = MFEMProblem
[]
[FESpaces]
[HDivFESpace]
type = MFEMVectorFESpace
fec_type = RT
fec_order = CONSTANT
ordering = "vdim"
[]
[L2FESpace]
type = MFEMScalarFESpace
fec_type = L2
fec_order = CONSTANT
[]
[]
[Variables<<<{"href": "../Variables/index.html"}>>>]
[F]
type = MFEMVariable
fespace = HDivFESpace
[]
[]
[AuxVariables<<<{"href": "../AuxVariables/index.html"}>>>]
[divF]
type = MFEMVariable
fespace = L2FESpace
[]
[]
[AuxKernels<<<{"href": "../AuxKernels/index.html"}>>>]
[div]
type = MFEMDivAux
variable = divF
source = F
execute_on = TIMESTEP_END
[]
[]
[Functions<<<{"href": "../Functions/index.html"}>>>]
[f]
type = ParsedVectorFunction<<<{"description": "Returns a vector function based on string descriptions for each component.", "href": "../../source/functions/MooseParsedVectorFunction.html"}>>>
expression_x<<<{"description": "x-component of function."}>>> = '(1. + 2*kappa * kappa) * cos(kappa * x) * sin(kappa * y)'
expression_y<<<{"description": "y-component of function."}>>> = '(1. + 2*kappa * kappa) * cos(kappa * y) * sin(kappa * x)'
expression_z<<<{"description": "z-component of function."}>>> = '0'
symbol_names<<<{"description": "Symbols (excluding t,x,y,z) that are bound to the values provided by the corresponding items in the vals vector."}>>> = kappa
symbol_values<<<{"description": "Constant numeric values, postprocessor names, function names, and scalar variables corresponding to the symbols in symbol_names."}>>> = 3.1415926535
[]
[F_exact]
type = ParsedVectorFunction<<<{"description": "Returns a vector function based on string descriptions for each component.", "href": "../../source/functions/MooseParsedVectorFunction.html"}>>>
expression_x<<<{"description": "x-component of function."}>>> = 'cos(kappa * x) * sin(kappa * y)'
expression_y<<<{"description": "y-component of function."}>>> = 'cos(kappa * y) * sin(kappa * x)'
expression_z<<<{"description": "z-component of function."}>>> = '0'
symbol_names<<<{"description": "Symbols (excluding t,x,y,z) that are bound to the values provided by the corresponding items in the vals vector."}>>> = kappa
symbol_values<<<{"description": "Constant numeric values, postprocessor names, function names, and scalar variables corresponding to the symbols in symbol_names."}>>> = 3.1415926535
[]
[]
[BCs<<<{"href": "../BCs/index.html"}>>>]
[dirichlet]
type = MFEMVectorNormalDirichletBC
variable = F
boundary = '1 2 3'
vector_coefficient = F_exact
[]
[]
[Kernels<<<{"href": "../Kernels/index.html"}>>>]
[divdiv]
type = MFEMDivDivKernel
variable = F
[]
[mass]
type = MFEMVectorFEMassKernel
variable = F
[]
[source]
type = MFEMVectorFEDomainLFKernel
variable = F
vector_coefficient = f
[]
[]
[Preconditioner]
[ADS]
type = MFEMHypreADS
fespace = HDivFESpace
[]
[]
[Solver]
type = MFEMCGSolver
preconditioner = ADS
l_tol = 1e-16
l_max_its = 1000
print_level = 2
[]
[Executioner<<<{"href": "../Executioner/index.html"}>>>]
type = MFEMSteady
device = "cpu"
[]
[Outputs<<<{"href": "../Outputs/index.html"}>>>]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/GradDiv
vtk_format = ASCII
[]
[]
(moose/test/tests/mfem/kernels/graddiv.i)