dolfinx.io.gmshio
Tools to extract data from Gmsh models.
Functions
|
The permutation array for permuting Gmsh ordering to DOLFINx ordering. |
|
Extract the mesh geometry from a Gmsh model. |
|
Extract all entities tagged with a physical marker in the gmsh model. |
|
Create a Mesh from a Gmsh model. |
|
Read a Gmsh .msh file and return a distributed |
|
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 thata_dolfinx[i] = a_gmsh[p[i]]
.
- dolfinx.io.gmshio.extract_geometry(model, name: Optional[str] = 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: Optional[str] = None)[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 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.
- dolfinx.io.gmshio.model_to_mesh(model, comm: ~mpi4py.MPI.Comm, rank: int, gdim: int = 3, partitioner: ~typing.Optional[~typing.Callable[[~mpi4py.MPI.Comm, int, int, ~dolfinx.cpp.graph.AdjacencyList_int32], ~dolfinx.cpp.graph.AdjacencyList_int32]] = None, dtype=<class 'numpy.float64'>) Tuple[Mesh, MeshTags_int32, MeshTags_int32] [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 thedolfinx.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:
A triplet
(mesh, cell_tags, facet_tags)
, where cell_tags hold markers for the cells and facet tags holds markers for facets (if tags are found in Gmsh model).
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
dolfinxio.XDMFFile
after creation for efficient access.
- dolfinx.io.gmshio.read_from_msh(filename: str, comm: Comm, rank: int = 0, gdim: int = 3, partitioner: Optional[Callable[[Comm, int, int, AdjacencyList_int32], AdjacencyList_int32]] = None) Tuple[Mesh, MeshTags_int32, MeshTags_int32] [source]
Read a Gmsh .msh file and return a distributed
dolfinx.mesh.Mesh
and 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:
A triplet
(mesh, cell_tags, facet_tags)
with meshtags for associated physical groups for cells and facets.
- dolfinx.io.gmshio.ufl_mesh(gmsh_cell: int, gdim: int) 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.