Boundary conditions#

The boundary conditions (BCs) are essential to FESTIM simulations. They describe the mathematical problem at the boundaries of the simulated domain. If no BC is set on a boundary, it is assumed that the flux is null. This is also called a symmetry BC.

All boundary conditions in FESTIM require a surface subdomain. The subdomain can be a surface, an edge or a point and can be defined with the festim.SurfaceSubdomain class. See the Surface Subdomains section for more information on how to define subdomains.

Hydrogen transport BCs#

Some BCs are specific to hydrogen transport. FESTIM provides a handful of convenience classes making things a bit easier for the users.

Imposing the concentration#

The concentration of a defined species can be imposed on boundaries with festim.FixedConcentrationBC.

from festim import FixedConcentrationBC, Species, SurfaceSubdomain

boundary = SurfaceSubdomain(id=1)
H = Species(name="Hydrogen")

my_bc = FixedConcentrationBC(subdomain=boundary, value=10, species=H)

The imposed concentration can be dependent on space, time and temperature:

from festim import FixedConcentrationBC, Species, SurfaceSubdomain

boundary = SurfaceSubdomain(id=1)
H = Species(name="Hydrogen")

my_custom_value = lambda x, t, T: 10 + x[0]**2 + t + T

my_bc = FixedConcentrationBC(subdomain=boundary, value=my_custom_value, species=H)

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.

Imposing a particle flux#

When a particle flux needs to be imposed on a boundary, use the festim.ParticleFluxBC class.

from festim import ParticleFluxBC, Species, SurfaceSubdomain

boundary = SurfaceSubdomain(id=1)
H = Species(name="Hydrogen")

my_flux_bc = ParticleFluxBC(subdomain=boundary, value=2, species=H)

As for the fixed concentration boundary condition, the flux can be dependent on space, time and temperature. But for particle fluxes, the values can also be dependent on a species’ concentration:

from festim import ParticleFluxBC, Species, SurfaceSubdomain

boundary = SurfaceSubdomain(id=1)
H = Species(name="Hydrogen")

my_custom_value = lambda t, c: 10*t**2 + 2*c

my_flux_bc = ParticleFluxBC(
    subdomain=boundary,
    value=my_custom_value,
    species=H,
    species_dependent_value={"c": H},
)

Note

The species_dependent_value argument requires a dictionary to be passed, mapping any arguments in the custom function given to value, to any species defined.

For instance with three species A, B and C, the dictionary can be defined as:

from festim import Species

A = Species(name="A")
B = Species(name="B")
C = Species(name="C")

my_custom_value = lambda c_A, c_B, c_C: 2*c_A + 3*c_B + 4*c_C

species_dependent_value = {"c_A": A, "c_B": B, "c_C": C}

Sievert’s law of solubility#

Impose the concentration of a species as \(c_\mathrm{m} = S(T) \sqrt{P}\) where \(S\) is the Sievert’s solubility and \(P\) is the partial pressure of the species on this surface (see festim.SievertsBC).

from festim import SievertsBC, SurfaceSubdomain, Species

boundary = SurfaceSubdomain(id=1)
H = Species(name="Hydrogen")

custom_pressure_value = lambda t: 2 + t

my_bc = SievertsBC(subdomain=3, S_0=2, E_S=0.1, species=H, pressure=custom_pressure_value)

Henry’s law of solubility#

Similarly, the the concentration of a species can be set from Henry’s law of solubility \(c_\mathrm{m} = K_H P\) where \(K_H\) is the Henry solubility (see festim.HenrysBC).

from festim import HenrysBC, SurfaceSubdomain, Species

boundary = SurfaceSubdomain(id=1)
H = Species(name="Hydrogen")

pressure_value = lambda t: 5 * t

my_bc = HenrysBC(subdomain=3, H_0=1.5, E_H=0.2, species=H, pressure=pressure_value)

Surface reactions#

Surface reactions on boundary can be defined with the festim.SurfaceReactionBC class.

The surface reaction class can be used to impose dissociation and recombination reactions on the surface of the material. A reaction is defined by specifying the reactants, products, and forward/backward rate constants. For example:

from festim import Species, SurfaceReactionBC, SurfaceSubdomain

boundary = SurfaceSubdomain(id=1)
A = Species("A")
B = Species("B")
C = Species("C")

my_bc = SurfaceReactionBC(
    reactant=[A, B],
    gas_pressure=1e5,
    k_r0=1,
    E_kr=0.1,
    k_d0=1e-5,
    E_kd=0.1,
    subdomain=boundary,
)

