Extra topic: Reading in meshes in parallel

Extra topic: Reading in meshes in parallel#

For more information about mesh-partitioning and mesh input, see: https://jsdokken.com/dolfinx_docs/meshes.html

Cell ownership#

For the remainder of this section we will consider a 3x3 unit square mesh:

from mpi4py import MPI
import dolfinx
import ipyparallel as ipp
domain = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 3, 3)

The mesh consists of cells, edges and vertices. A mesh is created by supplying the information regarding the connectivity between the cells and the mesh nodes. As a higher order mesh has more nodes than vertices, we can get the connectivity between the cells and mesh vertices through the mesh-topology.

We start by creating a simple function for inspecting these outputs in serial and parallel

def inspect_mesh(shared_facet: bool = False):
    from mpi4py import MPI

    import dolfinx

    ghost_mode = dolfinx.mesh.GhostMode.shared_facet if shared_facet else dolfinx.mesh.GhostMode.none
    domain = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 3, 3, ghost_mode=ghost_mode)
    topology = domain.topology
    tdim = topology.dim
    print(f"Number of cells in process: {topology.index_map(tdim).size_local}", flush=True)
    print(f"Number of shared cells: {topology.index_map(tdim).num_ghosts}", flush=True)
    print(f"Local range [{topology.index_map(tdim).local_range[0]}, {topology.index_map(tdim).local_range[1]})")
    cell_to_vertices = topology.connectivity(tdim, 0)
    print(cell_to_vertices)

We start by inspecting the outputs in serial

inspect_mesh()
Number of cells in process: 18
Number of shared cells: 0
Local range [0, 18)
<AdjacencyList> with 18 nodes
  0: [0 1 2 ]
  1: [0 3 2 ]
  2: [4 0 3 ]
  3: [3 2 5 ]
  4: [4 6 3 ]
  5: [3 7 5 ]
  6: [8 4 6 ]
  7: [6 3 7 ]
  8: [7 5 9 ]
  9: [8 10 6 ]
  10: [6 11 7 ]
  11: [7 12 9 ]
  12: [10 6 11 ]
  13: [11 7 12 ]
  14: [10 13 11 ]
  15: [11 14 12 ]
  16: [13 11 14 ]
  17: [13 15 14 ]

We observe that we have 18 cells in the mesh, that are connected to the 18 vertices. All indices start from 0, and can be mapped to its global owner by calling domain.topology.index_map(tdim).local_to_global([idx])

Parallel execution#

We use IPython-Parallel for running DOLFINx on multiple processes inside our script.

cluster = ipp.Cluster(engines="mpi", n=2)
rc = cluster.start_and_connect_sync()
Starting 2 engines with <class 'ipyparallel.cluster.launcher.MPIEngineSetLauncher'>

This has started to processes that can execute code with a MPI communicator. We run the mesh code on two processes, instructing DOLFINx not to share any cells between the processes

query = rc[:].apply_async(inspect_mesh, False)
query.wait()
query.display_outputs()
[stdout:0] 
Number of cells in process: 9
Number of shared cells: 0
Local range [0, 9)
<AdjacencyList> with 9 nodes
  0: [0 7 8 ]
  1: [0 1 8 ]
  2: [1 8 2 ]
  3: [1 3 2 ]
  4: [8 2 4 ]
  5: [3 2 5 ]
  6: [2 4 9 ]
  7: [3 6 5 ]
  8: [2 5 9 ]

[stdout:1] 
Number of cells in process: 9
Number of shared cells: 0
Local range [9, 18)
<AdjacencyList> with 9 nodes
  0: [0 1 2 ]
  1: [0 3 2 ]
  2: [4 0 3 ]
  3: [3 2 5 ]
  4: [4 6 3 ]
  5: [3 9 5 ]
  6: [6 3 9 ]
  7: [9 5 7 ]
  8: [9 8 7 ]

Next, we instruct DOLFINx to share cells between two processes if a facet is shared between the two processes.

new_query = rc[:].apply_async(inspect_mesh, True)
new_query.wait()
new_query.display_outputs()
[stdout:0] 
Number of cells in process: 9
Number of shared cells: 3
Local range [0, 9)
<AdjacencyList> with 12 nodes
  0: [0 7 8 ]
  1: [0 1 8 ]
  2: [1 8 2 ]
  3: [1 3 2 ]
  4: [8 2 4 ]
  5: [3 2 5 ]
  6: [2 4 9 ]
  7: [3 6 5 ]
  8: [2 5 9 ]
  9: [7 8 10 ]
  10: [8 10 4 ]
  11: [4 9 11 ]

[stdout:1] 
Number of cells in process: 9
Number of shared cells: 3
Local range [9, 18)
<AdjacencyList> with 12 nodes
  0: [0 1 2 ]
  1: [0 3 2 ]
  2: [4 0 3 ]
  3: [3 2 5 ]
  4: [4 6 3 ]
  5: [3 9 5 ]
  6: [6 3 9 ]
  7: [9 5 7 ]
  8: [9 8 7 ]
  9: [10 4 6 ]
  10: [6 11 9 ]
  11: [11 9 8 ]
cluster.stop_cluster_sync()
Hide code cell output
Stopping controller
Stopping engine(s): 1728596943