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]
type = MFEMMesh
file = gold/small_fichera.mesh
dim = 3
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
[]
[HDivFESpace]
type = MFEMVectorFESpace
fec_type = RT
fec_order = CONSTANT
[]
[]
[Variables]
[e_field]
type = MFEMVariable
fespace = HCurlFESpace
[]
[]
[AuxVariables]
[db_dt_field]
type = MFEMVariable
fespace = HDivFESpace
[]
[]
[AuxKernels]
[curl]
type = MFEMCurlAux
variable = db_dt_field
source = e_field
scale_factor = -1.0
execute_on = TIMESTEP_END
[]
[]
[Functions]
[exact_e_field]
type = ParsedVectorFunction
expression_x = 'sin(kappa * y)'
expression_y = 'sin(kappa * z)'
expression_z = 'sin(kappa * x)'
symbol_names = kappa
symbol_values = 3.1415926535
[]
[forcing_field]
type = ParsedVectorFunction
expression_x = '(1. + kappa * kappa) * sin(kappa * y)'
expression_y = '(1. + kappa * kappa) * sin(kappa * z)'
expression_z = '(1. + kappa * kappa) * sin(kappa * x)'
symbol_names = kappa
symbol_values = 3.1415926535
[]
[]
[BCs]
[tangential_E_bdr]
type = MFEMVectorFunctionTangentialDirichletBC
variable = e_field
function = exact_e_field
[]
[]
[Materials]
[Substance]
type = MFEMGenericConstantMaterial
prop_names = one
prop_values = 1.0
[]
[]
[Kernels]
[curlcurl]
type = MFEMCurlCurlKernel
variable = e_field
coefficient = one
[]
[mass]
type = MFEMVectorFEMassKernel
variable = e_field
coefficient = one
[]
[source]
type = MFEMVectorFEDomainLFKernel
variable = e_field
function = forcing_field
[]
[]
[Preconditioner]
[ams]
type = MFEMHypreAMS
fespace = HCurlFESpace
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = ams
l_tol = 1e-6
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/CurlCurl
vtk_format = ASCII
[]
[]
(test/tests/kernels/curlcurl.i)