- blockThe list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains.
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains.
- coefficient1.Name of property k to multiply the integrator by. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number.
Default:1.
C++ Type:MFEMScalarCoefficientName
Controllable:No
Description:Name of property k to multiply the integrator by. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number.
- variableVariable labelling the weak form this kernel is added to
C++ Type:VariableName
Unit:(no unit assumed)
Controllable:No
Description:Variable labelling the weak form this kernel is added to
MFEMVectorFEMassKernel
Overview
Adds the domain integrator for integrating the bilinear form
where or , and is a scalar coefficient.
This term arises from the weak form of the mass operator
Example Input File Syntax
[Kernels<<<{"href": "../../../syntax/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": "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": "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": "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
[]
[](test/tests/mfem/kernels/curlcurl.i)Input Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:No
Description:Set the enabled status of the MooseObject.
Advanced Parameters
Input Files
- (test/tests/mfem/vectorpostprocessors/point_value_sampler/point_value_sampler_curlcurl.i)
- (test/tests/mfem/kernels/curlcurl.i)
- (test/tests/mfem/complex/complex_waveguide.i)
- (test/tests/mfem/submeshes/av_magnetostatic.i)
- (test/tests/mfem/kernels/graddiv.i)
- (test/tests/mfem/kernels/darcy.i)
- (test/tests/mfem/vectorpostprocessors/line_value_sampler/line_value_sampler_curlcurl.i)
Child Objects
(test/tests/mfem/kernels/curlcurl.i)
# Definite Maxwell problem solved with Nedelec elements of the first kind
# based on MFEM Example 3.
[Mesh]
type = MFEMMesh
file = ../mesh/small_fichera.mesh
dim = 3
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
inactive = "L2FESpace"
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
[]
[HDivFESpace]
type = MFEMVectorFESpace
fec_type = RT
fec_order = CONSTANT
[]
[L2FESpace]
type = MFEMScalarFESpace
fec_type = L2
fec_order = CONSTANT
[]
[]
[Variables]
[e_field]
type = MFEMVariable
fespace = HCurlFESpace
[]
[]
[AuxVariables]
inactive = "joule_heating"
[db_dt_field]
type = MFEMVariable
fespace = HDivFESpace
[]
[joule_heating]
type = MFEMVariable
fespace = L2FESpace
[]
[]
[AuxKernels]
inactive = "joule_Q_aux"
[curl]
type = MFEMCurlAux
variable = db_dt_field
source = e_field
scale_factor = -1.0
execute_on = TIMESTEP_END
[]
[joule_Q_aux]
type = MFEMInnerProductAux
variable = joule_heating
first_source_vec = e_field
second_source_vec = e_field
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 = MFEMVectorTangentialDirichletBC
variable = e_field
vector_coefficient = exact_e_field
[]
[]
[Kernels]
[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-12
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/CurlCurl
vtk_format = ASCII
[]
[]
(test/tests/mfem/vectorpostprocessors/point_value_sampler/point_value_sampler_curlcurl.i)
# Definite Maxwell problem solved with Nedelec elements of the first kind
# based on MFEM Example 3. Sampled with MFEMPointValueSampler.
[Mesh]
type = MFEMMesh
file = ../../mesh/small_fichera.mesh
[]
[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 = MFEMVectorTangentialDirichletBC
variable = e_field
vector_coefficient = exact_e_field
[]
[]
[Kernels]
[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
[]
[VectorPostprocessors]
[point_sample]
type = MFEMPointValueSampler
variable = 'e_field'
points = '1 1 -0.5 1 1 0.5'
[]
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Outputs]
execute_on = 'timestep_end'
csv = true
[]
(test/tests/mfem/kernels/curlcurl.i)
# Definite Maxwell problem solved with Nedelec elements of the first kind
# based on MFEM Example 3.
[Mesh]
type = MFEMMesh
file = ../mesh/small_fichera.mesh
dim = 3
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
inactive = "L2FESpace"
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
[]
[HDivFESpace]
type = MFEMVectorFESpace
fec_type = RT
fec_order = CONSTANT
[]
[L2FESpace]
type = MFEMScalarFESpace
fec_type = L2
fec_order = CONSTANT
[]
[]
[Variables]
[e_field]
type = MFEMVariable
fespace = HCurlFESpace
[]
[]
[AuxVariables]
inactive = "joule_heating"
[db_dt_field]
type = MFEMVariable
fespace = HDivFESpace
[]
[joule_heating]
type = MFEMVariable
fespace = L2FESpace
[]
[]
[AuxKernels]
inactive = "joule_Q_aux"
[curl]
type = MFEMCurlAux
variable = db_dt_field
source = e_field
scale_factor = -1.0
execute_on = TIMESTEP_END
[]
[joule_Q_aux]
type = MFEMInnerProductAux
variable = joule_heating
first_source_vec = e_field
second_source_vec = e_field
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 = MFEMVectorTangentialDirichletBC
variable = e_field
vector_coefficient = exact_e_field
[]
[]
[Kernels]
[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-12
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/CurlCurl
vtk_format = ASCII
[]
[]
(test/tests/mfem/complex/complex_waveguide.i)
freq = 900e6
angfreq = ${fparse 2*pi*freq}
epsilon0 = 8.8541878176e-12
mu0 = ${fparse 4e-7*pi}
magnetic_reluctivity = ${fparse 1/mu0}
elec_cond_mouse = 0.97
elec_cond_air = 1e-323
[Mesh]
type = MFEMMesh
file = ../mesh/waveguide.g
dim = 3
[]
[Problem]
type = MFEMProblem
numeric_type = complex
[]
[FESpaces]
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
[]
[HDivFESpace]
type = MFEMVectorFESpace
fec_type = RT
fec_order = CONSTANT
[]
[]
[Variables]
[E]
type = MFEMComplexVariable
fespace = HCurlFESpace
[]
[]
[Functions]
[mass_coef_mouse]
type = ParsedFunction
expression = -43*${epsilon0}*${angfreq}^2
[]
[loss_coef_mouse]
type = ParsedFunction
expression = ${angfreq}*${elec_cond_mouse}
[]
[mass_coef_air]
type = ParsedFunction
expression = -${epsilon0}*${angfreq}^2
[]
[loss_coef_air]
type = ParsedFunction
expression = ${angfreq}*${elec_cond_air}
[]
[]
[FunctorMaterials]
[Mouse]
type = MFEMGenericFunctorMaterial
prop_names = 'massCoef lossCoef MagReluctivity'
prop_values = 'mass_coef_mouse loss_coef_mouse ${magnetic_reluctivity}'
block = 1
[]
[Air]
type = MFEMGenericFunctorMaterial
prop_names = 'massCoef lossCoef MagReluctivity'
prop_values = 'mass_coef_air loss_coef_air ${magnetic_reluctivity}'
block = 2
[]
[]
[BCs]
[tangential_E]
type = MFEMComplexVectorTangentialDirichletBC
variable = E
boundary = '2 3 4'
[]
[WaveguidePortIn]
type = MFEMRWTE10IntegratedBC
variable = E
boundary = '5'
input_port = true
port_length_vector = "24.76e-2 0.0 0.0"
port_width_vector = "0.0 12.38e-2 0.0"
frequency = ${freq}
epsilon = ${epsilon0}
mu = ${mu0}
[]
[WaveguidePortOut]
type = MFEMRWTE10IntegratedBC
variable = E
boundary = '6'
input_port = false
port_length_vector = "24.76e-2 0.0 0.0"
port_width_vector = "0.0 12.38e-2 0.0"
frequency = ${freq}
epsilon = ${epsilon0}
mu = ${mu0}
[]
[]
[Kernels]
[curlcurl]
type = MFEMComplexKernel
variable = E
[RealComponent]
type = MFEMCurlCurlKernel
coefficient = MagReluctivity
[]
[]
[mass_loss]
type = MFEMComplexKernel
variable = E
[RealComponent]
type = MFEMVectorFEMassKernel
coefficient = massCoef
[]
[ImagComponent]
type = MFEMVectorFEMassKernel
coefficient = lossCoef
[]
[]
[]
[Solver]
type = MFEMMUMPS
[]
[Executioner]
type = MFEMSteady
assembly_level = legacy
[]
[Postprocessors]
[ObstructionAbsorption]
type = MFEMComplexVectorPeriodAveragedPostprocessor
coefficient = ${elec_cond_mouse}
dual_variable = E
primal_variable = E
block = 1
[]
[]
[Outputs]
[ReportedPostprocessors]
type = CSV
file_base = OutputData/ComplexWaveguide
[]
[]
(test/tests/mfem/submeshes/av_magnetostatic.i)
# Magnetostatic problem solved on a closed conductor subject to
# global loop voltage constraint.
[Mesh]
type = MFEMMesh
file = ../mesh/embedded_concentric_torus.e
[]
[Problem]
type = MFEMProblem
[]
[SubMeshes]
inactive = 'fluxcut'
[fluxcut]
type = MFEMCutTransitionSubMesh
cut_boundary = 'MeasurementPlane'
block = 'TorusCore TorusSheath'
transition_subdomain = transition_dom
transition_subdomain_boundary = transition_bdr
closed_subdomain = coil_dom
[]
[]
[FESpaces]
inactive = 'FluxFESpace'
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
[]
[HDivFESpace]
type = MFEMVectorFESpace
fec_type = RT
fec_order = CONSTANT
[]
[FluxFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
submesh = fluxcut
[]
[]
[Variables]
[a_field]
type = MFEMVariable
fespace = HCurlFESpace
[]
[]
[AuxVariables]
inactive = 'flux_e_field'
[b_field]
type = MFEMVariable
fespace = HDivFESpace
[]
[e_field]
type = MFEMVariable
fespace = HCurlFESpace
[]
[flux_e_field]
type = MFEMVariable
fespace = FluxFESpace
[]
[]
[AuxKernels]
[curl]
type = MFEMCurlAux
variable = b_field
source = a_field
scale_factor = 1.0
execute_on = TIMESTEP_END
[]
[]
[Functions]
[exact_a_field]
type = ParsedVectorFunction
expression_x = '0'
expression_y = '0'
expression_z = '0'
[]
[]
[BCs]
[tangential_a_bdr]
type = MFEMVectorTangentialDirichletBC
variable = a_field
vector_coefficient = exact_a_field
boundary = 'Exterior'
[]
[]
[FunctorMaterials]
inactive = 'ConductorBoundary'
[Vacuum]
type = MFEMGenericFunctorMaterial
prop_names = reluctivity
prop_values = 1.0
[]
[Conductor]
type = MFEMGenericFunctorMaterial
prop_names = conductivity
prop_values = 1.0
block = 'TorusCore TorusSheath'
[]
[ConductorBoundary]
type = MFEMGenericFunctorMaterial
prop_names = conductivity_boundary
prop_values = 1.0
boundary = 'MeasurementPlane'
[]
[]
[Kernels]
[mass]
type = MFEMVectorFEMassKernel
variable = a_field
coefficient = 1e-10
[]
[curlcurl]
type = MFEMCurlCurlKernel
variable = a_field
coefficient = reluctivity
[]
[source]
type = MFEMMixedVectorMassKernel
variable = a_field
trial_variable = e_field
coefficient = conductivity
block = 'TorusCore TorusSheath'
[]
[]
[Preconditioner]
[ams]
type = MFEMHypreAMS
fespace = HCurlFESpace
[]
[]
[Solver]
type = MFEMHyprePCG
preconditioner = ams
l_tol = 1e-14
l_max_its = 1000
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[MultiApps]
[coil]
type = FullSolveMultiApp
input_files = cut_closed_coil.i
execute_on = INITIAL
[]
[]
[Transfers]
inactive = 'submesh_transfer_to_fluxsurface'
[from_coil]
type = MultiAppMFEMCopyTransfer
source_variables = e_field
variables = e_field
from_multi_app = coil
[]
[submesh_transfer_to_fluxsurface]
type = MFEMSubMeshTransfer
from_variable = e_field
to_variable = flux_e_field
execute_on = TIMESTEP_END
[]
[]
[Postprocessors]
inactive = 'CoilCurrent'
[CoilPower]
type = MFEMVectorFEInnerProductIntegralPostprocessor
coefficient = conductivity
dual_variable = e_field
primal_variable = e_field
block = 'TorusCore TorusSheath'
[]
[CoilCurrent]
type = MFEMVectorBoundaryFluxIntegralPostprocessor
coefficient = conductivity_boundary
variable = flux_e_field
boundary = 'MeasurementPlane'
[]
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/MagnetostaticClosedCoil
vtk_format = ASCII
[]
[ReportedPostprocessors]
type = CSV
file_base = OutputData/AVMagnetostaticClosedCoilCSV
[]
[]
(test/tests/mfem/kernels/graddiv.i)
# Grad-div problem using method of manufactured solutions,
# based on MFEM Example 4.
[Mesh]
type = MFEMMesh
file = ../mesh/beam-tet.mesh
dim = 3
uniform_refine = 1
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[HDivFESpace]
type = MFEMVectorFESpace
fec_type = RT
fec_order = CONSTANT
ordering = "vdim"
[]
[L2FESpace]
type = MFEMScalarFESpace
fec_type = L2
fec_order = CONSTANT
basis = GaussLegendre
[]
[]
[Variables]
[F]
type = MFEMVariable
fespace = HDivFESpace
[]
[]
[AuxVariables]
[divF]
type = MFEMVariable
fespace = L2FESpace
[]
[]
[AuxKernels]
[div]
type = MFEMDivAux
variable = divF
source = F
execute_on = TIMESTEP_END
[]
[]
[Functions]
[f]
type = ParsedVectorFunction
expression_x = '(1. + 2*kappa * kappa) * cos(kappa * x) * sin(kappa * y)'
expression_y = '(1. + 2*kappa * kappa) * cos(kappa * y) * sin(kappa * x)'
expression_z = '0'
symbol_names = kappa
symbol_values = 3.1415926535
[]
[F_exact]
type = ParsedVectorFunction
expression_x = 'cos(kappa * x) * sin(kappa * y)'
expression_y = 'cos(kappa * y) * sin(kappa * x)'
expression_z = '0'
symbol_names = kappa
symbol_values = 3.1415926535
[]
[]
[BCs]
[dirichlet]
type = MFEMVectorNormalDirichletBC
variable = F
boundary = '1 2 3'
vector_coefficient = F_exact
[]
[]
[Kernels]
[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]
type = MFEMSteady
device = "cpu"
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/GradDiv
vtk_format = ASCII
[]
[]
(test/tests/mfem/kernels/darcy.i)
[Mesh]
type = MFEMMesh
file = ../mesh/star.mesh
uniform_refine = 2
[]
[Problem]
type = MFEMProblem
[]
[Functions]
[exact_velocity]
type = ParsedVectorFunction
expression_x = '-exp(x) * sin(y)'
expression_y = '-exp(x) * cos(y)'
[]
[exact_pressure]
type = ParsedFunction
expression = 'exp(x) * sin(y)'
[]
[exact_pressure_rhs]
type = ParsedFunction
expression = '-exp(x) * sin(y)'
[]
[]
[FESpaces]
[HDivFESpace]
type = MFEMVectorFESpace
fec_type = RT
fec_order = SECOND
[]
[L2FESpace]
type = MFEMScalarFESpace
fec_type = L2
fec_order = SECOND
basis = GaussLegendre
[]
[]
[Variables]
[velocity]
type = MFEMVariable
fespace = HDivFESpace
[]
[pressure]
type = MFEMVariable
fespace = L2FESpace
[]
[]
[BCs]
[flux_boundaries]
type = MFEMVectorFEBoundaryFluxIntegratedBC
variable = velocity
coefficient = exact_pressure_rhs
[]
[]
[Kernels]
[VelocityMass]
type = MFEMVectorFEMassKernel
variable = velocity
[]
[PressureGrad]
type = MFEMVectorFEDivergenceKernel
trial_variable = pressure
variable = velocity
coefficient = -1
transpose = true
[]
[VelocityDiv]
type = MFEMVectorFEDivergenceKernel
trial_variable = velocity
variable = pressure
coefficient = -1
[]
[]
[Solver]
type = MFEMMUMPS
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Postprocessors]
[velocity_error]
type = MFEMVectorL2Error
variable = velocity
function = exact_velocity
[]
[pressure_error]
type = MFEML2Error
variable = pressure
function = exact_pressure
[]
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/Darcy
vtk_format = ASCII
[]
[DarcyErrorCSV]
type = CSV
file_base = OutputData/Darcy
[]
[]
(test/tests/mfem/vectorpostprocessors/line_value_sampler/line_value_sampler_curlcurl.i)
# Definite Maxwell problem solved with Nedelec elements of the first kind
# based on MFEM Example 3. Sampled with MFEMLineValueSampler.
[Mesh]
type = MFEMMesh
file = ../../mesh/small_fichera.mesh
[]
[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 = MFEMVectorTangentialDirichletBC
variable = e_field
vector_coefficient = exact_e_field
[]
[]
[Kernels]
[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
[]
[VectorPostprocessors]
[line_sample]
type = MFEMLineValueSampler
variable = 'e_field'
start_point = '1 1 -1'
end_point = '1 1 1'
num_points = 11
[]
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Outputs]
execute_on = 'timestep_end'
csv = true
[]
(framework/include/mfem/kernels/MFEMTimeDerivativeVectorFEMassKernel.h)
// This file is part of the MOOSE framework
// https://mooseframework.inl.gov
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#ifdef MOOSE_MFEM_ENABLED
#pragma once
#include "MFEMVectorFEMassKernel.h"
/**
* \f[
* (k du/dt, v)
* \f]
*/
class MFEMTimeDerivativeVectorFEMassKernel : public MFEMVectorFEMassKernel
{
public:
static InputParameters validParams();
MFEMTimeDerivativeVectorFEMassKernel(const InputParameters & parameters);
/// Get name of the trial variable (gridfunction) the kernel acts on.
/// Defaults to the name of the test variable labelling the weak form.
virtual const VariableName & getTrialVariableName() const override { return _var_dot_name; }
protected:
/// Name of variable (gridfunction) representing time derivative of variable.
const VariableName _var_dot_name;
};
#endif