Basic functionality#
Definition of a temperature field or problem is essential for hydrogen transport and FESTIM as a whole. Users can define it as a fixed value or as a function of space and time. Read more about temperature in FESTIM here. FESTIM’s unit for temperature is Kelvin.
This tutorial will go over the basics of defining the temperature in a FESTIM simulation.
Objectives:
Learn how to define homogeneous temperature
Learn how to define space and time dependent temperature functions
Learn how to define temperature as a
dolfinx.fem.Function
Defining temperatures in FESTIM#
The simplest way to define temperature is as a homogeneous value:
import festim as F
my_model = F.HydrogenTransportProblem()
my_model.temperature = 300 # K
/home/docs/checkouts/readthedocs.org/user_builds/festim-workshop/conda/latest/lib/python3.12/site-packages/festim/coupled_heat_hydrogen_problem.py:1: TqdmExperimentalWarning: Using `tqdm.autonotebook.tqdm` in notebook mode. Use `tqdm.tqdm` instead to force console mode (e.g. in jupyter console)
import tqdm.autonotebook
To define a space-dependent temperature function we can use lambda functions:
import festim as F
import ufl
T0 = 300
my_model = F.HydrogenTransportProblem()
my_model.temperature = lambda x: T0 * ufl.exp(-x[0])
For a time-dependent function such as:
my_model.temperature = lambda t: T0 * ufl.sin(t)
For a function in \( x, y, z \) space:
Note
When defining custom functions for values, only the arguments x, t and T can be defined. Spatial coordinates can be referred to by their indices, such as x[0], x[1], and x[2], regardless of the coordinate system used. Time dependence must use t, and T for temperature dependence.
We can assign a unit cube mesh and run the simulation to see how a spatially-dependent temperature field affects the transport:
from dolfinx.mesh import create_unit_cube
from mpi4py import MPI
import ufl
import festim as F
import numpy as np
mesh = F.Mesh(create_unit_cube(MPI.COMM_WORLD, 10, 10, 10))
my_model.mesh = mesh
mat = F.Material(D_0=1, E_D=0)
volume = F.VolumeSubdomain(id=1, material=mat)
top_surface = F.SurfaceSubdomain(id=1, locator=lambda x: np.isclose(x[2], 1.0))
bottom_surface = F.SurfaceSubdomain(id=2, locator=lambda x: np.isclose(x[2], 0.0))
my_model.subdomains = [top_surface, bottom_surface, volume]
H = F.Species("H")
my_model.species = [H]
my_model.boundary_conditions = [
F.FixedConcentrationBC(subdomain=top_surface, value=1.0, species=H),
F.FixedConcentrationBC(subdomain=bottom_surface, value=0.0, species=H),
]
my_model.settings = F.Settings(atol=1e-10, rtol=1e-10, transient=False)
my_model.initialise()
my_model.run()
2026-05-22 15:24:10.034 ( 0.980s) [ 75FD9A411740]vtkXOpenGLRenderWindow.:1458 WARN| bad X server connection. DISPLAY=
Defining temperature as a FEniCS function#
Users can pass a dolfinx Function to the temperature attribute, which can be helpful when using results from other programs (ie. converting an OpenFOAM output into a passable dolfinx function).
For example, we want to have a temperature function defined as:
To define the temperature as a function or expression, we can use ufl:
import dolfinx
from dolfinx.mesh import create_unit_square
from mpi4py import MPI
import festim as F
from basix.ufl import element
mesh_fenics = create_unit_square(MPI.COMM_WORLD, 20, 20)
my_model = F.HydrogenTransportProblem()
my_model.mesh = F.Mesh(mesh_fenics)
el = element("Lagrange", mesh_fenics.topology.cell_name(), 2)
V = dolfinx.fem.functionspace(my_model.mesh.mesh, el)
temperature = dolfinx.fem.Function(V)
coords = ufl.SpatialCoordinate(temperature.function_space.mesh)
x = coords[0]
y = coords[1]
interpolation = temperature.function_space.element.interpolation_points
expr = dolfinx.fem.Expression(300*ufl.exp(-((x-0.5)**2 + (y-0.5)**2)), interpolation)
temperature.interpolate(expr)
To set this temperature for our FESTIM simulation, we simply need to pass it to our model’s temperature attribute:
my_model.temperature = temperature