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<<<{"href": "../FESpaces/index.html"}>>>]
[HCurlFESpace]
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."}>>> = ND
fec_order<<<{"description": "Order of the FE shape function to use."}>>> = FIRST
[]
[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
[]
[]
[Variables<<<{"href": "../Variables/index.html"}>>>]
[e_field]
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."}>>> = HCurlFESpace
[]
[]
[AuxVariables<<<{"href": "../AuxVariables/index.html"}>>>]
[db_dt_field]
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
[]
[]
[AuxKernels<<<{"href": "../AuxKernels/index.html"}>>>]
[curl]
type = MFEMCurlAux<<<{"description": "Calculates the curl of an H(curl) conforming ND source variable and stores the result on an H(div) conforming RT result auxvariable", "href": "../../source/mfem/auxkernels/MFEMCurlAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = db_dt_field
source<<<{"description": "Vector H(curl) MFEMVariable to take the curl of."}>>> = e_field
scale_factor<<<{"description": "Factor to scale result auxvariable by."}>>> = -1.0
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"}>>>]
[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<<<{"description": "Applies a Dirichlet condition to the tangential components of a vector variable.", "href": "../../source/mfem/bcs/MFEMVectorTangentialDirichletBC.html"}>>>
variable<<<{"description": "Variable on which to apply the boundary condition"}>>> = e_field
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)."}>>> = exact_e_field
[]
[]
[Kernels<<<{"href": "../Kernels/index.html"}>>>]
[curlcurl]
type = MFEMCurlCurlKernel<<<{"description": "Adds the domain integrator to an MFEM problem for the bilinear form $(k\\vec\\nabla \\times \\vec u, \\vec\\nabla \\times \\vec v)_\\Omega$ arising from the weak form of the curl curl operator $k\\vec\\nabla \\times \\vec\\nabla \\times \\vec u$.", "href": "../../source/mfem/kernels/MFEMCurlCurlKernel.html"}>>>
variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = e_field
[]
[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"}>>> = e_field
[]
[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"}>>> = e_field
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)."}>>> = forcing_field
[]
[]
[Preconditioner<<<{"href": "../Preconditioner/index.html"}>>>]
[ams]
type = MFEMHypreAMS<<<{"description": "Hypre auxiliary-space Maxwell solver and preconditioner for the iterative solution of MFEM equation systems.", "href": "../../source/mfem/solvers/MFEMHypreAMS.html"}>>>
fespace<<<{"description": "H(curl) FESpace to use in HypreAMS setup."}>>> = HCurlFESpace
[]
[]
[Solver<<<{"href": "../Solver/index.html"}>>>]
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<<<{"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/CurlCurl
vtk_format<<<{"description": "Select VTK data format to use, choosing between BINARY, BINARY32, and ASCII."}>>> = ASCII
[]
[]
(test/tests/mfem/kernels/curlcurl.i)