Magnetostatic Problem on a Topologically Closed Conductor
Summary
Solves for the magnetic field around a topologically closed toroidal conductor carrying a net current, using a magnetic vector potential discretized using conforming Nédélec elements, following the method described in P. Dular, International Compumag Society Newsletter, 7, no. 2 (2000): 4-7. for enforcing global constraints.
Description
This problem solves the magnetostatic problem with strong form:
with the constitutive relations
where is the magnetic flux density, is the magnetic field, is the current density, is the electric field, is the externally applied electric field, is the material permeability, is the material resistivity, and is the material conductivity. To enforce the divergence free condition on , we shall also introduce the magnetic vector potential defined such that .
In addition, we wish to impose a global voltage constraint on a toroidal conductor domain of the form
where is any path traversing around the conductor lying entirely within with winding number equal to one, and is the loop voltage imposed on the conductor.
In a simply connected domain, satisfying can be achieved by introducing a continuous single-valued scalar electric potential such that . However, in order for the electric field to be nonzero circulating the closed conductor, this is no longer sufficient, since . We shall therefore represent the electric field by
With this, we can arrive at the weak form for we intend to solve, of the form
where
Discretizing the External Electric Field
We must now find a suitable choice of basis for the external electric field, such that
for some unique set of degrees of freedom . We shall also wish to restrict the support of to the smallest number of elements possible in to minimise the number of degrees of freedom required to represent .
To do this, we first identify a 'cut' surface of the toroidal conductor that is required to make it simply connected. Next, we define a one element wide 'transition' region on one side of the cut surface, using an MFEMCutTransitionSubMesh
object, comprised of all elements with at least one vertex lying on the cut surface that lie on one side of the cut.
Next, to enforce in this region, we define a scalar field variable defined on a (nodal) conforming finite element space on this submesh such that :
so .
Finally, to impose the global voltage constraint , we identify that
since outside . The loop voltage constraint is thus strongly imposed by Dirichlet conditions on , setting on and on . For lowest order conforming representations of , this fully constrains in the transition region.
This is shown explicitly in the example file:
# Solve for the electric field on a closed conductor subject to
# global loop voltage constraint.
initial_coil_domains = 'TorusCore TorusSheath'
coil_cut_surface = 'Cut'
coil_loop_voltage = -1.0
coil_conductivity = 1.0
[Problem<<<{"href": "../Problem/index.html"}>>>]
type = MFEMProblem
[]
[Mesh<<<{"href": "../Mesh/index.html"}>>>]
type = MFEMMesh
file = ../mesh/embedded_concentric_torus.e
[]
[FunctorMaterials<<<{"href": "../FunctorMaterials/index.html"}>>>]
[Conductor]
type = MFEMGenericFunctorMaterial<<<{"description": "Declares material scalar properties based on names and coefficients prescribed by input parameters.", "href": "../../source/mfem/functormaterials/MFEMGenericFunctorMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = conductivity
prop_values<<<{"description": "The corresponding names of coefficients associated with the named properties. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = ${coil_conductivity}
[]
[]
[ICs<<<{"href": "../ICs/index.html"}>>>]
[coil_external_potential_ic]
type = MFEMScalarBoundaryIC<<<{"description": "Sets the initial values of an MFEM scalar variable from a user-specified scalar coefficient.", "href": "../../source/mfem/ics/MFEMScalarBoundaryIC.html"}>>>
variable<<<{"description": "The variable to apply the initial condition on."}>>> = coil_external_potential
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies. Defaults to all boundaries."}>>> = ${coil_cut_surface}
coefficient<<<{"description": "The scalar coefficient. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = ${coil_loop_voltage}
[]
[]
[SubMeshes<<<{"href": "../SubMeshes/index.html"}>>>]
[cut]
type = MFEMCutTransitionSubMesh<<<{"description": "Class to construct an MFEMSubMesh formed from the set of elements that have least one vertex on the specified cut plane, that lie on one side of the plane, and that are restricted to the set of user-specified subdomains.", "href": "../../source/mfem/submeshes/MFEMCutTransitionSubMesh.html"}>>>
cut_boundary<<<{"description": "The boundary associated with the mesh cut."}>>> = ${coil_cut_surface}
block<<<{"description": "The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains."}>>> = ${initial_coil_domains}
transition_subdomain<<<{"description": "The name of the subdomain to be created on the mesh comprised of the set of elements adjacent to the cut surface on one side."}>>> = transition_dom
transition_subdomain_boundary<<<{"description": "Name to assign boundary of transition subdomain not shared with cut surface."}>>> = transition_bdr
closed_subdomain<<<{"description": "The name of the subdomain attribute to be created comprised of the set of all elements of the closed geometry, including the new transition region."}>>> = coil_dom
[]
[coil]
type = MFEMDomainSubMesh<<<{"description": "Class to construct an MFEMSubMesh formed from the subspace of the parent mesh restricted to the set of user-specified subdomains.", "href": "../../source/mfem/submeshes/MFEMDomainSubMesh.html"}>>>
block<<<{"description": "The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains."}>>> = coil_dom
execution_order_group<<<{"description": "Execution order groups are executed in increasing order (e.g., the lowest number is executed first). Note that negative group numbers may be used to execute groups before the default (0) group. Please refer to the user object documentation for ordering of user object execution within a group."}>>> = 2
[]
[]
[FESpaces<<<{"href": "../FESpaces/index.html"}>>>]
[H1FESpace]
type = MFEMScalarFESpace<<<{"description": "Convenience class to construct scalar finite element spaces.", "href": "../../source/mfem/fespaces/MFEMScalarFESpace.html"}>>>
fec_type<<<{"description": "Specifies the family of FE shape functions."}>>> = H1
fec_order<<<{"description": "Order of the FE shape function to use."}>>> = FIRST
[]
[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
[]
[CoilH1FESpace]
type = MFEMScalarFESpace<<<{"description": "Convenience class to construct scalar finite element spaces.", "href": "../../source/mfem/fespaces/MFEMScalarFESpace.html"}>>>
fec_type<<<{"description": "Specifies the family of FE shape functions."}>>> = H1
fec_order<<<{"description": "Order of the FE shape function to use."}>>> = FIRST
submesh<<<{"description": "Submesh to define the FESpace on. Leave blank to use base mesh."}>>> = coil
[]
[CoilHCurlFESpace]
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
submesh<<<{"description": "Submesh to define the FESpace on. Leave blank to use base mesh."}>>> = coil
[]
[TransitionH1FESpace]
type = MFEMScalarFESpace<<<{"description": "Convenience class to construct scalar finite element spaces.", "href": "../../source/mfem/fespaces/MFEMScalarFESpace.html"}>>>
fec_type<<<{"description": "Specifies the family of FE shape functions."}>>> = H1
fec_order<<<{"description": "Order of the FE shape function to use."}>>> = FIRST
submesh<<<{"description": "Submesh to define the FESpace on. Leave blank to use base mesh."}>>> = cut
[]
[TransitionHCurlFESpace]
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
submesh<<<{"description": "Submesh to define the FESpace on. Leave blank to use base mesh."}>>> = cut
[]
[]
[Variables<<<{"href": "../Variables/index.html"}>>>]
[coil_induced_potential]
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."}>>> = CoilH1FESpace
[]
[]
[AuxVariables<<<{"href": "../AuxVariables/index.html"}>>>]
[coil_external_potential]
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."}>>> = CoilH1FESpace
[]
[transition_external_potential]
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."}>>> = TransitionH1FESpace
[]
[transition_external_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."}>>> = TransitionHCurlFESpace
[]
[induced_potential]
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."}>>> = H1FESpace
[]
[induced_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
[]
[external_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
[]
[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
[]
[]
[AuxKernels<<<{"href": "../AuxKernels/index.html"}>>>]
[update_induced_e_field]
type = MFEMGradAux<<<{"description": "Calculates the gradient of an H1 conforming source variable and stores the result on an H(curl) conforming ND result auxvariable", "href": "../../source/mfem/auxkernels/MFEMGradAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = induced_e_field
source<<<{"description": "Scalar H1 MFEMVariable to take the gradient of."}>>> = induced_potential
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
[]
[update_external_e_field]
type = MFEMGradAux<<<{"description": "Calculates the gradient of an H1 conforming source variable and stores the result on an H(curl) conforming ND result auxvariable", "href": "../../source/mfem/auxkernels/MFEMGradAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = transition_external_e_field
source<<<{"description": "Scalar H1 MFEMVariable to take the gradient of."}>>> = transition_external_potential
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
[]
[update_total_e_field]
type = MFEMSumAux<<<{"description": "Calculates the sum of two variables sharing an FE space, each optionally scaled by a real constant, and stores the result in a third.", "href": "../../source/mfem/auxkernels/MFEMSumAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = e_field
first_source_variable<<<{"description": "First variable to sum."}>>> = induced_e_field
second_source_variable<<<{"description": "Second variable to sum."}>>> = external_e_field
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
execution_order_group<<<{"description": "Execution order groups are executed in increasing order (e.g., the lowest number is executed first). Note that negative group numbers may be used to execute groups before the default (0) group. Please refer to the user object documentation for ordering of user object execution within a group."}>>> = 3
[]
[]
[Kernels<<<{"href": "../Kernels/index.html"}>>>]
[diff]
type = MFEMDiffusionKernel<<<{"description": "Adds the domain integrator to an MFEM problem for the bilinear form $(k\\vec\\nabla u, \\vec\\nabla v)_\\Omega$ arising from the weak form of the Laplacian operator $- \\vec\\nabla \\cdot \\left( k \\vec \\nabla u \\right)$.", "href": "../../source/mfem/kernels/MFEMDiffusionKernel.html"}>>>
variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = coil_induced_potential
coefficient<<<{"description": "Name of property for diffusion coefficient k. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = conductivity
[]
[source]
type = MFEMMixedGradGradKernel<<<{"description": "Adds the domain integrator to an MFEM problem for the mixed bilinear form $(k\\vec\\nabla u, \\vec\\nabla v)_\\Omega$.", "href": "../../source/mfem/kernels/MFEMMixedGradGradKernel.html"}>>>
trial_variable<<<{"description": "The trial variable this kernel is acting on and which will be solved for. If empty (default), it will be the same as the test variable."}>>> = coil_external_potential
variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = coil_induced_potential
coefficient<<<{"description": "Name of property k to use. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = conductivity
block<<<{"description": "The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains."}>>> = 'transition_dom'
[]
[]
[Solver<<<{"href": "../Solver/index.html"}>>>]
type = MFEMSuperLU
[]
[Executioner<<<{"href": "../Executioner/index.html"}>>>]
type = MFEMSteady
[]
[Transfers<<<{"href": "../Transfers/index.html"}>>>]
[submesh_transfer_from_coil]
type = MFEMSubMeshTransfer<<<{"description": "Class to transfer MFEM variable data to or from a restricted copy of the variable defined on an a subspace of an MFEMMesh, represented as an MFEMSubMesh.", "href": "../../source/mfem/transfers/MFEMSubMeshTransfer.html"}>>>
from_variable<<<{"description": "MFEM variable to transfer data from. Can be defined on either the parent mesh or a submesh of it."}>>> = coil_induced_potential
to_variable<<<{"description": "MFEM variable to transfer data into. Can be defined on either the parent mesh or a submesh of it."}>>> = induced_potential
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
[]
[submesh_transfer_to_transition]
type = MFEMSubMeshTransfer<<<{"description": "Class to transfer MFEM variable data to or from a restricted copy of the variable defined on an a subspace of an MFEMMesh, represented as an MFEMSubMesh.", "href": "../../source/mfem/transfers/MFEMSubMeshTransfer.html"}>>>
from_variable<<<{"description": "MFEM variable to transfer data from. Can be defined on either the parent mesh or a submesh of it."}>>> = coil_external_potential
to_variable<<<{"description": "MFEM variable to transfer data into. Can be defined on either the parent mesh or a submesh of it."}>>> = transition_external_potential
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
[]
[submesh_transfer_from_transition]
type = MFEMSubMeshTransfer<<<{"description": "Class to transfer MFEM variable data to or from a restricted copy of the variable defined on an a subspace of an MFEMMesh, represented as an MFEMSubMesh.", "href": "../../source/mfem/transfers/MFEMSubMeshTransfer.html"}>>>
from_variable<<<{"description": "MFEM variable to transfer data from. Can be defined on either the parent mesh or a submesh of it."}>>> = transition_external_e_field
to_variable<<<{"description": "MFEM variable to transfer data into. Can be defined on either the parent mesh or a submesh of it."}>>> = external_e_field
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
execution_order_group<<<{"description": "Execution order groups are executed in increasing order (e.g., the lowest number is executed first). Note that negative group numbers may be used to execute groups before the default (0) group. Please refer to the user object documentation for ordering of user object execution within a group."}>>> = 2
[]
[]
[Outputs<<<{"href": "../Outputs/index.html"}>>>]
[GlobalParaViewDataCollection]
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/WholePotentialCoil
vtk_format<<<{"description": "Select VTK data format to use, choosing between BINARY, BINARY32, and ASCII."}>>> = ASCII
[]
[TransitionParaViewDataCollection]
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/CutPotentialCoil
vtk_format<<<{"description": "Select VTK data format to use, choosing between BINARY, BINARY32, and ASCII."}>>> = ASCII
submesh<<<{"description": "Submesh to output variables on. Leave blank to use base mesh."}>>> = cut
[]
[CoilParaViewDataCollection]
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/Coil
vtk_format<<<{"description": "Select VTK data format to use, choosing between BINARY, BINARY32, and ASCII."}>>> = ASCII
submesh<<<{"description": "Submesh to output variables on. Leave blank to use base mesh."}>>> = coil
[]
[]
(test/tests/mfem/submeshes/cut_closed_coil.i)Solving for the Total Electric Field
With now constrained, we now must solve for the remaining degrees of freedom in arising due to surface charges in the coil as required for continuity of current in the conductor (), by solving
for . We then can calculate the total electric field from the sum of its induced and external components.
Solving for the Magnetic Flux Density
Finally, with the electric field in the conductor, we can now find the magnetic flux density around the current density flowing in the coil in the magnetostatic limit.
This is typically solved using HypreAMS
as a preconditioner. Note that the above magnetostatic problem is singular in non-conductive regions, and thus either singular = true
should be passed to the preconditioner or a small mass term with should be added for stability.
Example File
The full closed coil magnetostatic example detailed above can be found below:
# Magnetostatic problem solved on a closed conductor subject to
# global loop voltage constraint.
[Mesh<<<{"href": "../Mesh/index.html"}>>>]
type = MFEMMesh
file = ../mesh/embedded_concentric_torus.e
[]
[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"}>>>]
[a_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"}>>>]
[b_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
[]
[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
[]
[]
[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"}>>> = b_field
source<<<{"description": "Vector H(curl) MFEMVariable to take the curl of."}>>> = a_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
execution_order_group<<<{"description": "Execution order groups are executed in increasing order (e.g., the lowest number is executed first). Note that negative group numbers may be used to execute groups before the default (0) group. Please refer to the user object documentation for ordering of user object execution within a group."}>>> = 3
[]
[]
[Functions<<<{"href": "../Functions/index.html"}>>>]
[exact_a_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."}>>> = '0'
expression_y<<<{"description": "y-component of function."}>>> = '0'
expression_z<<<{"description": "z-component of function."}>>> = '0'
[]
[]
[BCs<<<{"href": "../BCs/index.html"}>>>]
[tangential_a_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"}>>> = a_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_a_field
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies. Defaults to all boundaries."}>>> = 'Exterior'
[]
[]
[FunctorMaterials<<<{"href": "../FunctorMaterials/index.html"}>>>]
[Vacuum]
type = MFEMGenericFunctorMaterial<<<{"description": "Declares material scalar properties based on names and coefficients prescribed by input parameters.", "href": "../../source/mfem/functormaterials/MFEMGenericFunctorMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = reluctivity
prop_values<<<{"description": "The corresponding names of coefficients associated with the named properties. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = 1.0
[]
[Conductor]
type = MFEMGenericFunctorMaterial<<<{"description": "Declares material scalar properties based on names and coefficients prescribed by input parameters.", "href": "../../source/mfem/functormaterials/MFEMGenericFunctorMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = conductivity
prop_values<<<{"description": "The corresponding names of coefficients associated with the named properties. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = 1.0
block<<<{"description": "The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains."}>>> = 'TorusCore TorusSheath'
[]
[]
[Kernels<<<{"href": "../Kernels/index.html"}>>>]
[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"}>>> = a_field
coefficient<<<{"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."}>>> = 1e-10
[]
[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"}>>> = a_field
coefficient<<<{"description": "Name of scalar coefficient 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."}>>> = reluctivity
[]
[source]
type = MFEMMixedVectorMassKernel<<<{"description": "Adds the domain integrator to an MFEM problem for the mixed bilinear form $(k\\vec u, \\vec v)_\\Omega$.", "href": "../../source/mfem/kernels/MFEMMixedVectorMassKernel.html"}>>>
variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = a_field
trial_variable<<<{"description": "The trial variable this kernel is acting on and which will be solved for. If empty (default), it will be the same as the test variable."}>>> = e_field
coefficient<<<{"description": "Name of property k to use. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = conductivity
block<<<{"description": "The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains."}>>> = 'TorusCore TorusSheath'
[]
[]
[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 = MFEMHyprePCG
preconditioner = ams
l_tol = 1e-14
l_max_its = 1000
[]
[Executioner<<<{"href": "../Executioner/index.html"}>>>]
type = MFEMSteady
device = cpu
[]
[MultiApps<<<{"href": "../MultiApps/index.html"}>>>]
[coil]
type = FullSolveMultiApp<<<{"description": "Performs a complete simulation during each execution.", "href": "../../source/multiapps/FullSolveMultiApp.html"}>>>
input_files<<<{"description": "The input file for each App. If this parameter only contains one input file it will be used for all of the Apps. When using 'positions_from_file' it is also admissable to provide one input_file per file."}>>> = cut_closed_coil.i
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."}>>> = INITIAL
[]
[]
[Transfers<<<{"href": "../Transfers/index.html"}>>>]
[from_coil]
type = MultiAppMFEMCopyTransfer<<<{"description": "Copies variable values from one MFEM application to another", "href": "../../source/mfem/transfers/MultiAppMFEMCopyTransfer.html"}>>>
source_variable<<<{"description": "Variable to transfer from"}>>> = e_field
variable<<<{"description": "AuxVariable to store transferred value in."}>>> = e_field
from_multi_app<<<{"description": "The name of the MultiApp to receive data from"}>>> = coil
[]
[]
[Postprocessors<<<{"href": "../Postprocessors/index.html"}>>>]
[CoilPower]
type = MFEMVectorFEInnerProductIntegralPostprocessor<<<{"description": "Calculates the total flux of a vector field through an interface", "href": "../../source/mfem/postprocessors/MFEMVectorFEInnerProductIntegralPostprocessor.html"}>>>
coefficient<<<{"description": "Name of optional scalar coefficient to scale integrand by. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = conductivity
dual_variable<<<{"description": "Name of the second vector variable in the inner product."}>>> = e_field
primal_variable<<<{"description": "Name of the first vector variable in the inner product."}>>> = e_field
block<<<{"description": "The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains."}>>> = 'TorusCore TorusSheath'
[]
[]
[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/MagnetostaticClosedCoil
vtk_format<<<{"description": "Select VTK data format to use, choosing between BINARY, BINARY32, and ASCII."}>>> = ASCII
[]
[ReportedPostprocessors]
type = CSV<<<{"description": "Output for postprocessors, vector postprocessors, and scalar variables using comma seperated values (CSV).", "href": "../../source/outputs/CSV.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/AVMagnetostaticClosedCoilCSV
[]
[]
(test/tests/mfem/submeshes/av_magnetostatic.i)Global Current Constraints
Strongly constraining the total current through the conductor instead of the loop voltage is also possible, by defining a cut plane in the domain outside the conductor region using the H-phi formulation instead of the A-V formulation detailed above; the details are similar, using a magnetic scalar potential instead of an electric scalar potential. The constraint to be enforced is
where is any path encircling the conductor, and I is the total current passing through the path . Outside the conductor, the magnetic field is curl-free, and thus can be expressed by the gradient of a scalar magnetic potential (up to a sign).
An example file for this approach can be found below:
# Solve for the magnetic field around a closed conductor subject to
# global current constraint.
initial_vacuum_domains = 'Exterior'
vacuum_cut_surface = 'Cut'
conductor_current = 1.0
vacuum_permeability = 1.0
[Problem<<<{"href": "../Problem/index.html"}>>>]
type = MFEMProblem
[]
[Mesh<<<{"href": "../Mesh/index.html"}>>>]
type = MFEMMesh
file = ../mesh/split_embedded_concentric_torus.e
[]
[FunctorMaterials<<<{"href": "../FunctorMaterials/index.html"}>>>]
[Conductor]
type = MFEMGenericFunctorMaterial<<<{"description": "Declares material scalar properties based on names and coefficients prescribed by input parameters.", "href": "../../source/mfem/functormaterials/MFEMGenericFunctorMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = permeability
prop_values<<<{"description": "The corresponding names of coefficients associated with the named properties. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = ${vacuum_permeability}
[]
[]
[ICs<<<{"href": "../ICs/index.html"}>>>]
[vacuum_cut_potential_ic]
type = MFEMScalarBoundaryIC<<<{"description": "Sets the initial values of an MFEM scalar variable from a user-specified scalar coefficient.", "href": "../../source/mfem/ics/MFEMScalarBoundaryIC.html"}>>>
variable<<<{"description": "The variable to apply the initial condition on."}>>> = vacuum_cut_potential
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies. Defaults to all boundaries."}>>> = ${vacuum_cut_surface}
coefficient<<<{"description": "The scalar coefficient. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = ${conductor_current}
[]
[]
[SubMeshes<<<{"href": "../SubMeshes/index.html"}>>>]
[cut]
type = MFEMCutTransitionSubMesh<<<{"description": "Class to construct an MFEMSubMesh formed from the set of elements that have least one vertex on the specified cut plane, that lie on one side of the plane, and that are restricted to the set of user-specified subdomains.", "href": "../../source/mfem/submeshes/MFEMCutTransitionSubMesh.html"}>>>
cut_boundary<<<{"description": "The boundary associated with the mesh cut."}>>> = ${vacuum_cut_surface}
block<<<{"description": "The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains."}>>> = ${initial_vacuum_domains}
transition_subdomain<<<{"description": "The name of the subdomain to be created on the mesh comprised of the set of elements adjacent to the cut surface on one side."}>>> = transition_dom
transition_subdomain_boundary<<<{"description": "Name to assign boundary of transition subdomain not shared with cut surface."}>>> = transition_bdr
closed_subdomain<<<{"description": "The name of the subdomain attribute to be created comprised of the set of all elements of the closed geometry, including the new transition region."}>>> = vacuum_dom
[]
[vacuum]
type = MFEMDomainSubMesh<<<{"description": "Class to construct an MFEMSubMesh formed from the subspace of the parent mesh restricted to the set of user-specified subdomains.", "href": "../../source/mfem/submeshes/MFEMDomainSubMesh.html"}>>>
block<<<{"description": "The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains."}>>> = vacuum_dom
execution_order_group<<<{"description": "Execution order groups are executed in increasing order (e.g., the lowest number is executed first). Note that negative group numbers may be used to execute groups before the default (0) group. Please refer to the user object documentation for ordering of user object execution within a group."}>>> = 2
[]
[]
[FESpaces<<<{"href": "../FESpaces/index.html"}>>>]
[VacuumH1FESpace]
type = MFEMScalarFESpace<<<{"description": "Convenience class to construct scalar finite element spaces.", "href": "../../source/mfem/fespaces/MFEMScalarFESpace.html"}>>>
fec_type<<<{"description": "Specifies the family of FE shape functions."}>>> = H1
fec_order<<<{"description": "Order of the FE shape function to use."}>>> = FIRST
submesh<<<{"description": "Submesh to define the FESpace on. Leave blank to use base mesh."}>>> = vacuum
[]
[VacuumHCurlFESpace]
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
submesh<<<{"description": "Submesh to define the FESpace on. Leave blank to use base mesh."}>>> = vacuum
[]
[TransitionH1FESpace]
type = MFEMScalarFESpace<<<{"description": "Convenience class to construct scalar finite element spaces.", "href": "../../source/mfem/fespaces/MFEMScalarFESpace.html"}>>>
fec_type<<<{"description": "Specifies the family of FE shape functions."}>>> = H1
fec_order<<<{"description": "Order of the FE shape function to use."}>>> = FIRST
submesh<<<{"description": "Submesh to define the FESpace on. Leave blank to use base mesh."}>>> = cut
[]
[TransitionHCurlFESpace]
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
submesh<<<{"description": "Submesh to define the FESpace on. Leave blank to use base mesh."}>>> = cut
[]
[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
[]
[]
[Variables<<<{"href": "../Variables/index.html"}>>>]
[vacuum_magnetic_potential]
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."}>>> = VacuumH1FESpace
[]
[]
[AuxVariables<<<{"href": "../AuxVariables/index.html"}>>>]
[vacuum_cut_potential]
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."}>>> = VacuumH1FESpace
[]
[transition_cut_potential]
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."}>>> = TransitionH1FESpace
[]
[transition_cut_function_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."}>>> = TransitionHCurlFESpace
[]
[background_h_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."}>>> = VacuumHCurlFESpace
[]
[cut_function_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."}>>> = VacuumHCurlFESpace
[]
[vacuum_h_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."}>>> = VacuumHCurlFESpace
[]
[h_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
[]
[]
[AuxKernels<<<{"href": "../AuxKernels/index.html"}>>>]
[update_background_h_field]
type = MFEMGradAux<<<{"description": "Calculates the gradient of an H1 conforming source variable and stores the result on an H(curl) conforming ND result auxvariable", "href": "../../source/mfem/auxkernels/MFEMGradAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = background_h_field
source<<<{"description": "Scalar H1 MFEMVariable to take the gradient of."}>>> = vacuum_magnetic_potential
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
[]
[update_transition_cut_function_field]
type = MFEMGradAux<<<{"description": "Calculates the gradient of an H1 conforming source variable and stores the result on an H(curl) conforming ND result auxvariable", "href": "../../source/mfem/auxkernels/MFEMGradAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = transition_cut_function_field
source<<<{"description": "Scalar H1 MFEMVariable to take the gradient of."}>>> = transition_cut_potential
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
[]
[update_total_h_field]
type = MFEMSumAux<<<{"description": "Calculates the sum of two variables sharing an FE space, each optionally scaled by a real constant, and stores the result in a third.", "href": "../../source/mfem/auxkernels/MFEMSumAux.html"}>>>
variable<<<{"description": "The name of the variable that this object applies to"}>>> = vacuum_h_field
first_source_variable<<<{"description": "First variable to sum."}>>> = background_h_field
second_source_variable<<<{"description": "Second variable to sum."}>>> = cut_function_field
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
execution_order_group<<<{"description": "Execution order groups are executed in increasing order (e.g., the lowest number is executed first). Note that negative group numbers may be used to execute groups before the default (0) group. Please refer to the user object documentation for ordering of user object execution within a group."}>>> = 3
[]
[]
[BCs<<<{"href": "../BCs/index.html"}>>>]
# Set zero of magnetic potential on symmetry plane
[Exterior]
type = MFEMScalarDirichletBC<<<{"description": "Applies a Dirichlet condition to a scalar variable.", "href": "../../source/mfem/bcs/MFEMScalarDirichletBC.html"}>>>
variable<<<{"description": "Variable on which to apply the boundary condition"}>>> = vacuum_magnetic_potential
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies. Defaults to all boundaries."}>>> = 'Cut'
coefficient<<<{"description": "The coefficient setting the values on the essential boundary. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = 0.0
[]
[]
[Kernels<<<{"href": "../Kernels/index.html"}>>>]
[diff]
type = MFEMDiffusionKernel<<<{"description": "Adds the domain integrator to an MFEM problem for the bilinear form $(k\\vec\\nabla u, \\vec\\nabla v)_\\Omega$ arising from the weak form of the Laplacian operator $- \\vec\\nabla \\cdot \\left( k \\vec \\nabla u \\right)$.", "href": "../../source/mfem/kernels/MFEMDiffusionKernel.html"}>>>
variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = vacuum_magnetic_potential
coefficient<<<{"description": "Name of property for diffusion coefficient k. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = permeability
[]
[source]
type = MFEMMixedGradGradKernel<<<{"description": "Adds the domain integrator to an MFEM problem for the mixed bilinear form $(k\\vec\\nabla u, \\vec\\nabla v)_\\Omega$.", "href": "../../source/mfem/kernels/MFEMMixedGradGradKernel.html"}>>>
trial_variable<<<{"description": "The trial variable this kernel is acting on and which will be solved for. If empty (default), it will be the same as the test variable."}>>> = vacuum_cut_potential
variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = vacuum_magnetic_potential
coefficient<<<{"description": "Name of property k to use. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = permeability
block<<<{"description": "The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains."}>>> = 'transition_dom'
[]
[]
[Preconditioner<<<{"href": "../Preconditioner/index.html"}>>>]
[boomeramg]
type = MFEMHypreBoomerAMG<<<{"description": "Hypre BoomerAMG solver and preconditioner for the iterative solution of MFEM equation systems.", "href": "../../source/mfem/solvers/MFEMHypreBoomerAMG.html"}>>>
[]
[]
[Solver<<<{"href": "../Solver/index.html"}>>>]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 100
[]
[Executioner<<<{"href": "../Executioner/index.html"}>>>]
type = MFEMSteady
[]
[Transfers<<<{"href": "../Transfers/index.html"}>>>]
[submesh_transfer_to_transition]
type = MFEMSubMeshTransfer<<<{"description": "Class to transfer MFEM variable data to or from a restricted copy of the variable defined on an a subspace of an MFEMMesh, represented as an MFEMSubMesh.", "href": "../../source/mfem/transfers/MFEMSubMeshTransfer.html"}>>>
from_variable<<<{"description": "MFEM variable to transfer data from. Can be defined on either the parent mesh or a submesh of it."}>>> = vacuum_cut_potential
to_variable<<<{"description": "MFEM variable to transfer data into. Can be defined on either the parent mesh or a submesh of it."}>>> = transition_cut_potential
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
[]
[submesh_transfer_from_transition]
type = MFEMSubMeshTransfer<<<{"description": "Class to transfer MFEM variable data to or from a restricted copy of the variable defined on an a subspace of an MFEMMesh, represented as an MFEMSubMesh.", "href": "../../source/mfem/transfers/MFEMSubMeshTransfer.html"}>>>
from_variable<<<{"description": "MFEM variable to transfer data from. Can be defined on either the parent mesh or a submesh of it."}>>> = transition_cut_function_field
to_variable<<<{"description": "MFEM variable to transfer data into. Can be defined on either the parent mesh or a submesh of it."}>>> = cut_function_field
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
execution_order_group<<<{"description": "Execution order groups are executed in increasing order (e.g., the lowest number is executed first). Note that negative group numbers may be used to execute groups before the default (0) group. Please refer to the user object documentation for ordering of user object execution within a group."}>>> = 2
[]
[submesh_transfer_from_vacuum]
type = MFEMSubMeshTransfer<<<{"description": "Class to transfer MFEM variable data to or from a restricted copy of the variable defined on an a subspace of an MFEMMesh, represented as an MFEMSubMesh.", "href": "../../source/mfem/transfers/MFEMSubMeshTransfer.html"}>>>
from_variable<<<{"description": "MFEM variable to transfer data from. Can be defined on either the parent mesh or a submesh of it."}>>> = vacuum_h_field
to_variable<<<{"description": "MFEM variable to transfer data into. Can be defined on either the parent mesh or a submesh of it."}>>> = h_field
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
execution_order_group<<<{"description": "Execution order groups are executed in increasing order (e.g., the lowest number is executed first). Note that negative group numbers may be used to execute groups before the default (0) group. Please refer to the user object documentation for ordering of user object execution within a group."}>>> = 4
[]
[]
[Postprocessors<<<{"href": "../Postprocessors/index.html"}>>>]
[MagneticEnergy]
type = MFEMVectorFEInnerProductIntegralPostprocessor<<<{"description": "Calculates the total flux of a vector field through an interface", "href": "../../source/mfem/postprocessors/MFEMVectorFEInnerProductIntegralPostprocessor.html"}>>>
coefficient<<<{"description": "Name of optional scalar coefficient to scale integrand by. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = ${fparse 0.5*vacuum_permeability}
dual_variable<<<{"description": "Name of the second vector variable in the inner product."}>>> = vacuum_h_field
primal_variable<<<{"description": "Name of the first vector variable in the inner product."}>>> = vacuum_h_field
execution_order_group<<<{"description": "Execution order groups are executed in increasing order (e.g., the lowest number is executed first). Note that negative group numbers may be used to execute groups before the default (0) group. Please refer to the user object documentation for ordering of user object execution within a group."}>>> = 4
block<<<{"description": "The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains."}>>> = 'Exterior'
[]
[]
[Outputs<<<{"href": "../Outputs/index.html"}>>>]
[ReportedPostprocessors]
type = CSV<<<{"description": "Output for postprocessors, vector postprocessors, and scalar variables using comma seperated values (CSV).", "href": "../../source/outputs/CSV.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/HPhiMagnetostaticClosedCoilCSV
[]
[VacuumParaViewDataCollection]
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/HPhiMagnetostaticClosedCoil
vtk_format<<<{"description": "Select VTK data format to use, choosing between BINARY, BINARY32, and ASCII."}>>> = ASCII
submesh<<<{"description": "Submesh to output variables on. Leave blank to use base mesh."}>>> = vacuum
[]
[]
(test/tests/mfem/submeshes/hphi_magnetostatic.i)