phiFEM v0.7.0

Author

Raphaël Bulle

Published

April 17, 2026

compute_tags_measures API

The main phiFEM function is compute_tags_measures (source).

compute_tags_measures(mesh: dolfinx.mesh.Mesh,
                      discrete_levelset: dolfinx.fem.Function,
                      detection_degree: int,
                      box_mode: bool = False,
                      single_layer_cut: bool = False) ->
                      Tuple[
                        MeshTags,
                        MeshTags,
                        Mesh | None,
                        ufl.Measure,
                        list[npt.NDArray[np.int32]] | None]

Args:

  • mesh: dolfinx.mesh.Mesh: the background mesh.
  • discrete_levelset: dolfinx.fem.Function: the interpolation of the levelset function.
  • detection_degree: int: an integer determining the number of quadrature points used on the boundary of the cells to detect cut cells (detection_degree should not be larger than the degree of discrete_levelset).
  • box_mode: bool: if False a submesh is computed by discarding the outside cells and the tags are computed on this submesh. If True, no submesh is computed and the tags are computed on the background mesh.
  • single_layer_cut: bool: if True, the layer of cut cells is forced to be only of one cell width. If False, the layer of cut cells might be larger (this can lead to unwanted behavior of the numerical scheme).

Returns:

  • The cells tags of the background mesh or submesh depending on box_mode (see Cells tags section).
  • The facets tags of the background mesh or submesh depending on box_mode (see Facets tags section).
  • The submesh or None depending on box_mode.
  • The boundaries measure a dS measure with properly oriented facet normal if box_mode=True or a ds measure if box_mode=False (see Facets tags section).
  • The list of maps returned by dolfinx.mesh.create_submesh (see dolfinx documentation) if box_mode=False. None if box_mode=True.

phiFEM MeshTags and Measures

(a) phiFEM cells tags on a background mesh.
(b) phiFEM facets tags on a background mesh.
Figure 1: The exact and discrete boundaries \(\class{levelset}{\{\varphi = 0\}}\) and \(\class{discrete-levelset}{\{\varphi_h = 0\}}\) are represented light and dark blue respectively.

Cells tags (Figure 1 (a))

Consider dx = ufl.Measure("dx", domain=background_mesh, subdomain_data=cells_tags).

 

Tag n\(^{\circ}\) Position Usage
\(\mathbf{\class{inside-cells}{1}}\) Inside \(\{\varphi_h < 0\}\) dx(1)
\(\mathbf{\class{cut-cells}{2}}\) Cut by \(\{\varphi_h = 0\}\) dx(2)
\(\mathbf{\class{outside-cells}{3}}\) Ouside \(\{\varphi_h < 0\}\) dx(3)

 

Facets tags (Figure 1 (b))

Consider dS = ufl.Measure("dS", domain=background_mesh, subdomain_data=facets_tags).

 

Tag n\(^{\circ}\) Position Usage
\(\mathbf{1}\) Inside \(\{\varphi_h < 0\}\) dS(1)
\(\mathbf{\class{cut-facets}{2}}\) Cut by \(\{\varphi_h = 0\}\) dS(2)
\(\mathbf{\class{inside-bd-facets}{3}}\) Shared by an interior and a cut cell dS(3)
\(\mathbf{\class{outside-bd-facets}{4}}\) Shared by an exterior and a cut cell dS(4)
\(\mathbf{\class{outside-facets}{5}}\) Outside \(\{\varphi_h <0 \}\) dS(5)
\(\mathbf{\class{direct-facets}{6}}\) Fitted facet of \(\{\varphi = 0 \}\) ⚠️

 

Notes:

  • ⚠️ The fitted facets with tag \(\mathbf{\class{direct-facets}{6}}\) are not represented on Figure 1 (b) and correspond to a limiting case where some facets of the background mesh are subsets of \(\{\varphi_h = 0\}\). We do not consider this as a reliable tag.
  • The facets of tags \(\mathbf{\class{inside-bd-facets}{3}}\) and \(\mathbf{\class{outside-bd-facets}{4}}\) can not be used to compute integrals involving oriented normals (see Boundaries measure section).

Boundaries measure

Let us denote \(\Omega_h\) the union of inside and cut cells (red and yellow cells in Figure 1 (a)). In \(\varphi\)-FEM schemes, we often need to compute integrals over \(\partial \Omega_h\) (orange facets in Figure 1 (b)) involving the outward from \(\Omega_h\) pointing unit normal, like e.g. \[ \int_{\partial \Omega_h} \nabla u_h \cdot n v_h. \tag{1}\]

The integral in (1) can not be computed using the above dS(4) since the facet normal n = ufl.FacetNormal(background_mesh) computed by FEniCSx is oriented randomly and not always pointing outward \(\Omega_h\).

To compute (1) we need to use the boundaries measure returned by compute_tags_measures:

import ufl
from phifem.mesh_scripts import compute_tags_measures

d_bdy = compute_tags_measures(background_mesh, discrete_levelset, detection_degree, box_mode=True)[3]
n = ufl.FacetNormal(background_mesh)
Omega_h_bdy_integral = ufl.inner(ufl.inner(ufl.grad(u_h), n), v_h) * d_bdy(100)

 

Tag n\(^{\circ}\) Position Normal orientation Usage
\(\mathbf{100}\) Shared by an exterior and a cut cell Outward \(\Omega_h\) d_bdy(100)
\(\mathbf{101}\) Shared by an interior and a cut cell Inward \(\Omega_h\) d_bdy(101)