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<<<{"href": "../FESpaces/index.html"}>>>]
  [HDivFESpace]
    type = MFEMVectorFESpace<<<{"description": "Convenience class to construct vector finite element spaces, abstracting away some of the mathematical complexity of specifying the dimensions.", "href": "../../source/mfem/fespaces/MFEMVectorFESpace.html"}>>>
    fec_type<<<{"description": "Specifies the family of FE shape functions."}>>> = RT
    fec_order<<<{"description": "Order of the FE shape function to use."}>>> = CONSTANT
    ordering<<<{"description": "Ordering style to use for vector DoFs."}>>> = "vdim"
  []
  [L2FESpace]
    type = MFEMScalarFESpace<<<{"description": "Convenience class to construct scalar finite element spaces.", "href": "../../source/mfem/fespaces/MFEMScalarFESpace.html"}>>>
    fec_type<<<{"description": "Specifies the family of FE shape functions."}>>> = L2
    fec_order<<<{"description": "Order of the FE shape function to use."}>>> = CONSTANT
  []
[]

[Variables<<<{"href": "../Variables/index.html"}>>>]
  [F]
    type = MFEMVariable<<<{"description": "Class for adding MFEM variables to the problem (`mfem::ParGridFunction`s).", "href": "../../source/mfem/variables/MFEMVariable.html"}>>>
    fespace<<<{"description": "The finite element space this variable is defined on."}>>> = HDivFESpace
  []
[]

[AuxVariables<<<{"href": "../AuxVariables/index.html"}>>>]
  [divF]
    type = MFEMVariable<<<{"description": "Class for adding MFEM variables to the problem (`mfem::ParGridFunction`s).", "href": "../../source/mfem/variables/MFEMVariable.html"}>>>
    fespace<<<{"description": "The finite element space this variable is defined on."}>>> = L2FESpace
  []
[]

[AuxKernels<<<{"href": "../AuxKernels/index.html"}>>>]
  [div]
    type = MFEMDivAux<<<{"description": "Calculates the divergence of an H(div) conforming RT source variable and stores the result on an L2 conforming result auxvariable", "href": "../../source/mfem/auxkernels/MFEMDivAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = divF
    source<<<{"description": "Vector H(div) MFEMVariable to take the divergence of."}>>> = F
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 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<<<{"description": "Applies a Dirichlet condition to the normal components of a vector variable.", "href": "../../source/mfem/bcs/MFEMVectorNormalDirichletBC.html"}>>>
    variable<<<{"description": "Variable on which to apply the boundary condition"}>>> = F
    boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies. Defaults to all boundaries."}>>> = '1 2 3'
    vector_coefficient<<<{"description": "Vector coefficient specifying the values the variable takes on the boundary. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a numeric vector value (enclosed in curly braces)."}>>> = F_exact
  []
[]

[Kernels<<<{"href": "../Kernels/index.html"}>>>]
  [divdiv]
    type = MFEMDivDivKernel<<<{"description": "Adds the domain integrator to an MFEM problem for the bilinear form $(k\\vec\\nabla \\cdot \\vec u, \\vec\\nabla \\cdot \\vec v)_\\Omega$ arising from the weak form of the grad-div operator $-\\vec\\nabla \\left( k \\vec\\nabla \\cdot \\vec u \\right)$.", "href": "../../source/mfem/kernels/MFEMDivDivKernel.html"}>>>
    variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = F
  []
  [mass]
    type = MFEMVectorFEMassKernel<<<{"description": "Adds the domain integrator to an MFEM problem for the bilinear form $(k \\vec u, \\vec v)_\\Omega$ arising from the weak form of the mass operator $k \\vec u$.", "href": "../../source/mfem/kernels/MFEMVectorFEMassKernel.html"}>>>
    variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = F
  []
  [source]
    type = MFEMVectorFEDomainLFKernel<<<{"description": "Adds the domain integrator to an MFEM problem for the linear form $(\\vec f, \\vec v)_\\Omega$ arising from the weak form of the forcing term $\\vec f$.", "href": "../../source/mfem/kernels/MFEMVectorFEDomainLFKernel.html"}>>>
    variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = F
    vector_coefficient<<<{"description": "The name of the vector coefficient f. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a numeric vector value (enclosed in curly braces)."}>>> = f
  []
[]

[Preconditioner<<<{"href": "../Preconditioner/index.html"}>>>]
  [ADS]
    type = MFEMHypreADS<<<{"description": "Hypre auxiliary-space divergence solver and preconditioner for the iterative solution of MFEM equation systems.", "href": "../../source/mfem/solvers/MFEMHypreADS.html"}>>>
    fespace<<<{"description": "H(div) FESpace to use in HypreADS setup."}>>> = HDivFESpace
  []
[]

[Solver<<<{"href": "../Solver/index.html"}>>>]
  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<<<{"description": "Output for controlling export to an mfem::ParaViewDataCollection.", "href": "../../source/mfem/outputs/MFEMParaViewDataCollection.html"}>>>
    file_base<<<{"description": "The desired solution output name without an extension. If not provided, MOOSE sets it with Outputs/file_base when available. Otherwise, MOOSE uses input file name and this object name for a master input or uses master file_base, the subapp name and this object name for a subapp input to set it."}>>> = OutputData/GradDiv
    vtk_format<<<{"description": "Select VTK data format to use, choosing between BINARY, BINARY32, and ASCII."}>>> = ASCII
  []
[]
(test/tests/mfem/kernels/graddiv.i)