- blockThe list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains.
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of subdomains (names or ids) that this object will be restricted to. Leave empty to apply to all subdomains.
- coefficient1.Name of property for the mass coefficient k. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number.
Default:1.
C++ Type:MFEMScalarCoefficientName
Controllable:No
Description:Name of property for the mass coefficient k. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number.
- variableVariable labelling the weak form this kernel is added to
C++ Type:VariableName
Unit:(no unit assumed)
Controllable:No
Description:Variable labelling the weak form this kernel is added to
MFEMTimeDerivativeMassKernel
Overview
Adds the domain integrator for integrating the bilinear form
where and is a scalar coefficient.
This term arises from the weak form of the operator
Example Input File Syntax
[Kernels<<<{"href": "../../../syntax/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": "MFEMDiffusionKernel.html"}>>>
variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = temperature
[]
[dT_dt]
type = MFEMTimeDerivativeMassKernel<<<{"description": "Adds the domain integrator to an MFEM problem for the bilinear form $(k \\dot{u}, v)_\\Omega$ arising from the weak form of the operator $k \\dot{u}$.", "href": "MFEMTimeDerivativeMassKernel.html"}>>>
variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = temperature
[]
[](test/tests/mfem/kernels/heattransfer.i)Input Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Controllable:No
Description:Set the enabled status of the MooseObject.
Advanced Parameters
Input Files
- (test/tests/mfem/multiapps/full_solve_parent.i)
- (test/tests/mfem/multiapps/dt_from_parent.i)
- (test/tests/mfem/timesteppers/mfem_multiple_timesequences.i)
- (test/tests/mfem/multiapps/dt_from_parent_sub.i)
- (test/tests/mfem/multiapps/full_solve_sub.i)
- (test/tests/mfem/kernels/heattransfer.i)
- (test/tests/mfem/kernels/nl_heatconduction.i)
- (test/tests/mfem/kernels/nl_heattransfer.i)
- (test/tests/mfem/ics/transient_scalar_ic.i)
- (test/tests/mfem/multiapps/sub_cycling_parent.i)
- (test/tests/mfem/timesteppers/mfem_multiple_timesteppers.i)
- (test/tests/mfem/kernels/mixed_heattransfer.i)
- (test/tests/mfem/multiapps/sub_cycling_sub.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/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.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'
[]
[]
(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/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/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/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/kernels/nl_heatconduction.i)
# Implementation of MFEM Example 16, for a time dependent nonlinear heat equation problem of the
# form
# dT/dt = \nabla \cdot (\kappa + \alpha T) \nabla T
kappa = 0.5
alpha = 1e-2
[Mesh]
type = MFEMMesh
file = ../mesh/star.mesh
uniform_refine = 1
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = SECOND
[]
[]
[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
[]
[]
[Functions]
[initial]
type = ParsedFunction
expression = 'if((x*x + y*y > 0.251), 1.0, 2.0)'
[]
[diffusivity_temperature_dependence]
type = MFEMParsedFunction
expression = 'alpha * temperature'
symbol_names = 'alpha temperature'
symbol_values = '${alpha} temperature'
[]
[]
[ICs]
[diffused_ic]
type = MFEMScalarIC
coefficient = initial
variable = temperature
[]
[]
[Kernels]
[nl_diffusion]
type = MFEMNLDiffusionKernel
variable = temperature
k_coefficient = diffusivity_temperature_dependence
dk_du_coefficient = ${alpha}
[]
[linear_diffusion]
type = MFEMDiffusionKernel
variable = temperature
coefficient = ${kappa}
[]
[dT_dt]
type = MFEMTimeDerivativeMassKernel
variable = temperature
[]
[]
[Solver]
type = MFEMMUMPS
print_level = 0
[]
[Executioner]
type = MFEMTransient
device = cpu
assembly_level = legacy
dt = 1e-2
start_time = 0.0
end_time = 0.5
nl_max_its = 30
nl_abs_tol = 1.0e-5
nl_rel_tol = 1.0e-5
print_level = 1
[]
[VectorPostprocessors]
[centre_temperature]
type = MFEMPointValueSampler
variable = 'temperature'
points = '0.0 0.0 0.0'
execute_on = TIMESTEP_END
[]
[]
[Outputs]
file_base = OutputData/NLHeatConduction
csv = true
time_step_interval = 10
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/NLHeatConduction
vtk_format = ASCII
[]
[]
(test/tests/mfem/kernels/nl_heattransfer.i)
[Problem]
type = MFEMProblem
[]
[Mesh]
type = MFEMMesh
file = ../mesh/stacked_hexes.e
[]
[FESpaces]
[H1FESpace]
type = MFEMScalarFESpace
fec_type = H1
fec_order = FIRST
[]
[]
[Variables]
[temperature]
type = MFEMVariable
fespace = H1FESpace
[]
[]
[ICs]
[temperature_ic]
type = MFEMScalarIC
coefficient = 200.0
variable = temperature
[]
[]
[Functions]
[T_inf]
type = MFEMParsedFunction
expression = 'temperature + 1'
symbol_names = 'temperature'
symbol_values = 'temperature'
[]
[htc]
type = MFEMParsedFunction
expression = 'temperature/100 + 1'
symbol_names = 'temperature'
symbol_values = 'temperature'
[]
[dhtc_dT]
type = MFEMParsedFunction
expression = '1 / 100'
symbol_names = 'temperature'
symbol_values = 'temperature'
[]
[]
[Kernels]
[dT_dt]
type = MFEMTimeDerivativeMassKernel
variable = temperature
[]
[diffusion]
type = MFEMDiffusionKernel
variable = temperature
[]
[]
[BCs]
active = nonlinear
[nonlinear]
type = MFEMNLConvectiveHeatFluxBC
variable = temperature
boundary = 'right'
T_infinity = T_inf
heat_transfer_coefficient = htc
d_heat_transfer_dT_coefficient = dhtc_dT
[]
[linearized]
type = MFEMNLConvectiveHeatFluxBC
variable = temperature
boundary = 'right'
T_infinity = 201.0
heat_transfer_coefficient = 3.0
d_heat_transfer_dT_coefficient = 0.0
[]
[]
[VectorPostprocessors]
[line_sample]
type = MFEMLineValueSampler
variable = 'temperature'
start_point = '0.0 0.5 0.5'
end_point = '1.0 0.5 0.5'
num_points = 3
execute_on = TIMESTEP_END
[]
[]
[Solver]
type = MFEMMUMPS
print_level = 0
[]
[Executioner]
type = MFEMTransient
device = cpu
assembly_level = legacy
dt = 1
num_steps = 3
nl_max_its = 150
nl_abs_tol = 1e-12
print_level = 1
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/NLHeatTransfer
vtk_format = ASCII
[]
[CSV]
type = CSV
file_base = OutputData/NLHeatTransfer
time_step_interval = 3
[]
[]
(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
basis = GaussLegendre
[]
[]
[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/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/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/kernels/mixed_heattransfer.i)
# Mixed heat transfer problem.
# Based on Firedrake Irksome demo_mixed_heat example:
# https://www.firedrakeproject.org/Irksome/demos/demo_mixed_heat.py.html
[Mesh]
type = MFEMMesh
file = ../mesh/square.e
[]
[Problem]
type = MFEMProblem
[]
[FESpaces]
[HDivFESpace]
type = MFEMVectorFESpace
fec_type = RT
fec_order = FIRST
[]
[L2FESpace]
type = MFEMScalarFESpace
fec_type = L2
fec_order = FIRST
[]
[]
[Variables]
[time_integrated_heat_flux]
type = MFEMVariable
fespace = HDivFESpace
time_derivative = heat_flux
[]
[temperature]
type = MFEMVariable
fespace = L2FESpace
[]
[]
[Kernels]
[dT_dt,T']
type = MFEMTimeDerivativeMassKernel
variable = temperature
[]
[divh,T']
type = MFEMVectorFEDivergenceKernel
trial_variable = heat_flux
variable = temperature
[]
[h,h']
type = MFEMTimeDerivativeVectorFEMassKernel
variable = time_integrated_heat_flux
[]
[-T,div.h']
type = MFEMVectorFEDivergenceKernel
trial_variable = temperature
variable = time_integrated_heat_flux
coefficient = -1.0
transpose = true
[]
[]
[BCs]
[gamma_T_right]
type = MFEMVectorFEBoundaryFluxIntegratedBC
variable = time_integrated_heat_flux
coefficient = 0.0
boundary = 2
[]
[gamma_T_left]
type = MFEMVectorFEBoundaryFluxIntegratedBC
variable = time_integrated_heat_flux
coefficient = -1.0
boundary = 4
[]
[gamma_h_topbottom]
type = MFEMVectorNormalDirichletBC
variable = time_integrated_heat_flux
vector_coefficient = '0.0 0.0'
boundary = '1 3'
[]
[]
[Solver]
type = MFEMSuperLU
[]
[Executioner]
type = MFEMTransient
device = cpu
assembly_level = legacy
dt = 0.03
start_time = 0.0
end_time = 0.09
[]
[Outputs]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/MixedHeatTransfer
vtk_format = ASCII
[]
[]
(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
[]
[]