- fespaceH(curl) FESpace to use in HypreAMS setup.
C++ Type:MFEMFESpaceName
Controllable:No
Description:H(curl) FESpace to use in HypreAMS setup.
- low_order_refinedFalseSet usage of Low-Order Refined solver.
Default:False
C++ Type:bool
Controllable:No
Description:Set usage of Low-Order Refined solver.
- print_level2Set the solver verbosity.
Default:2
C++ Type:int
Controllable:No
Description:Set the solver verbosity.
- singularFalseDeclare that the system is singular; use when solving curl-curl problem if mass term is zero
Default:False
C++ Type:bool
Controllable:No
Description:Declare that the system is singular; use when solving curl-curl problem if mass term is zero
MFEMHypreAMS
Overview
Defines and builds an mfem::HypreAMS solver to use as a preconditioner or solver to solve the MFEM equation system. Most effective for preconditioning and solving a curl-curl problem when using Nédélec elements, in which case the FE space should be passed to the mfem::HypreAMS solver during construction.
If the system of equations is singular - commonly arising, for example, when solving for the magnetic vector potential in magnetostatic systems in the steady state - users should set the singular parameter to true to add a small mass term to ensure solvability.
A Low-Order-Refined (LOR) version of this solver may be used instead by setting the parameter "low_order_refined" to true. Using an LOR solver improves performance for high polynomial order systems.
Example Input File Syntax
[FESpaces<<<{"href": "../../../syntax/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": "../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
[]
[]
[Preconditioner<<<{"href": "../../../syntax/Preconditioner/index.html"}>>>]
[ams]
type = MFEMHypreAMS<<<{"description": "Hypre auxiliary-space Maxwell solver and preconditioner for the iterative solution of MFEM equation systems.", "href": "MFEMHypreAMS.html"}>>>
fespace<<<{"description": "H(curl) FESpace to use in HypreAMS setup."}>>> = HCurlFESpace
[]
[]
[Solver<<<{"href": "../../../syntax/Solver/index.html"}>>>]
type = MFEMHypreGMRES
preconditioner = ams
l_tol = 1e-12
[](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/submeshes/magnetostatic.i)
- (test/tests/mfem/submeshes/hphi_magnetodynamic.i)
- (test/tests/mfem/vectorpostprocessors/line_value_sampler/line_value_sampler_curlcurl.i)
- (test/tests/mfem/kernels/maxwell_eigenproblem.i)
- (test/tests/mfem/kernels/curlcurl.i)
- (test/tests/mfem/vectorpostprocessors/point_value_sampler/point_value_sampler_curlcurl.i)
- (test/tests/mfem/submeshes/av_magnetostatic.i)
low_order_refined
Default:False
C++ Type:bool
Controllable:No
Description:Set usage of Low-Order Refined solver.
(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/submeshes/magnetostatic.i)
# Definite Maxwell problem solved with Nedelec elements of the first kind
# based on MFEM Example 3.
[Mesh]
type = MFEMMesh
file = ../mesh/cylinder-hex-q2.gen
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
[]
[HDivFESpace]
type = MFEMVectorFESpace
fec_type = RT
fec_order = CONSTANT
[]
[]
[Variables]
[a_field]
type = MFEMVariable
fespace = HCurlFESpace
[]
[electric_potential]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[AuxVariables]
[b_field]
type = MFEMVariable
fespace = HDivFESpace
[]
[]
[AuxKernels]
[curl]
type = MFEMCurlAux
variable = b_field
source = a_field
execute_on = TIMESTEP_END
[]
[]
[BCs]
[tangential_a_bdr]
type = MFEMVectorTangentialDirichletBC
variable = a_field
boundary = '1 2 3'
[]
[]
[Kernels]
inactive = coefficient_source
[curlcurl]
type = MFEMCurlCurlKernel
variable = a_field
[]
[auxvar_source]
type = MFEMMixedVectorGradientKernel
trial_variable = electric_potential
variable = a_field
block = 1
[]
[coefficient_source]
type = MFEMVectorFEDomainLFKernel
vector_coefficient = electric_potential_grad
variable = a_field
block = 1
[]
[]
[Preconditioner]
[ams]
type = MFEMHypreAMS
fespace = HCurlFESpace
singular = true
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = ams
l_tol = 1e-12
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[MultiApps]
[subapp]
type = FullSolveMultiApp
input_files = open_coil_source.i
execute_on = INITIAL
[]
[]
[Transfers]
[from_sub]
type = MultiAppMFEMCopyTransfer
source_variables = electric_potential
variables = electric_potential
from_multi_app = subapp
[]
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/Magnetostatic
vtk_format = ASCII
[]
[]
(test/tests/mfem/submeshes/hphi_magnetodynamic.i)
# Solve for the magnetic field around a closed conductor subject to
# global current constraint.
conductor_domains = 'TorusCore TorusSheath'
conductor_resistivity = 1.0
vacuum_permeability = 1.0
[Problem]
type = MFEMProblem
[]
[Mesh]
type = MFEMMesh
file = ../mesh/split_embedded_concentric_torus.e
[]
[FunctorMaterials]
[Conductor]
type = MFEMGenericFunctorMaterial
prop_names = 'resistivity'
prop_values = ${conductor_resistivity}
block = ${conductor_domains}
[]
[Vacuum]
type = MFEMGenericFunctorMaterial
prop_names = 'permeability'
prop_values = '${vacuum_permeability}'
[]
[]
[SubMeshes]
[conductor]
type = MFEMDomainSubMesh
block = ${conductor_domains}
submesh_boundary = conductor_surface
[]
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
[]
[HDivFESpace]
type = MFEMVectorFESpace
fec_type = RT
fec_order = CONSTANT
[]
[CoilHCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
submesh = conductor
[]
[]
[Variables]
[coil_induced_h_field]
type = MFEMVariable
fespace = CoilHCurlFESpace
[]
[]
[AuxVariables]
[h_field]
type = MFEMVariable
fespace = HCurlFESpace
[]
[coil_external_h_field]
type = MFEMVariable
fespace = CoilHCurlFESpace
[]
[j_field]
type = MFEMVariable
fespace = HDivFESpace
[]
[]
[AuxKernels]
[update_j_field]
type = MFEMCurlAux
variable = j_field
source = h_field
scale_factor = 1.0
execute_on = TIMESTEP_END
[]
[]
[BCs]
[conductor_bdr]
type = MFEMVectorTangentialDirichletBC
variable = coil_induced_h_field
vector_coefficient = coil_external_h_field
boundary = conductor_surface
[]
[]
[Kernels]
[dBdt]
type = MFEMTimeDerivativeVectorFEMassKernel
variable = coil_induced_h_field
coefficient = permeability
[]
[curlE]
type = MFEMCurlCurlKernel
variable = coil_induced_h_field
coefficient = resistivity
[]
[]
[Preconditioner]
[ams]
type = MFEMHypreAMS
fespace = CoilHCurlFESpace
[]
[]
[Solver]
type = MFEMHyprePCG
preconditioner = ams
l_tol = 1e-9
l_max_its = 100
[]
[Executioner]
type = MFEMTransient
dt = 0.5
start_time = 0.0
end_time = 2.0
[]
[MultiApps]
[hphi_magnetostatic]
type = FullSolveMultiApp
input_files = hphi_magnetostatic.i
execute_on = INITIAL
[]
[]
[Transfers]
[from_external_field]
type = MultiAppMFEMCopyTransfer
source_variables = h_field
variables = h_field
from_multi_app = hphi_magnetostatic
[]
[submesh_transfer_to_coil]
type = MFEMSubMeshTransfer
from_variable = h_field
to_variable = coil_external_h_field
execute_on = TIMESTEP_BEGIN
[]
[submesh_transfer_from_coil]
type = MFEMSubMeshTransfer
from_variable = coil_induced_h_field
to_variable = h_field
execute_on = TIMESTEP_END
[]
[]
[Postprocessors]
[CoilPower]
type = MFEMVectorFEInnerProductIntegralPostprocessor
coefficient = resistivity
dual_variable = j_field
primal_variable = j_field
block = 'TorusCore TorusSheath'
[]
[]
[Outputs]
[ReportedPostprocessors]
type = CSV
file_base = OutputData/HPhiMagnetodynamicClosedCoilCSV
[]
[VacuumParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/HPhiMagnetodynamicClosedCoil
vtk_format = ASCII
[]
[]
(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
[]
(test/tests/mfem/kernels/maxwell_eigenproblem.i)
[Mesh]
type = MFEMMesh
file = ../mesh/beam-tet.mesh
dim = 3
[]
[Problem]
type = MFEMEigenproblem
num_modes = 5
[]
[FESpaces]
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
[]
[]
[Variables]
[E]
type = MFEMVariable
fespace = HCurlFESpace
[]
[]
[BCs]
[all]
type = MFEMVectorTangentialDirichletBC
variable = E
vector_coefficient = '0 0 0'
[]
[]
[Kernels]
[diff]
type = MFEMCurlCurlKernel
variable = E
[]
[]
[Preconditioner]
[ams]
type = MFEMHypreAMS
fespace = HCurlFESpace
print_level = 0
singular = true
[]
[]
[Solver]
type = MFEMHypreAME
preconditioner = ams
print_level = 0
l_tol = 1e-8
l_max_its = 100
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[VectorPostprocessors]
[eigenvalues]
type = MFEMEigenvaluesPostprocessor
[]
[]
[Outputs]
execute_on = 'timestep_end'
csv = true
file_base = OutputData/MaxwellEigenproblem
[]
(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/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
[]
[]