Transient Heat Transfer
Summary
Solves a transient heat conduction problem with a boundary parameterized by a heat transfer coefficient that exchanges heat with a thermal reservoir.
Description
This problem solves the transient heat equation with strong form:
where , , , and , subject to the initial condition .
In this example, we solve this using the weak form
where
Example File
[Mesh<<<{"href": "../Mesh/index.html"}>>>]
type = MFEMMesh
file = ../mesh/mug.e
dim = 3
[]
[Problem<<<{"href": "../Problem/index.html"}>>>]
type = MFEMProblem
[]
[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
[]
[]
[Variables<<<{"href": "../Variables/index.html"}>>>]
[temperature]
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
[]
[]
[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"}>>> = 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": "../../source/mfem/kernels/MFEMTimeDerivativeMassKernel.html"}>>>
variable<<<{"description": "Variable labelling the weak form this kernel is added to"}>>> = temperature
[]
[]
[BCs<<<{"href": "../BCs/index.html"}>>>]
active<<<{"description": "If specified only the blocks named will be visited and made active"}>>> = 'bottom top_convective'
[bottom]
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"}>>> = temperature
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies. Defaults to all boundaries."}>>> = '1'
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."}>>> = 1.0
[]
[top_convective]
type = MFEMConvectiveHeatFluxBC<<<{"description": "Convective heat transfer boundary condition with temperature and heat transfer coefficent given by material properties to add to MFEM problems.", "href": "../../source/mfem/bcs/MFEMConvectiveHeatFluxBC.html"}>>>
variable<<<{"description": "Variable on which to apply the boundary condition"}>>> = temperature
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies. Defaults to all boundaries."}>>> = '2'
T_infinity<<<{"description": "Name of a coefficient specifying the far-field temperature. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = .5
heat_transfer_coefficient<<<{"description": "Name of the coefficient specifying the heat transfer coefficient. A functor is any of the following: a variable, an MFEM material property, a function, a postprocessor or a number."}>>> = 5
[]
[top_dirichlet]
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"}>>> = temperature
boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies. Defaults to all boundaries."}>>> = '2'
[]
[]
[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"}>>>
[]
[jacobi]
type = MFEMOperatorJacobiSmoother<<<{"description": "MFEM solver for performing Jacobi smoothing of the equation system.", "href": "../../source/mfem/solvers/MFEMOperatorJacobiSmoother.html"}>>>
[]
[]
[Solver<<<{"href": "../Solver/index.html"}>>>]
type = MFEMHypreGMRES
preconditioner = boomeramg
l_tol = 1e-16
l_max_its = 1000
[]
[Executioner<<<{"href": "../Executioner/index.html"}>>>]
type = MFEMTransient
device = cpu
assembly_level = legacy
dt = 2.0
start_time = 0.0
end_time = 6.0
[]
[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/HeatTransfer
vtk_format<<<{"description": "Select VTK data format to use, choosing between BINARY, BINARY32, and ASCII."}>>> = ASCII
[]
[]
(test/tests/mfem/kernels/heattransfer.i)