Transient Heat Transfer (Mixed Form)
Summary
Solves a transient heat conduction problem using a mixed weak form, representing temperature on piecewise constant conforming finite elements, and heat fluxes on conforming Raviart-Thomas elements.
This example is based on the Firedrake Irksome demo_mixed_heat example.
Description
This problem solves the transient heat equation with strong form:
where and , subject to the initial condition .
In this example, we solve this by introducing the heat flux to generate the mixed weak form of the transient heat equation
where
and the polynomial order of the FEs used to represent and is the same as and .
Notably, in contrast to the primal weak form for the transient heat equation solved here, the heat flux on is strongly imposed via a Dirichlet condition, and temperatures on boundaries are imposed weakly via boundary integrals.
Example File
# Mixed heat transfer problem.
# Based on Firedrake Irksome demo_mixed_heat example:
# https://www.firedrakeproject.org/Irksome/demos/demo_mixed_heat.py.html
[Mesh<<<{"href": "../Mesh/index.html"}>>>]
type = MFEMMesh
file = ../mesh/square.e
[]
[Problem<<<{"href": "../Problem/index.html"}>>>]
type = MFEMProblem
[]
[FESpaces]
[HDivFESpace]
type = MFEMVectorFESpace
fec_type = RT
fec_order = FIRST
[]
[L2FESpace]
type = MFEMScalarFESpace
fec_type = L2
fec_order = FIRST
[]
[]
[Variables<<<{"href": "../Variables/index.html"}>>>]
[time_integrated_heat_flux]
type = MFEMVariable
fespace = HDivFESpace
time_derivative = heat_flux
[]
[temperature]
type = MFEMVariable
fespace = L2FESpace
[]
[]
[Kernels<<<{"href": "../Kernels/index.html"}>>>]
[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<<<{"href": "../BCs/index.html"}>>>]
[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<<<{"href": "../Executioner/index.html"}>>>]
type = MFEMTransient
device = cpu
assembly_level = legacy
dt = 0.03
start_time = 0.0
end_time = 0.09
[]
[Outputs<<<{"href": "../Outputs/index.html"}>>>]
[ParaViewDataCollection]
type = MFEMParaViewDataCollection
file_base = OutputData/MixedHeatTransfer
vtk_format = ASCII
[]
[](moose/test/tests/mfem/kernels/mixed_heattransfer.i)