Transfer meshtags to submesh in dolfinx

Hello everyone,

I’m currently trying to solve an equation on only a subdomain. For this, I created a submesh using dolfinx.mesh.create_submesh. This works perfectly, however, I want to transfer the meshtags from the parent mesh to the submesh, to use them for the marking of the measure. Is there a way of achieving this in the current dolfinx version?

Best regards
Max

Are the mesh-tags of cells or facets?

I’m interested in the tags for the facets.

Hi,

I am also modelling a problem involving parent and submesh. Does anyone know how the problem described in the original question is done in dolfinx? I am interested in both cell and facet tages.

Thanks!

Hi, is this already solved in the meantime?

I would also like to create a submesh and transfer the facet tags from the parent mesh to the submesh. But it seems that this is still not done automatically.

Here is a minimal example that does this (mapping facet tags in the parent mesh to the submesh):

from IPython import embed
import dolfinx
import ufl
from mpi4py import MPI
import numpy as np


def facets(x):
    return np.isclose(x[0], 0)


mesh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 10, 10)
tdim = mesh.topology.dim
fdim = tdim - 1

left_facets = dolfinx.mesh.locate_entities_boundary(mesh, fdim, facets)


def some_facets(x):
    return (x[1] <= 0.51) & np.isclose(x[0], 0)


some_marked_facets = dolfinx.mesh.locate_entities_boundary(
    mesh, fdim, some_facets)
values = np.full_like(some_marked_facets, 1, dtype=np.int32)

facet_tags = dolfinx.mesh.meshtags(mesh, fdim, some_marked_facets, values)

with dolfinx.io.XDMFFile(mesh.comm, "facet.xdmf", "w") as xdmf:
    xdmf.write_mesh(mesh)
    xdmf.write_meshtags(facet_tags)


f_map = mesh.topology.index_map(fdim)
all_facets = f_map.size_local + f_map.num_ghosts
# Create array with zeros for all facets that are not marked
all_values = np.zeros(all_facets, dtype=np.int32)
all_values[facet_tags.indices] = facet_tags.values

submesh, entity_map, _, _ = dolfinx.mesh.create_submesh(
    mesh, fdim, left_facets)
sub_values = np.empty(len(entity_map), dtype=np.int32)
for i, entity in enumerate(entity_map):
    sub_values[i] = all_values[entity]
sub_meshtag = dolfinx.mesh.meshtags(submesh, submesh.topology.dim, np.arange(
    len(entity_map), dtype=np.int32), sub_values)

with dolfinx.io.XDMFFile(mesh.comm, "sub_tag.xdmf", "w") as xdmf:
    xdmf.write_mesh(submesh)
    xdmf.write_meshtags(sub_meshtag)
1 Like

And here is a similar version transfering facet tags from a parent mesh to a submesh of the same co-dimension:

from IPython import embed
import dolfinx
import ufl
from mpi4py import MPI
import numpy as np


def facets(x):
    return np.isclose(x[0], 0)


mesh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 10, 10)
tdim = mesh.topology.dim
fdim = tdim - 1

left_facets = dolfinx.mesh.locate_entities_boundary(mesh, fdim, facets)


def some_facets(x):
    return (x[1] <= 0.51) & np.isclose(x[0], 0)


def some_cells(x):
    return (x[0] <= 0.51)


some_marked_facets = dolfinx.mesh.locate_entities_boundary(
    mesh, fdim, some_facets)
values = np.full_like(some_marked_facets, 1, dtype=np.int32)

facet_tags = dolfinx.mesh.meshtags(mesh, fdim, some_marked_facets, values)
c_to_f = mesh.topology.connectivity(tdim, fdim)
with dolfinx.io.XDMFFile(mesh.comm, "facet.xdmf", "w") as xdmf:
    xdmf.write_mesh(mesh)
    xdmf.write_meshtags(facet_tags)


f_map = mesh.topology.index_map(fdim)
all_facets = f_map.size_local + f_map.num_ghosts
# Create array with zeros for all facets that are not marked
all_values = np.zeros(all_facets, dtype=np.int32)
all_values[facet_tags.indices] = facet_tags.values


submesh, entity_map, _, _ = dolfinx.mesh.create_submesh(
    mesh, tdim, dolfinx.mesh.locate_entities(mesh, mesh.topology.dim, some_cells))

submesh.topology.create_entities(fdim)
subf_map = submesh.topology.index_map(fdim)
submesh.topology.create_connectivity(tdim, fdim)
c_to_f_sub = submesh.topology.connectivity(tdim, fdim)
num_sub_facets = subf_map.size_local + subf_map.size_global
sub_values = np.empty(num_sub_facets, dtype=np.int32)
for i, entity in enumerate(entity_map):
    parent_facets = c_to_f.links(entity)
    child_facets = c_to_f_sub.links(i)
    for child, parent in zip(child_facets, parent_facets):
        sub_values[child] = all_values[parent]
sub_meshtag = dolfinx.mesh.meshtags(submesh, submesh.topology.dim-1, np.arange(
    num_sub_facets, dtype=np.int32), sub_values)

with dolfinx.io.XDMFFile(mesh.comm, "sub_tag.xdmf", "w") as xdmf:
    xdmf.write_mesh(submesh)
    submesh.topology.create_connectivity(fdim, tdim)
    xdmf.write_meshtags(sub_meshtag)

2 Likes

That’s great! Thank you!

Is there a way to create submeshes with the help of cell tags/facet tags instead of maker functions?

Yes, the input to create submesh is just the list of entities that you want to have in your submesh.
https://docs.fenicsproject.org/dolfinx/v0.7.2/cpp/mesh.html#_CPPv4I_NSt14floating_pointEEN7dolfinx4mesh14create_submeshENSt5tupleI4MeshI1TENSt6vectorINSt7int32_tEEENSt6vectorINSt7int32_tEEENSt6vectorINSt7int32_tEEEEERK4MeshI1TEiNSt4spanIKNSt7int32_tEEE

1 Like