dolfinx.io.gmshio

Tools to extract data from Gmsh models.

Functions

cell_perm_array(cell_type, num_nodes)

The permutation array for permuting Gmsh ordering to DOLFINx ordering.

extract_geometry(model[, name])

Extract the mesh geometry from a Gmsh model.

extract_topology_and_markers(model[, name])

Extract all entities tagged with a physical marker in the gmsh model.

model_to_mesh(model, comm, rank[, gdim, ...])

Create a Mesh from a Gmsh model.

read_from_msh(filename, comm[, rank, gdim, ...])

Read a Gmsh .msh file and return a dolfinx.mesh.Mesh and cell facet markers.

ufl_mesh(gmsh_cell, gdim, dtype)

Create a UFL mesh from a Gmsh cell identifier and geometric dimension.

dolfinx.io.gmshio.cell_perm_array(cell_type: CellType, num_nodes: int) list[int][source]

The permutation array for permuting Gmsh ordering to DOLFINx ordering.

Parameters:
  • cell_type – DOLFINx cell type.

  • num_nodes – Number of nodes in the cell.

Returns:

An array p such that a_dolfinx[i] = a_gmsh[p[i]].

dolfinx.io.gmshio.extract_geometry(model, name: str | None = None) ndarray[Any, dtype[float64]][source]

Extract the mesh geometry from a Gmsh model.

Returns an array of shape (num_nodes, 3), where the i-th row corresponds to the i-th node in the mesh.

Parameters:
  • model – Gmsh model

  • name – Name of the Gmsh model. If not set the current model will be used.

Returns:

The mesh geometry as an array of shape (num_nodes, 3).

dolfinx.io.gmshio.extract_topology_and_markers(model, name: str | None = None) tuple[dict[int, TopologyDict], dict[str, tuple[int, int]]][source]

Extract all entities tagged with a physical marker in the gmsh model.

Returns a nested dictionary where the first key is the gmsh MSH element type integer. Each element type present in the model contains the cell topology of the elements and corresponding markers.

Parameters:
  • model – Gmsh model.

  • name – Name of the gmsh model. If not set the current model will be used.

Returns:

A tuple (topologies, physical_groups), where topologies is a nested dictionary where each key corresponds to a gmsh cell type. Each cell type found in the mesh has a 2D array containing the topology of the marked cell and a list with the corresponding markers. physical_groups is a dictionary where the key is the physical name and the value is a tuple with the dimension and tag.

dolfinx.io.gmshio.model_to_mesh(model, comm: ~mpi4py.MPI.Comm, rank: int, gdim: int = 3, partitioner: ~typing.Callable[[~mpi4py.MPI.Comm, int, int, ~dolfinx.cpp.graph.AdjacencyList_int32], ~dolfinx.cpp.graph.AdjacencyList_int32] | None = None, dtype=<class 'numpy.float64'>) MeshData[source]

Create a Mesh from a Gmsh model.

Creates a dolfinx.mesh.Mesh from the physical entities of the highest topological dimension in the Gmsh model. In parallel, the gmsh model is processed on one MPI rank, and the dolfinx.mesh.Mesh is distributed across ranks.

Parameters:
  • model – Gmsh model.

  • comm – MPI communicator to use for mesh creation.

  • rank – MPI rank that the Gmsh model is initialized on.

  • gdim – Geometrical dimension of the mesh.

  • partitioner – Function that computes the parallel distribution of cells across MPI ranks.

Returns:

MeshData with mesh, cell tags, facet tags, edge tags, vertex tags and physical groups.

Note

For performance, this function should only be called once for large problems. For re-use, it is recommended to save the mesh and corresponding tags using dolfinx.io.XDMFFile after creation for efficient access.

dolfinx.io.gmshio.read_from_msh(filename: str | Path, comm: Comm, rank: int = 0, gdim: int = 3, partitioner: Callable[[Comm, int, int, AdjacencyList], AdjacencyList_int32] | None = None) MeshData[source]

Read a Gmsh .msh file and return a dolfinx.mesh.Mesh and cell facet markers.

Note

This function requires the Gmsh Python module.

Parameters:
  • filename – Name of .msh file.

  • comm – MPI communicator to create the mesh on.

  • rank – Rank of comm responsible for reading the .msh file.

  • gdim – Geometric dimension of the mesh

Returns:

Meshdata with mesh, cell tags, facet tags, edge tags, vertex tags and physical groups.

dolfinx.io.gmshio.ufl_mesh(gmsh_cell: int, gdim: int, dtype: dtype[Any] | None | type[Any] | _SupportsDType[dtype[Any]] | str | tuple[Any, int] | tuple[Any, SupportsIndex | Sequence[SupportsIndex]] | list[Any] | _DTypeDict | tuple[Any, Any]) Mesh[source]

Create a UFL mesh from a Gmsh cell identifier and geometric dimension.

See https://gmsh.info//doc/texinfo/gmsh.html#MSH-file-format.

Parameters:
  • gmsh_cell – Gmsh cell identifier.

  • gdim – Geometric dimension of the mesh.

Returns:

UFL Mesh using Lagrange elements (equispaced) of the corresponding DOLFINx cell.