Definite Maxwell Problem

Summary

Solves a 3D electromagnetic diffusion problem for the electric field on a cube missing an octant, discretized using conforming Nédélec elements. This example is based on MFEM Example 3.

Description

This problem solves a definite Maxwell equation with strong form:

where

In this example, we solve this using the weak form

where

Example File

# Definite Maxwell problem solved with Nedelec elements of the first kind
# based on MFEM Example 3.

[Mesh<<<{"href": "../Mesh/index.html"}>>>]
  type = MFEMMesh
  file = ../mesh/small_fichera.mesh
  dim = 3
[]

[Problem<<<{"href": "../Problem/index.html"}>>>]
  type = MFEMProblem
[]

[FESpaces]
  [HCurlFESpace]
    type = MFEMVectorFESpace
    fec_type = ND
    fec_order = FIRST
  []
  [HDivFESpace]
    type = MFEMVectorFESpace
    fec_type = RT
    fec_order = CONSTANT
  []
[]

[Variables<<<{"href": "../Variables/index.html"}>>>]
  [e_field]
    type = MFEMVariable
    fespace = HCurlFESpace
  []
[]

[AuxVariables<<<{"href": "../AuxVariables/index.html"}>>>]
  [db_dt_field]
    type = MFEMVariable
    fespace = HDivFESpace
  []
[]

[AuxKernels<<<{"href": "../AuxKernels/index.html"}>>>]
  [curl]
    type = MFEMCurlAux
    variable = db_dt_field
    source = e_field
    scale_factor = -1.0
    execute_on = TIMESTEP_END
  []
[]

[Functions<<<{"href": "../Functions/index.html"}>>>]
  [exact_e_field]
    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."}>>> = 'sin(kappa * y)'
    expression_y<<<{"description": "y-component of function."}>>> = 'sin(kappa * z)'
    expression_z<<<{"description": "z-component of function."}>>> = 'sin(kappa * x)'

    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
  []

  [forcing_field]
    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. + kappa * kappa) * sin(kappa * y)'
    expression_y<<<{"description": "y-component of function."}>>> = '(1. + kappa * kappa) * sin(kappa * z)'
    expression_z<<<{"description": "z-component of function."}>>> = '(1. + kappa * kappa) * sin(kappa * x)'

    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"}>>>]
  [tangential_E_bdr]
    type = MFEMVectorTangentialDirichletBC
    variable = e_field
    vector_coefficient = exact_e_field
  []
[]

[Kernels<<<{"href": "../Kernels/index.html"}>>>]
  [curlcurl]
    type = MFEMCurlCurlKernel
    variable = e_field
  []
  [mass]
    type = MFEMVectorFEMassKernel
    variable = e_field
  []
  [source]
    type = MFEMVectorFEDomainLFKernel
    variable = e_field
    vector_coefficient = forcing_field
  []
[]

[Preconditioner]
  [ams]
    type = MFEMHypreAMS
    fespace = HCurlFESpace
  []
[]

[Solver]
  type = MFEMHypreGMRES
  preconditioner = ams
  l_tol = 1e-6
[]

[Executioner<<<{"href": "../Executioner/index.html"}>>>]
  type = MFEMSteady
  device = cpu
[]

[Outputs<<<{"href": "../Outputs/index.html"}>>>]
  [ParaViewDataCollection]
    type = MFEMParaViewDataCollection
    file_base = OutputData/CurlCurl
    vtk_format = ASCII
  []
[]
(moose/test/tests/mfem/kernels/curlcurl.i)