- basisGaussLobattoSpecifies the quadrature basis used for scalar elements.
Default:GaussLobatto
C++ Type:std::string
Controllable:No
Description:Specifies the quadrature basis used for scalar elements.
- fec_orderFIRSTOrder of the FE shape function to use.
Default:FIRST
C++ Type:MooseEnum
Controllable:No
Description:Order of the FE shape function to use.
- fec_typeH1Specifies the family of FE shape functions.
Default:H1
C++ Type:MooseEnum
Controllable:No
Description:Specifies the family of FE shape functions.
- orderingVDIMOrdering style to use for vector DoFs.
Default:VDIM
C++ Type:MooseEnum
Controllable:No
Description:Ordering style to use for vector DoFs.
- submeshSubmesh to define the FESpace on. Leave blank to use base mesh.
C++ Type:std::string
Controllable:No
Description:Submesh to define the FESpace on. Leave blank to use base mesh.
- vdim1The number of degrees of freedom per basis function.
Default:1
C++ Type:int
Controllable:No
Description:The number of degrees of freedom per basis function.
MFEMScalarFESpace
Summary
Convenience class to construct scalar finite element spaces.
Overview
This is a convenience class for building finite element spaces to represent scalar variables. The family of shape functions is selected from the fec_type
parameter, and the order is controlled using the fec_order
parameter.
If you need a finite element space that can't be constructed using the options available in this class, you can use MFEMGenericFESpace instead.
Example Input File Syntax
[FESpaces<<<{"href": "../../../syntax/FESpaces/index.html"}>>>]
[H1FESpace]
type = MFEMScalarFESpace<<<{"description": "Convenience class to construct scalar finite element spaces.", "href": "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
[]
[]
[Variables<<<{"href": "../../../syntax/Variables/index.html"}>>>]
[temperature]
type = MFEMVariable<<<{"description": "Class for adding MFEM variables to the problem (`mfem::ParGridFunction`s).", "href": "../variables/MFEMVariable.html"}>>>
fespace<<<{"description": "The finite element space this variable is defined on."}>>> = H1FESpace
[]
[]
(test/tests/mfem/kernels/heattransfer.i)Input Parameters
- allow_duplicate_execution_on_initialFalseIn the case where this UserObject is depended upon by an initial condition, allow it to be executed twice during the initial setup (once before the IC and again after mesh adaptivity (if applicable).
Default:False
C++ Type:bool
Controllable:No
Description:In the case where this UserObject is depended upon by an initial condition, allow it to be executed twice during the initial setup (once before the IC and again after mesh adaptivity (if applicable).
- execute_onTIMESTEP_ENDThe 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.
Default:TIMESTEP_END
C++ Type:ExecFlagEnum
Options:XFEM_MARK, FORWARD, ADJOINT, HOMOGENEOUS_FORWARD, ADJOINT_TIMESTEP_BEGIN, ADJOINT_TIMESTEP_END, NONE, INITIAL, LINEAR, LINEAR_CONVERGENCE, NONLINEAR, NONLINEAR_CONVERGENCE, POSTCHECK, TIMESTEP_END, TIMESTEP_BEGIN, MULTIAPP_FIXED_POINT_END, MULTIAPP_FIXED_POINT_BEGIN, MULTIAPP_FIXED_POINT_CONVERGENCE, FINAL, CUSTOM, TRANSFER
Controllable:No
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.
- execution_order_group0Execution 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.
Default:0
C++ Type:int
Controllable:No
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.
- force_postauxFalseForces the UserObject to be executed in POSTAUX
Default:False
C++ Type:bool
Controllable:No
Description:Forces the UserObject to be executed in POSTAUX
- force_preauxFalseForces the UserObject to be executed in PREAUX
Default:False
C++ Type:bool
Controllable:No
Description:Forces the UserObject to be executed in PREAUX
- force_preicFalseForces the UserObject to be executed in PREIC during initial setup
Default:False
C++ Type:bool
Controllable:No
Description:Forces the UserObject to be executed in PREIC during initial setup
Execution Scheduling 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:Yes
Description:Set the enabled status of the MooseObject.
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
- use_interpolated_stateFalseFor the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Default:False
C++ Type:bool
Controllable:No
Description:For the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Material Property Retrieval Parameters
Input Files
- (test/tests/mfem/transfers/h1_mfem_sub_mfem_sub/sub_send.i)
- (test/tests/mfem/submeshes/domain_submesh.i)
- (test/tests/mfem/kernels/diffusion.i)
- (test/tests/mfem/multiapps/full_solve_sub.i)
- (test/tests/mfem/submeshes/magnetostatic.i)
- (test/tests/mfem/kernels/darcy.i)
- (test/tests/mfem/transfers/h1_mfem_parent_mfem_sub/parent.i)
- (test/tests/mfem/submeshes/cut_closed_coil.i)
- (test/tests/mfem/transfers/h1_mfem_sub_mfem_sub/parent.i)
- (test/tests/mfem/ics/transient_scalar_ic.i)
- (test/tests/mfem/timesteppers/mfem_multiple_timesequences.i)
- (test/tests/mfem/auxkernels/projection.i)
- (test/tests/mfem/transfers/h1_mfem_sub_mfem_sub/sub_recv.i)
- (test/tests/mfem/multiapps/full_solve_parent.i)
- (test/tests/mfem/multiapps/dt_from_parent_sub.i)
- (test/tests/mfem/kernels/graddiv.i)
- (test/tests/mfem/multiapps/sub_cycling_parent.i)
- (test/tests/mfem/kernels/irrotational.i)
- (test/tests/mfem/submeshes/domain_submesh_transfer.i)
- (test/tests/mfem/timesteppers/mfem_multiple_timesteppers.i)
- (test/tests/mfem/transfers/h1_mfem_parent_mfem_sub/sub.i)
- (test/tests/mfem/ics/scalar_ic.i)
- (test/tests/mfem/submeshes/hphi_magnetostatic.i)
- (test/tests/mfem/multiapps/sub_cycling_sub.i)
- (test/tests/mfem/submeshes/open_coil_source.i)
- (test/tests/mfem/kernels/heattransfer.i)
- (test/tests/mfem/submeshes/hphi_magnetodynamic.i)
- (test/tests/mfem/multiapps/dt_from_parent.i)
(test/tests/mfem/kernels/heattransfer.i)
[Mesh]
type = MFEMMesh
file = ../mesh/mug.e
dim = 3
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[temperature]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[AuxVariables]
inactive = average_temperature
[average_temperature]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[AuxKernels]
inactive = average_field
[average_field]
type = MFEMScalarTimeAverageAux
variable = average_temperature
source = temperature
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = temperature
[]
[dT_dt]
type = MFEMTimeDerivativeMassKernel
variable = temperature
[]
[]
[BCs]
active = 'bottom top_convective'
[bottom]
type = MFEMScalarDirichletBC
variable = temperature
boundary = '1'
coefficient = 1.0
[]
[top_convective]
type = MFEMConvectiveHeatFluxBC
variable = temperature
boundary = '2'
T_infinity = .5
heat_transfer_coefficient = 5
[]
[top_dirichlet]
type = MFEMScalarDirichletBC
variable = temperature
boundary = '2'
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[jacobi]
type = MFEMOperatorJacobiSmoother
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-16
l_max_its = 1000
[]
[Executioner]
type = MFEMTransient
device = cpu
assembly_level = legacy
dt = 2.0
start_time = 0.0
end_time = 6.0
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/HeatTransfer
vtk_format = ASCII
[]
[]
(test/tests/mfem/transfers/h1_mfem_sub_mfem_sub/sub_send.i)
[Mesh]
type = MFEMMesh
file = ../../mesh/square.msh
dim = 3
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[send]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[BCs]
[back]
type = MFEMScalarDirichletBC
variable = send
boundary = 1
coefficient = 1.0
[]
[bottom]
type = MFEMScalarDirichletBC
variable = send
boundary = 2
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = send
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-16
l_max_its = 1000
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/DiffusionSendApp
vtk_format = ASCII
[]
[]
[Executioner]
type = MFEMSteady
[]
(test/tests/mfem/submeshes/domain_submesh.i)
[Mesh]
type = MFEMMesh
file = ../mesh/cylinder-hex-q2.gen
[]
[Problem]
type = MFEMProblem
[]
[SubMeshes]
[wire]
type = MFEMDomainSubMesh
block = interior
[]
[]
[FESpaces]
[SubMeshH1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
submesh = wire
[]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[submesh_potential]
type = MFEMVariable
fespace = SubMeshH1FESpace
[]
[]
[BCs]
[top]
type = MFEMScalarDirichletBC
variable = submesh_potential
boundary = front
coefficient = 1.0
[]
[bottom]
type = MFEMScalarDirichletBC
variable = submesh_potential
boundary = back
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = submesh_potential
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 1000
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/DomainPotential
vtk_format = ASCII
submesh = wire
[]
[]
(test/tests/mfem/kernels/diffusion.i)
[Mesh]
type = MFEMMesh
file = ../mesh/mug.e
dim = 3
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
[]
[]
[Variables]
[concentration]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[AuxVariables]
[concentration_gradient]
type = MFEMVariable
fespace = HCurlFESpace
[]
[]
[AuxKernels]
[grad]
type = MFEMGradAux
variable = concentration_gradient
source = concentration
execute_on = TIMESTEP_END
[]
[]
[BCs]
[bottom]
type = MFEMScalarDirichletBC
variable = concentration
boundary = 'bottom'
coefficient = 1.0
[]
[top]
type = MFEMScalarDirichletBC
variable = concentration
boundary = 'top'
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = concentration
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[jacobi]
type = MFEMOperatorJacobiSmoother
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-16
l_max_its = 1000
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Outputs]
active = ParaViewDataCollection
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/Diffusion
vtk_format = ASCII
[]
[VisItDataCollection]
type = MFEMVisItDataCollection
file_base = OutputData/VisItDataCollection
[]
[ConduitDataCollection]
type = MFEMConduitDataCollection
file_base = OutputData/ConduitDataCollection/Run
protocol = conduit_bin
[]
[]
(test/tests/mfem/multiapps/full_solve_sub.i)
[Problem]
type = MFEMProblem
verbose_multiapps = true
[]
[Mesh]
type = MFEMMesh
file = ../mesh/square.e
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[u]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = u
[]
[td]
type = MFEMTimeDerivativeMassKernel
variable = u
[]
[]
[BCs]
[left]
type = MFEMScalarDirichletBC
variable = u
boundary = left
coefficient = 0
[]
[right]
type = MFEMScalarDirichletBC
variable = u
boundary = right
coefficient = 1
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
print_level = 0
[]
[]
[Solver]
type = MFEMHyprePCG
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 100
print_level = 0
[]
[Executioner]
type = MFEMTransient
num_steps = 2
dt = 0.01
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/full_solve_sub
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
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 = '1 2 3'
[]
[]
[FunctorMaterials]
[Vacuum]
type = MFEMGenericFunctorMaterial
prop_names = reluctivity
prop_values = 1.0
[]
[]
[Kernels]
[curlcurl]
type = MFEMCurlCurlKernel
variable = a_field
coefficient = reluctivity
[]
[source]
type = MFEMMixedVectorGradientKernel
trial_variable = electric_potential
variable = a_field
coefficient = 1.0
block = 1
[]
[]
[Preconditioner]
[ams]
type = MFEMHypreAMS
fespace = HCurlFESpace
singular = true
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = ams
l_tol = 1e-9
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[MultiApps]
[subapp]
type = FullSolveMultiApp
input_files = open_coil_source.i
execute_on = INITIAL
[]
[]
[Transfers]
[from_sub]
type = MultiAppMFEMCopyTransfer
source_variable = electric_potential
variable = electric_potential
from_multi_app = subapp
[]
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/Magnetostatic
vtk_format = ASCII
[]
[]
(test/tests/mfem/kernels/darcy.i)
[Mesh]
type = MFEMMesh
file = ../mesh/star.mesh
uniform_refine = 4
[]
[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
[]
[]
[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 = MFEMSuperLU
[]
[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/transfers/h1_mfem_parent_mfem_sub/parent.i)
[Mesh]
type = MFEMMesh
file = ../../mesh/square.msh
dim = 3
[]
[Problem]
type = MFEMProblem
solve = false
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[AuxVariables]
[u]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[Executioner]
type = MFEMSteady
[]
[MultiApps]
[subapp]
type = FullSolveMultiApp
input_files = sub.i
execute_on = INITIAL
[]
[]
[Transfers]
[from_sub]
type = MultiAppMFEMCopyTransfer
source_variable = u
variable = u
from_multi_app = subapp
[]
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/Diffusion
vtk_format = ASCII
[]
[]
(test/tests/mfem/submeshes/cut_closed_coil.i)
# 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]
type = MFEMProblem
[]
[Mesh]
type = MFEMMesh
file = ../mesh/embedded_concentric_torus.e
[]
[FunctorMaterials]
[Conductor]
type = MFEMGenericFunctorMaterial
prop_names = conductivity
prop_values = ${coil_conductivity}
[]
[]
[ICs]
[coil_external_potential_ic]
type = MFEMScalarBoundaryIC
variable = coil_external_potential
boundary = ${coil_cut_surface}
coefficient = ${coil_loop_voltage}
[]
[]
[SubMeshes]
[cut]
type = MFEMCutTransitionSubMesh
cut_boundary = ${coil_cut_surface}
block = ${initial_coil_domains}
transition_subdomain = transition_dom
transition_subdomain_boundary = transition_bdr
closed_subdomain = coil_dom
[]
[coil]
type = MFEMDomainSubMesh
block = coil_dom
execution_order_group = 2
[]
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
[]
[CoilH1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
submesh = coil
[]
[CoilHCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
submesh = coil
[]
[TransitionH1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
submesh = cut
[]
[TransitionHCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
submesh = cut
[]
[]
[Variables]
[coil_induced_potential]
type = MFEMVariable
fespace = CoilH1FESpace
[]
[]
[AuxVariables]
[coil_external_potential]
type = MFEMVariable
fespace = CoilH1FESpace
[]
[transition_external_potential]
type = MFEMVariable
fespace = TransitionH1FESpace
[]
[transition_external_e_field]
type = MFEMVariable
fespace = TransitionHCurlFESpace
[]
[induced_potential]
type = MFEMVariable
fespace = H1FESpace
[]
[induced_e_field]
type = MFEMVariable
fespace = HCurlFESpace
[]
[external_e_field]
type = MFEMVariable
fespace = HCurlFESpace
[]
[e_field]
type = MFEMVariable
fespace = HCurlFESpace
[]
[]
[AuxKernels]
[update_induced_e_field]
type = MFEMGradAux
variable = induced_e_field
source = induced_potential
scale_factor = -1.0
execute_on = TIMESTEP_END
[]
[update_external_e_field]
type = MFEMGradAux
variable = transition_external_e_field
source = transition_external_potential
scale_factor = -1.0
execute_on = TIMESTEP_END
[]
[update_total_e_field]
type = MFEMSumAux
variable = e_field
first_source_variable = induced_e_field
second_source_variable = external_e_field
execute_on = TIMESTEP_END
execution_order_group = 3
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = coil_induced_potential
coefficient = conductivity
[]
[source]
type = MFEMMixedGradGradKernel
trial_variable = coil_external_potential
variable = coil_induced_potential
coefficient = conductivity
block = 'transition_dom'
[]
[]
[Solver]
type = MFEMSuperLU
[]
[Executioner]
type = MFEMSteady
[]
[Transfers]
[submesh_transfer_from_coil]
type = MFEMSubMeshTransfer
from_variable = coil_induced_potential
to_variable = induced_potential
execute_on = TIMESTEP_END
[]
[submesh_transfer_to_transition]
type = MFEMSubMeshTransfer
from_variable = coil_external_potential
to_variable = transition_external_potential
execute_on = TIMESTEP_END
[]
[submesh_transfer_from_transition]
type = MFEMSubMeshTransfer
from_variable = transition_external_e_field
to_variable = external_e_field
execute_on = TIMESTEP_END
execution_order_group = 2
[]
[]
[Outputs]
[GlobalParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/WholePotentialCoil
vtk_format = ASCII
[]
[TransitionParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/CutPotentialCoil
vtk_format = ASCII
submesh = cut
[]
[CoilParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/Coil
vtk_format = ASCII
submesh = coil
[]
[]
(test/tests/mfem/transfers/h1_mfem_sub_mfem_sub/parent.i)
[Mesh]
type = MFEMMesh
file = ../../mesh/square.msh
dim = 3
[]
[Problem]
type = MFEMProblem
solve = false
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
#[AuxVariables]
# [recv]
# type = MFEMVariable
# fespace = H1FESpace
# []
#[]
[MultiApps]
[recv_app]
type = FullSolveMultiApp
input_files = sub_recv.i
execute_on = FINAL
[]
[send_app]
type = FullSolveMultiApp
input_files = sub_send.i
execute_on = INITIAL
[]
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Transfers]
[to_sub]
type = MultiAppMFEMCopyTransfer
source_variable = send
variable = recv
from_multi_app = send_app
to_multi_app = recv_app
[]
[]
(test/tests/mfem/ics/transient_scalar_ic.i)
[Mesh]
type = MFEMMesh
file = ../mesh/cylinder-hex-q2.gen
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[L2FESpace]
type = MFEMScalarFESpace
fec_type = L2
fec_order = CONSTANT
[]
[]
[Variables]
[h1_scalar]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[AuxVariables]
[l2_scalar]
type = MFEMVariable
fespace = L2FESpace
[]
[]
[Functions]
[height]
type = ParsedFunction
expression = 'z'
[]
[]
[ICs]
[l2_scalar_ic]
type = MFEMScalarIC
variable = l2_scalar
coefficient = 2.0
[]
[h1_scalar_ic]
type = MFEMScalarIC
variable = h1_scalar
coefficient = height
[]
[]
[Kernels]
[h1_laplacian]
type = MFEMDiffusionKernel
variable = h1_scalar
coefficient = 1.0
[]
[dh1_dt]
type = MFEMTimeDerivativeMassKernel
variable = h1_scalar
coefficient = 1.0
[]
[]
[BCs]
[bottom]
type = MFEMScalarDirichletBC
variable = h1_scalar
boundary = '1'
coefficient = height
[]
[top_dirichlet]
type = MFEMScalarDirichletBC
variable = h1_scalar
boundary = '2'
coefficient = height
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 100
[]
[Executioner]
type = MFEMTransient
device = cpu
dt = 2.0
start_time = 0.0
end_time = 2.0
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/TransientScalarIC
vtk_format = ASCII
[]
[]
(test/tests/mfem/timesteppers/mfem_multiple_timesequences.i)
[Problem]
type = MFEMProblem
[]
[Mesh]
type = MFEMMesh
file = ../mesh/square.e
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[u]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = u
coefficient = 0.1
[]
[time]
type = MFEMTimeDerivativeMassKernel
variable = u
[]
[]
[BCs]
[left]
type = MFEMScalarDirichletBC
variable = u
boundary = left
coefficient = 0
[]
[right]
type = MFEMScalarDirichletBC
variable = u
boundary = right
coefficient = 1
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHyprePCG
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 100
[]
[Executioner]
type = MFEMTransient
device = cpu
end_time = 0.8
# Use as many different time sequence steppers as we could to test the compositionDT
[TimeSteppers]
[ConstDT1]
type = ConstantDT
dt = 0.2
[]
[ConstDT2]
type = ConstantDT
dt = 0.1
[]
[LogConstDT]
type = LogConstantDT
log_dt = 0.2
first_dt = 0.1
[]
[Timesequence1]
type = TimeSequenceStepper
time_sequence = '0 0.25 0.3 0.5 0.8'
[]
[Timesequence2]
type = CSVTimeSequenceStepper
file_name = timesequence.csv
column_name = time
[]
[Timesequence3]
type = ExodusTimeSequenceStepper
mesh = timesequence.e
[]
[]
[]
[Postprocessors]
[timestep]
type = TimePostprocessor
execute_on = 'timestep_end'
[]
[]
[Outputs]
csv = true
file_base='mfem_multiple_timesequences'
[]
(test/tests/mfem/auxkernels/projection.i)
[Mesh]
type = MFEMMesh
file = ../mesh/hinomaru.e
dim = 2
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
[]
[L2FESpace]
type = MFEMScalarFESpace
fec_type = L2
fec_order = CONSTANT
[]
[]
[Variables]
[Az]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[AuxVariables]
[J]
type = MFEMVariable
fespace = L2FESpace
[]
[GAz]
type = MFEMVariable
fespace = HCurlFESpace
[]
[GAz(copy)]
type = MFEMVariable
fespace = HCurlFESpace
[]
[]
[Kernels]
[diffusion]
type = MFEMDiffusionKernel
variable = Az
[]
[source]
type = MFEMDomainLFKernel
variable = Az
coefficient = J_source
[]
[]
[AuxKernels]
[J]
type = MFEMScalarProjectionAux
variable = J
coefficient = J_source
[]
[GAz]
type = MFEMGradAux
variable = GAz
source = Az
[]
[GAz(copy)]
type = MFEMVectorProjectionAux
variable = GAz(copy)
vector_coefficient = GAz
[]
[]
[BCs]
[essential]
type = MFEMScalarDirichletBC
variable = Az
boundary = 1
coefficient = 1
[]
[]
[FunctorMaterials]
[J_wire]
type = MFEMGenericFunctorMaterial
prop_names = J_source
prop_values = 8.0
block = wire
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHyprePCG
preconditioner = boomeramg
l_tol = 1e-16
[]
[Executioner]
type = MFEMSteady
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/Projection
vtk_format = ASCII
[]
[]
(test/tests/mfem/transfers/h1_mfem_sub_mfem_sub/sub_recv.i)
[Mesh]
type = MFEMMesh
file = ../../mesh/square.msh
dim = 3
[]
[Problem]
type = MFEMProblem
solve = false
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[AuxVariables]
[recv]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[Executioner]
type = MFEMSteady
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/DiffusionRecvApp
vtk_format = ASCII
[]
[]
(test/tests/mfem/multiapps/full_solve_parent.i)
[Problem]
type = MFEMProblem
verbose_multiapps = true
[]
[Mesh]
type = MFEMMesh
file = ../mesh/square.e
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[u]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = u
[]
[td]
type = MFEMTimeDerivativeMassKernel
variable = u
[]
[]
[BCs]
[left]
type = MFEMScalarDirichletBC
variable = u
boundary = left
coefficient = 0
[]
[right]
type = MFEMScalarDirichletBC
variable = u
boundary = right
coefficient = 1
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
print_level = 0
[]
[]
[Solver]
type = MFEMHyprePCG
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 100
print_level = 0
[]
[Executioner]
type = MFEMTransient
device = cpu
num_steps = 2
dt = 0.1
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/full_solve_parent
vtk_format = ASCII
[]
[]
[MultiApps]
[full_solve]
type = FullSolveMultiApp
app_type = MooseTestApp
execute_on = timestep_begin
input_files = full_solve_sub.i
keep_full_output_history = true
[]
[]
(test/tests/mfem/multiapps/dt_from_parent_sub.i)
[Problem]
type = MFEMProblem
[]
[Mesh]
type = MFEMMesh
file = ../mesh/square.e
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[u]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = u
[]
[td]
type = MFEMTimeDerivativeMassKernel
variable = u
[]
[]
[BCs]
[left]
type = MFEMScalarDirichletBC
variable = u
boundary = left
coefficient = 0
[]
[right]
type = MFEMScalarDirichletBC
variable = u
boundary = right
coefficient = 1
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHyprePCG
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 100
[]
[Executioner]
type = MFEMTransient
num_steps = 10
dt = 1 # This will be constrained by the parent solve
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/dt_from_parent_sub
vtk_format = ASCII
[]
[]
(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
[]
[]
[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/multiapps/sub_cycling_parent.i)
[Problem]
type = MFEMProblem
verbose_multiapps = true
[]
[Mesh]
type = MFEMMesh
file = ../mesh/square.e
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[u]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = u
[]
[td]
type = MFEMTimeDerivativeMassKernel
variable = u
[]
[]
[BCs]
[left]
type = MFEMScalarDirichletBC
variable = u
boundary = left
coefficient = 0
[]
[right]
type = MFEMScalarDirichletBC
variable = u
boundary = right
coefficient = 1
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHyprePCG
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 100
[]
[Executioner]
type = MFEMTransient
device = cpu
num_steps = 2
dt = 0.1
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/sub_cycling_parent
vtk_format = ASCII
[]
[]
[MultiApps]
[sub]
type = TransientMultiApp
app_type = MooseTestApp
execute_on = timestep_end
positions = '0 0 0'
input_files = sub_cycling_sub.i
sub_cycling = true
[]
[]
(test/tests/mfem/kernels/irrotational.i)
# 2D irrotational vortex with Nedelec elements of the first kind.
centre_x = -0.75
centre_y = 0.1
[Mesh]
type = MFEMMesh
file = ../mesh/vortex.msh
dim = 2
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = SEVENTH
[]
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = SEVENTH
[]
[]
[Variables]
[velocity_potential]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[AuxVariables]
[velocity]
type = MFEMVariable
fespace = HCurlFESpace
[]
[]
[Functions]
[speed]
type = ParsedFunction
expression = '1 / sqrt((x-x0)^2 + (y-y0)^2)'
symbol_names = 'x0 y0'
symbol_values = '${centre_x} ${centre_y}'
[]
[theta]
type = ParsedFunction
expression = 'atan2(y-y0, x-x0)'
symbol_names = 'x0 y0'
symbol_values = '${centre_x} ${centre_y}'
[]
[exact_velocity]
type = ParsedVectorFunction
expression_x = '-v * sin(th)'
expression_y = 'v * cos(th)'
symbol_names = 'v th'
symbol_values = 'speed theta'
[]
[]
[BCs]
[potential_velocity_boundary]
type = MFEMScalarDirichletBC
variable = velocity_potential
boundary = '1'
coefficient = theta
[]
[]
[Kernels]
[laplacian]
type = MFEMDiffusionKernel
variable = velocity_potential
[]
[]
[AuxKernels]
[grad]
type = MFEMGradAux
variable = velocity
source = velocity_potential
execute_on = TIMESTEP_END
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-16
l_max_its = 1000
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Postprocessors]
[potential_error]
type = MFEML2Error
variable = velocity_potential
function = theta
[]
[velocity_error]
type = MFEMVectorL2Error
variable = velocity
function = exact_velocity
execution_order_group = 1
[]
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/Irrotational
vtk_format = ASCII
[]
[L2CSV]
type = CSV
file_base = OutputData/Irrotational
[]
[]
(test/tests/mfem/submeshes/domain_submesh_transfer.i)
[Mesh]
type = MFEMMesh
file = ../mesh/cylinder-hex-q2.gen
[]
[Problem]
type = MFEMProblem
[]
[SubMeshes]
[wire]
type = MFEMDomainSubMesh
block = interior
[]
[]
[FESpaces]
[SubMeshH1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
submesh = wire
[]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[submesh_potential]
type = MFEMVariable
fespace = SubMeshH1FESpace
[]
[]
[AuxVariables]
[potential]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[BCs]
[top]
type = MFEMScalarDirichletBC
variable = submesh_potential
boundary = front
coefficient = 1.0
[]
[bottom]
type = MFEMScalarDirichletBC
variable = submesh_potential
boundary = back
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = submesh_potential
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 1000
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Transfers]
[submesh_transfer]
type = MFEMSubMeshTransfer
from_variable = submesh_potential
to_variable = potential
[]
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/DomainPotentialTransfer
vtk_format = ASCII
[]
[]
(test/tests/mfem/timesteppers/mfem_multiple_timesteppers.i)
[Problem]
type = MFEMProblem
[]
[Mesh]
type = MFEMMesh
file = ../mesh/square.e
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[u]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = u
coefficient = 0.1
[]
[time]
type = MFEMTimeDerivativeMassKernel
variable = u
[]
[]
[Functions]
[dts]
type = PiecewiseLinear
x = '0 0.85 2'
y = '0.2 0.15 0.2'
[]
[]
[BCs]
[left]
type = MFEMScalarDirichletBC
variable = u
boundary = left
coefficient = 0
[]
[right]
type = MFEMScalarDirichletBC
variable = u
boundary = right
coefficient = 1
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHyprePCG
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 100
[]
[Executioner]
type = MFEMTransient
device = cpu
end_time = 0.8
# Use as many different time steppers as we could to test the compositionDT,
# SolutionTimeAdaptiveDT give slightly different dt per run, set rel_err = 1e-2
# to ensure the test won't fail due to the small difference in the high-digit.
[TimeSteppers]
[ConstDT1]
type = ConstantDT
dt = 0.2
[]
[FunctionDT]
type = FunctionDT
function = dts
[]
[LogConstDT]
type = LogConstantDT
log_dt = 0.2
first_dt = 0.1
[]
[IterationAdapDT]
type = IterationAdaptiveDT
dt = 0.5
[]
[Timesequence]
type = TimeSequenceStepper
time_sequence = '0 0.25 0.3 0.5 0.8'
[]
[]
[]
[Postprocessors]
[timestep]
type = TimePostprocessor
execute_on = 'timestep_end'
[]
[]
[Outputs]
csv = true
file_base='mfem_multiple_timesteppers'
[]
(test/tests/mfem/transfers/h1_mfem_parent_mfem_sub/sub.i)
[Mesh]
type = MFEMMesh
file = ../../mesh/square.msh
dim = 3
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[u]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[BCs]
[bottom]
type = MFEMScalarDirichletBC
variable = u
boundary = 2
coefficient = 1.0
[]
[top]
type = MFEMScalarDirichletBC
variable = u
boundary = 4
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = u
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-16
l_max_its = 1000
[]
[Executioner]
type = MFEMSteady
[]
[MultiApps]
active = ''
[subapp]
type = FullSolveMultiApp
input_files = parent.i
execute_on = FINAL
[]
[]
[Transfers]
active = ''
[to_sub]
type = MultiAppMFEMCopyTransfer
source_variable = u
variable = u
to_multi_app = subapp
[]
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/DiffusionSub
vtk_format = ASCII
[]
[]
(test/tests/mfem/ics/scalar_ic.i)
[Mesh]
type = MFEMMesh
file = ../mesh/cylinder-hex-q2.gen
[]
[Problem]
type = MFEMProblem
solve = false
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[L2FESpace]
type = MFEMScalarFESpace
fec_type = L2
fec_order = CONSTANT
[]
[]
[Variables]
[h1_scalar]
type = MFEMVariable
fespace = H1FESpace
[]
[l2_scalar]
type = MFEMVariable
fespace = L2FESpace
[]
[]
[Functions]
[height]
type = ParsedFunction
expression = 'z'
[]
[]
[ICs]
[l2_scalar_ic]
type = MFEMScalarIC
variable = l2_scalar
coefficient = 2.0
[]
[h1_scalar_ic]
type = MFEMScalarIC
variable = h1_scalar
coefficient = height
[]
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/ScalarIC
vtk_format = ASCII
[]
[]
(test/tests/mfem/submeshes/hphi_magnetostatic.i)
# 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]
type = MFEMProblem
[]
[Mesh]
type = MFEMMesh
file = ../mesh/split_embedded_concentric_torus.e
[]
[FunctorMaterials]
[Conductor]
type = MFEMGenericFunctorMaterial
prop_names = permeability
prop_values = ${vacuum_permeability}
[]
[]
[ICs]
[vacuum_cut_potential_ic]
type = MFEMScalarBoundaryIC
variable = vacuum_cut_potential
boundary = ${vacuum_cut_surface}
coefficient = ${conductor_current}
[]
[]
[SubMeshes]
[cut]
type = MFEMCutTransitionSubMesh
cut_boundary = ${vacuum_cut_surface}
block = ${initial_vacuum_domains}
transition_subdomain = transition_dom
transition_subdomain_boundary = transition_bdr
closed_subdomain = vacuum_dom
[]
[vacuum]
type = MFEMDomainSubMesh
block = vacuum_dom
execution_order_group = 2
[]
[]
[FESpaces]
[VacuumH1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
submesh = vacuum
[]
[VacuumHCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
submesh = vacuum
[]
[TransitionH1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
submesh = cut
[]
[TransitionHCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
submesh = cut
[]
[HCurlFESpace]
type = MFEMVectorFESpace
fec_type = ND
fec_order = FIRST
[]
[]
[Variables]
[vacuum_magnetic_potential]
type = MFEMVariable
fespace = VacuumH1FESpace
[]
[]
[AuxVariables]
[vacuum_cut_potential]
type = MFEMVariable
fespace = VacuumH1FESpace
[]
[transition_cut_potential]
type = MFEMVariable
fespace = TransitionH1FESpace
[]
[transition_cut_function_field]
type = MFEMVariable
fespace = TransitionHCurlFESpace
[]
[background_h_field]
type = MFEMVariable
fespace = VacuumHCurlFESpace
[]
[cut_function_field]
type = MFEMVariable
fespace = VacuumHCurlFESpace
[]
[vacuum_h_field]
type = MFEMVariable
fespace = VacuumHCurlFESpace
[]
[h_field]
type = MFEMVariable
fespace = HCurlFESpace
[]
[]
[AuxKernels]
[update_background_h_field]
type = MFEMGradAux
variable = background_h_field
source = vacuum_magnetic_potential
scale_factor = -1.0
execute_on = TIMESTEP_END
[]
[update_transition_cut_function_field]
type = MFEMGradAux
variable = transition_cut_function_field
source = transition_cut_potential
scale_factor = -1.0
execute_on = TIMESTEP_END
[]
[update_total_h_field]
type = MFEMSumAux
variable = vacuum_h_field
first_source_variable = background_h_field
second_source_variable = cut_function_field
execute_on = TIMESTEP_END
execution_order_group = 3
[]
[]
[BCs]
# Set zero of magnetic potential on symmetry plane
[Exterior]
type = MFEMScalarDirichletBC
variable = vacuum_magnetic_potential
boundary = 'Cut'
coefficient = 0.0
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = vacuum_magnetic_potential
coefficient = permeability
[]
[source]
type = MFEMMixedGradGradKernel
trial_variable = vacuum_cut_potential
variable = vacuum_magnetic_potential
coefficient = permeability
block = 'transition_dom'
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 100
[]
[Executioner]
type = MFEMSteady
[]
[Transfers]
[submesh_transfer_to_transition]
type = MFEMSubMeshTransfer
from_variable = vacuum_cut_potential
to_variable = transition_cut_potential
execute_on = TIMESTEP_END
[]
[submesh_transfer_from_transition]
type = MFEMSubMeshTransfer
from_variable = transition_cut_function_field
to_variable = cut_function_field
execute_on = TIMESTEP_END
execution_order_group = 2
[]
[submesh_transfer_from_vacuum]
type = MFEMSubMeshTransfer
from_variable = vacuum_h_field
to_variable = h_field
execute_on = TIMESTEP_END
execution_order_group = 4
[]
[]
[Postprocessors]
[MagneticEnergy]
type = MFEMVectorFEInnerProductIntegralPostprocessor
coefficient = ${fparse 0.5*vacuum_permeability}
dual_variable = vacuum_h_field
primal_variable = vacuum_h_field
execution_order_group = 4
block = 'Exterior'
[]
[]
[Outputs]
[ReportedPostprocessors]
type = CSV
file_base = OutputData/HPhiMagnetostaticClosedCoilCSV
[]
[VacuumParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/HPhiMagnetostaticClosedCoil
vtk_format = ASCII
submesh = vacuum
[]
[]
(test/tests/mfem/multiapps/sub_cycling_sub.i)
[Problem]
type = MFEMProblem
[]
[Mesh]
type = MFEMMesh
file = ../mesh/square.e
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[u]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = u
[]
[td]
type = MFEMTimeDerivativeMassKernel
variable = u
[]
[]
[BCs]
[left]
type = MFEMScalarDirichletBC
variable = u
boundary = left
coefficient = 0
[]
[right]
type = MFEMScalarDirichletBC
variable = u
boundary = right
coefficient = 1
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHyprePCG
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 100
[]
[Executioner]
type = MFEMTransient
num_steps = 2
dt = 0.01
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/sub_cycling_sub
vtk_format = ASCII
[]
[]
(test/tests/mfem/submeshes/open_coil_source.i)
[Mesh]
type = MFEMMesh
file = ../mesh/cylinder-hex-q2.gen
[]
[Problem]
type = MFEMProblem
[]
[SubMeshes]
[wire]
type = MFEMDomainSubMesh
block = 1
[]
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[SubMeshH1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
submesh = wire
[]
[]
[Variables]
[electric_potential]
type = MFEMVariable
fespace = H1FESpace
[]
[submesh_potential]
type = MFEMVariable
fespace = SubMeshH1FESpace
[]
[]
[BCs]
[high_terminal]
type = MFEMScalarDirichletBC
variable = submesh_potential
boundary = '1'
coefficient = 1.0
[]
[low_terminal]
type = MFEMScalarDirichletBC
variable = submesh_potential
boundary = '2'
coefficient = 0.0
[]
[]
[FunctorMaterials]
[Substance]
type = MFEMGenericFunctorMaterial
prop_names = conductivity
prop_values = 1.0
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = submesh_potential
coefficient = conductivity
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-10
l_max_its = 1000
[]
[Executioner]
type = MFEMSteady
device = cpu
[]
[Transfers]
[submesh_potential_transfer]
type = MFEMSubMeshTransfer
from_variable = submesh_potential
to_variable = electric_potential
[]
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/OpenCoilSourceSubMesh
vtk_format = ASCII
submesh = wire
[]
[]
(test/tests/mfem/kernels/heattransfer.i)
[Mesh]
type = MFEMMesh
file = ../mesh/mug.e
dim = 3
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[temperature]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[AuxVariables]
inactive = average_temperature
[average_temperature]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[AuxKernels]
inactive = average_field
[average_field]
type = MFEMScalarTimeAverageAux
variable = average_temperature
source = temperature
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = temperature
[]
[dT_dt]
type = MFEMTimeDerivativeMassKernel
variable = temperature
[]
[]
[BCs]
active = 'bottom top_convective'
[bottom]
type = MFEMScalarDirichletBC
variable = temperature
boundary = '1'
coefficient = 1.0
[]
[top_convective]
type = MFEMConvectiveHeatFluxBC
variable = temperature
boundary = '2'
T_infinity = .5
heat_transfer_coefficient = 5
[]
[top_dirichlet]
type = MFEMScalarDirichletBC
variable = temperature
boundary = '2'
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[jacobi]
type = MFEMOperatorJacobiSmoother
[]
[]
[Solver]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-16
l_max_its = 1000
[]
[Executioner]
type = MFEMTransient
device = cpu
assembly_level = legacy
dt = 2.0
start_time = 0.0
end_time = 6.0
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/HeatTransfer
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
execution_order_group = 4
[]
[]
[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_variable = h_field
variable = 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
execution_order_group = 5
block = 'TorusCore TorusSheath'
[]
[]
[Outputs]
[ReportedPostprocessors]
type = CSV
file_base = OutputData/HPhiMagnetodynamicClosedCoilCSV
[]
[VacuumParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/HPhiMagnetodynamicClosedCoil
vtk_format = ASCII
[]
[]
(test/tests/mfem/multiapps/dt_from_parent.i)
[Problem]
type = MFEMProblem
[]
[Mesh]
type = MFEMMesh
file = ../mesh/square.e
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[u]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[Kernels]
[diff]
type = MFEMDiffusionKernel
variable = u
[]
[td]
type = MFEMTimeDerivativeMassKernel
variable = u
[]
[]
[BCs]
[left]
type = MFEMScalarDirichletBC
variable = u
boundary = left
coefficient = 0
[]
[right]
type = MFEMScalarDirichletBC
variable = u
boundary = right
coefficient = 1
[]
[]
[Preconditioner]
[boomeramg]
type = MFEMHypreBoomerAMG
[]
[]
[Solver]
type = MFEMHyprePCG
preconditioner = boomeramg
l_tol = 1e-8
l_max_its = 100
[]
[Executioner]
type = MFEMTransient
device = cpu
num_steps = 10
dt = 0.2
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/dt_from_parent
vtk_format = ASCII
[]
[]
[MultiApps]
[sub_app]
type = TransientMultiApp
app_type = MooseTestApp
input_files = 'dt_from_parent_sub.i'
[]
[]