The net reaction rate is:

\[R = K_r c_A c_B - K_d c_C\]

where \(K_r\) and \(K_d\) are the temperature-dependent forward and backward rate constants, respectively, and \(c_A\), \(c_B\), and \(c_C\) are the concentrations of species A, B, and C at the surface. From this FESTIM automatically applies the corresponding fluxes as Neumann boundary conditions (see festim.ParticleFluxBC).

Recombination and Dissociation#

Hydrogen recombination/dissociation can be modelled as:

\[\mathrm{H + H} \overset{K_r}{\underset{K_d}{\rightleftharpoons}} \mathrm{H_2}\]

If the partial pressure of \(\mathrm{H}_2\) is known or assumed constant, the product species does not need to be included explicitly, and the flux simplifies to:

\[\mathbf{J}_{\mathrm{H}} \cdot \mathbf{n} = K_r c_{\mathrm{H}}^2 - K_d P_{\mathrm{H_2}}\]

Both rate coefficients can follow Arrhenius laws.

from festim import Species, SurfaceReactionBC, SurfaceSubdomain

boundary = SurfaceSubdomain(id=1)
H = Species("H")

my_bc = SurfaceReactionBC(
    reactant=[H, H],
    gas_pressure=1e5,
    k_r0=1,
    E_kr=0.1,
    k_d0=1e-5,
    E_kd=0.1,
    subdomain=boundary,
)

Multiple isoptopes#

Surface reactions can involve multiple hydrogen isotopes, enabling the modelling of more complex interactions between species. For example, in a system with both mobile hydrogen and tritium, various molecular recombination pathways may occur at the surface, resulting in the formation of \(\mathrm{H_2}\), \(\mathrm{T_2}\), and \(\mathrm{HT}\):

\[\mathrm{H + H} \rightleftharpoons \mathrm{H_2}, \quad \mathrm{T + T} \rightleftharpoons \mathrm{T_2}, \quad \mathrm{H + T} \rightleftharpoons \mathrm{HT}\]
from festim import Species, SurfaceReactionBC, SurfaceSubdomain

boundary = SurfaceSubdomain(id=1)
H = Species("H")
T = Species("T")

reac1 = SurfaceReactionBC(
    reactant=[H, H],
    gas_pressure=1e6,
    k_r0=1.0,
    E_kr=0.1,
    k_d0=0.5,
    E_kd=0.1,
    subdomain=boundary,
)
reac2 = SurfaceReactionBC(
    reactant=[T, T],
    gas_pressure=1e6,
    k_r0=1.0,
    E_kr=0.1,
    k_d0=0.5,
    E_kd=0.1,
    subdomain=boundary,
)
reac3 = SurfaceReactionBC(
    reactant=[H, T],
    gas_pressure=1e6,
    k_r0=1.0,
    E_kr=0.1,
    k_d0=0.5,
    E_kd=0.1,
    subdomain=boundary,
)

Heat transfer BCs#

Some BCs are specific to heat transfer. FESTIM provides a handful of convenience classes making things a bit easier for the users.

Imposing the temperature#

The temperature can be imposed on boundaries with festim.FixedTemperatureBC.

from festim import FixedTemperatureBC, SurfaceSubdomain

boundary = SurfaceSubdomain(id=1)

my_bc = FixedTemperatureBC(subdomain=boundary, value=10)

To define the temperature as space or time dependent, a function can be passed to the value argument:

from festim import FixedTemperatureBC, SurfaceSubdomain

boundary = SurfaceSubdomain(id=1)

my_custom_value = lambda x, t: 10 + x[0]**2 + t

my_bc = FixedTemperatureBC(subdomain=boundary, value=my_custom_value)

Imposing a heat flux#

When a heat flux needs to be imposed on a boundary, use the festim.HeatFluxBC class.

from festim import HeatFluxBC, SurfaceSubdomain

boundary = SurfaceSubdomain(id=1)

my_flux_bc = HeatFluxBC(subdomain=boundary, value=5)

As for the fixed temperature boundary condition, the flux can be dependent on space and time. But for heat fluxes, the values can also be dependent on a temperature:

from festim import HeatFluxBC, SurfaceSubdomain

boundary = SurfaceSubdomain(id=1)

my_custom_value = lambda x, t, T: 2 * x[0] + 10 * t + T

my_flux_bc = HeatFluxBC(subdomain=boundary, value=my_custom_value)