Steady State Diffusion

Summary

Solves a steady state diffusion problem for the temperature profile across a mug with fixed temperatures on two boundaries.

Description

This problem solves the steady state heat equation with strong form:

where , on the base of the mug and at the top of the mug.

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
  []
  [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"}>>>]
  [concentration]
    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
  []
[]

[AuxVariables<<<{"href": "../AuxVariables/index.html"}>>>]
  [concentration_gradient]
    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"}>>>]
  [grad]
    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"}>>> = concentration_gradient
    source<<<{"description": "Scalar H1 MFEMVariable to take the gradient of."}>>> = concentration
    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
  []
[]

[ICs<<<{"href": "../ICs/index.html"}>>>]
  [diffused_ic]
    type = MFEMScalarIC<<<{"description": "Sets the initial values of an MFEM scalar variable from a user-specified scalar coefficient.", "href": "../../source/mfem/ics/MFEMScalarIC.html"}>>>
    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."}>>> = one
    variable<<<{"description": "The variable to apply the initial condition on."}>>> = concentration
  []
[]

[Functions<<<{"href": "../Functions/index.html"}>>>]
  [one]
    type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../source/functions/MooseParsedFunction.html"}>>>
    expression<<<{"description": "The user defined function."}>>> = 1.0
  []
[]

[BCs<<<{"href": "../BCs/index.html"}>>>]
  [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"}>>> = concentration
    boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies. Defaults to all boundaries."}>>> = 'bottom'
    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]
    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"}>>> = concentration
    boundary<<<{"description": "The list of boundaries (ids or names) from the mesh where this object applies. Defaults to all boundaries."}>>> = 'top'
  []
[]

[FunctorMaterials<<<{"href": "../FunctorMaterials/index.html"}>>>]
  [Substance]
    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"}>>> = diffusivity
    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."}>>> = 'the_domain'
  []
[]

[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"}>>> = concentration
    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."}>>> = diffusivity
  []
[]

[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 = MFEMSteady
  device = cpu
[]

[Outputs<<<{"href": "../Outputs/index.html"}>>>]
  active<<<{"description": "If specified only the blocks named will be visited and made active"}>>> = ParaViewDataCollection
  [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/Diffusion
    vtk_format<<<{"description": "Select VTK data format to use, choosing between BINARY, BINARY32, and ASCII."}>>> = ASCII
  []
  [VisItDataCollection]
    type = MFEMVisItDataCollection<<<{"description": "Output for controlling export to an mfem::VisItDataCollection.", "href": "../../source/mfem/outputs/MFEMVisItDataCollection.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/VisItDataCollection
  []
  [ConduitDataCollection]
    type = MFEMConduitDataCollection<<<{"description": "Output for controlling MFEMConduitDataCollection inherited data.", "href": "../../source/mfem/outputs/MFEMConduitDataCollection.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/ConduitDataCollection/Run
    protocol<<<{"description": "Conduit relay I/O protocol to use. Options: hdf5 (default), json, conduit_json, conduit_bin."}>>> = conduit_bin
  []
[]
(test/tests/mfem/kernels/diffusion.i)