dolfinx.mesh
Creation, refining and marking of meshes.
Functions
|
Create a :class:dolfinx.mesh.MeshTags` object that associates data with a subset of mesh entities, where the entities are defined by their vertices. |
|
Compute mesh entities satisfying a geometric marking function. |
|
Compute mesh entities that are connected to an owned boundary facet and satisfy a geometric marking function. |
|
Refine a mesh. |
|
Create a mesh from topology and geometry arrays. |
|
Create a mesh with specified entities from another mesh. |
|
Create a MeshTags object that associates data with a subset of mesh entities. |
|
Compute the midpoints of a set of mesh entities. |
|
Compute the indices of all exterior facets that are owned by the caller. |
|
Compute all entities of |
|
Create an interval mesh. |
|
Create a mesh on the unit interval. |
|
Create a rectangle mesh. |
|
Create a mesh of a unit square. |
|
Create a box mesh. |
|
Create a mesh of a unit cube. |
|
Generate cell mesh tags on a refined mesh from the mesh tags on the coarse parent mesh. |
|
Compute the geometric DOFs associated with the closure of the given mesh entities. |
|
Create a Geometry object. |
Classes
|
A mesh. |
|
Mesh tags associate data (markers) with a subset of mesh entities of a given dimension. |
|
|
|
|
|
The geometry of a |
|
Topology for a |
- class dolfinx.mesh.CellType(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)
Bases:
Enum
- hexahedron = -8
- interval = 2
- property name
(self) -> object
- point = 1
- prism = -6
- pyramid = -5
- quadrilateral = -4
- tetrahedron = 4
- triangle = 3
- class dolfinx.mesh.Geometry(geometry: Geometry_float32 | Geometry_float64)[source]
Bases:
object
The geometry of a
dolfinx.mesh.Mesh
Initialize a geometry from a C++ geometry.
- Parameters:
geometry – The C++ geometry object.
Note
Geometry objects should usually be constructed with the
create_geometry()
and not using this class initializer. This class is combined with different base classes that depend on the scalar type used in the Geometry.- property cmap: CoordinateElement
Element that describes the geometry map.
- property dim
Dimension of the Euclidean coordinate system.
- property dofmap: ndarray[Any, dtype[int32]]
Dofmap for the geometry, shape
(num_cells, dofs_per_cell)
.
- property input_global_indices: ndarray[Any, dtype[int64]]
Global input indices of the geometry nodes.
- property x: ndarray[Any, dtype[float32]] | ndarray[Any, dtype[float64]]
Geometry coordinate points,
shape=(num_points, 3)
.
- class dolfinx.mesh.GhostMode(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)
Bases:
Enum
- none = 0
- class dolfinx.mesh.Mesh(msh, domain: Mesh | None)[source]
Bases:
object
A mesh.
Initialize mesh from a C++ mesh.
- Parameters:
msh – A C++ mesh object.
domain – A UFL domain.
Note
Mesh objects should usually be constructed using
create_mesh()
and not using this class initializer. This class is combined with different base classes that depend on the scalar type used in the Mesh.- property comm
- h(dim: int, entities: ndarray[Any, dtype[int32]]) ndarray[Any, dtype[float64]] [source]
Geometric size measure of cell entities.
- Parameters:
dim – Topological dimension of the entities to compute the size measure of.
entities – Indices of entities of dimension
dim
to compute size measure of.
- Returns:
Size measure for each requested entity.
- property name
- class dolfinx.mesh.MeshTags(meshtags)[source]
Bases:
object
Mesh tags associate data (markers) with a subset of mesh entities of a given dimension.
Initialize tags from a C++ MeshTags object.
- Parameters:
meshtags – C++ mesh tags object.
Note
MeshTags objects should not usually be created using this initializer directly.
A Python mesh is passed to the initializer as it may have UFL data attached that is not attached the C + + Mesh that is associated with the C + +
meshtags
object. If mesh is passed,mesh
andmeshtags
must share the same C + + mesh.- property dim: int
Topological dimension of the tagged entities.
- find(value) ndarray[Any, dtype[int32]] [source]
Get a list of all entity indices with a given value.
- Parameters:
value – Tag value to search for.
- Returns:
Indices of entities with tag
value
.
- property indices: ndarray[Any, dtype[int32]]
Indices of tagged mesh entities.
- property name: str
Name of the mesh tags object.
- property values
Values associated with tagged mesh entities.
- class dolfinx.mesh.Topology(topology: Topology)[source]
Bases:
object
Topology for a
dolfinx.mesh.Mesh
Initialize a topology from a C++ topology. :param topology: The C++ topology object
Note
Topology objects should usually be constructed with the
dolfinx.cpp.mesh.create_topology()
and not this class initializer.- property comm
- connectivity(d0: int, d1: int) AdjacencyList_int32 [source]
Return connectivity from entities of dimension
d0
to entities of dimensiond1
.- Parameters:
d0 – Dimension of entity one is mapping from
d1 – Dimension of entity one is mapping to
- create_connectivity(d0: int, d1: int)[source]
Create connectivity between given pair of dimensions,
d0
andd1
.- Parameters:
d0 – Dimension of entities one is mapping from
d1 – Dimension of entities one is mapping to
- create_entities(dim: int) int [source]
Create entities of given topological dimension.
- Parameters:
dim – Topological dimension
- Returns:
Number of newly created entities, returns -1 if entities already existed
- property dim: int
Return the topological dimension of the mesh.
- property entity_types: list[list[CellType]]
Entity types in the topology for all topological dimensions.
- get_cell_permutation_info() ndarray[Any, dtype[uint32]] [source]
Get permutation information.
The returned data is used for packing coefficients and assembling of tensors. The bits of each integer encodes the number of reflections and permutations for each sub-entity of the cell to be able to map it to the reference element.
- get_facet_permutations() ndarray[Any, dtype[uint8]] [source]
Get the permutation integer to apply to facets.
The bits of each integer describes the number of reflections and rotations that has to be applied to a facet to map between a facet in the mesh (relative to a cell) and the corresponding facet on the reference element. The data has the shape
(num_cells, num_facets_per_cell)
, flattened row-wise. The number of cells include potential ghost cells.Note
The data can be unpacked with
numpy.unpack_bits
.
- index_map(dim: int) IndexMap [source]
Get the IndexMap that described the parallel distribution of the mesh entities.
- Parameters:
dim – Topological dimension
- Returns:
Index map for the entities of dimension
dim
.
- interprocess_facets() ndarray[Any, dtype[int32]] [source]
List of inter-process facets, if facet topology has been computed.
- property original_cell_index: ndarray[Any, dtype[int64]]
Get the original cell index
- set_connectivity(graph: AdjacencyList_int32, d0: int, d1: int)[source]
Set connectivity for given pair of topological dimensions.
- Parameters:
graph – Connectivity graph
d0 – Topological dimension mapping from
d1 – Topological dimension mapping to
- dolfinx.mesh.compute_incident_entities(topology: Topology, entities: ndarray[Any, dtype[int32]], d0: int, d1: int) ndarray[Any, dtype[int32]] [source]
Compute all entities of
d1
connected toentities
of dimensiond0
.- Parameters:
topology – The topology.
entities – List of entities fo dimension
d0
.d0 – Topological dimension.
d1 – Topological dimension to map to.
- Returns:
Incident entity indices.
- dolfinx.mesh.compute_midpoints(msh: Mesh, dim: int, entities: ndarray[Any, dtype[int32]])[source]
Compute the midpoints of a set of mesh entities.
- Parameters:
msh – The mesh.
dim – Topological dimension of the mesh entities to consider.
entities – Indices of entities in
mesh
to consider.
- Returns:
Midpoints of the entities, shape
(num_entities, 3)
.
- dolfinx.mesh.create_box(comm: ~mpi4py.MPI.Comm, points: list[~collections.abc.Buffer | ~numpy._typing._array_like._SupportsArray[~numpy.dtype[~typing.Any]] | ~numpy._typing._nested_sequence._NestedSequence[~numpy._typing._array_like._SupportsArray[~numpy.dtype[~typing.Any]]] | bool | int | float | complex | str | bytes | ~numpy._typing._nested_sequence._NestedSequence[bool | int | float | complex | str | bytes]], n: list, cell_type=CellType.tetrahedron, dtype: ~numpy.dtype[~typing.Any] | None | type[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | str | tuple[~typing.Any, int] | tuple[~typing.Any, ~typing.SupportsIndex | ~collections.abc.Sequence[~typing.SupportsIndex]] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | tuple[~typing.Any, ~typing.Any] = <class 'numpy.float64'>, ghost_mode=GhostMode.shared_facet, partitioner=None) Mesh [source]
Create a box mesh.
- Parameters:
comm – MPI communicator.
points – Coordinates of the ‘lower-left’ and ‘upper-right’ corners of the box.
n – List of cells in each direction
cell_type – The cell type.
dtype – Float type for the mesh geometry(
numpy.float32
ornumpy.float64
).ghost_mode – The ghost mode used in the mesh partitioning.
partitioner – Function that computes the parallel distribution of cells across MPI ranks.
- Returns:
A mesh of a box domain.
- dolfinx.mesh.create_geometry(index_map: IndexMap, dofmap: ndarray[Any, dtype[int32]], element: CoordinateElement, x: ndarray, input_global_indices: ndarray[Any, dtype[int64]]) Geometry [source]
Create a Geometry object.
- Parameters:
index_map – Index map describing the layout of the geometry points (nodes).
dofmap – The geometry (point) dofmap. For a cell, it gives the row in the point coordinates
x
of each local geometry node.shape=(num_cells, num_dofs_per_cell)
.element – Element that describes the cell geometry map.
x – The point coordinates. The shape is
(num_points, geometric_dimension).
input_global_indices – The ‘global’ input index of each point, commonly from a mesh input file.
- dolfinx.mesh.create_interval(comm: ~mpi4py.MPI.Comm, nx: int, points: ~collections.abc.Buffer | ~numpy._typing._array_like._SupportsArray[~numpy.dtype[~typing.Any]] | ~numpy._typing._nested_sequence._NestedSequence[~numpy._typing._array_like._SupportsArray[~numpy.dtype[~typing.Any]]] | bool | int | float | complex | str | bytes | ~numpy._typing._nested_sequence._NestedSequence[bool | int | float | complex | str | bytes], dtype: ~numpy.dtype[~typing.Any] | None | type[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | str | tuple[~typing.Any, int] | tuple[~typing.Any, ~typing.SupportsIndex | ~collections.abc.Sequence[~typing.SupportsIndex]] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | tuple[~typing.Any, ~typing.Any] = <class 'numpy.float64'>, ghost_mode=GhostMode.shared_facet, partitioner=None) Mesh [source]
Create an interval mesh.
- Parameters:
comm – MPI communicator.
nx – Number of cells.
points – Coordinates of the end points.
dtype – Float type for the mesh geometry(
numpy.float32
ornumpy.float64
).ghost_mode – Ghost mode used in the mesh partitioning. Options are
GhostMode.none
andGhostMode.shared_facet
.partitioner – Partitioning function to use for determining the parallel distribution of cells across MPI ranks.
- Returns:
An interval mesh.
- dolfinx.mesh.create_mesh(comm: Comm, cells: ndarray[Any, dtype[int64]], x: ndarray[Any, dtype[floating]], e: Mesh | FiniteElement | _BasixElement | CoordinateElement, partitioner: Callable | None = None) Mesh [source]
Create a mesh from topology and geometry arrays.
- Parameters:
comm – MPI communicator to define the mesh on.
cells – Cells of the mesh.
cells[i]
are the ‘nodes’ of celli
.x – Mesh geometry (‘node’ coordinates), with shape
(num_nodes, gdim)
.e – UFL mesh. The mesh scalar type is determined by the scalar type of
e
.partitioner – Function that determines the parallel distribution of cells across MPI ranks.
Note
If required, the coordinates
x
will be cast to the same scalar type as the domain/elemente
.- Returns:
A mesh.
- dolfinx.mesh.create_rectangle(comm: ~mpi4py.MPI.Comm, points: ~collections.abc.Buffer | ~numpy._typing._array_like._SupportsArray[~numpy.dtype[~typing.Any]] | ~numpy._typing._nested_sequence._NestedSequence[~numpy._typing._array_like._SupportsArray[~numpy.dtype[~typing.Any]]] | bool | int | float | complex | str | bytes | ~numpy._typing._nested_sequence._NestedSequence[bool | int | float | complex | str | bytes], n: ~collections.abc.Buffer | ~numpy._typing._array_like._SupportsArray[~numpy.dtype[~typing.Any]] | ~numpy._typing._nested_sequence._NestedSequence[~numpy._typing._array_like._SupportsArray[~numpy.dtype[~typing.Any]]] | bool | int | float | complex | str | bytes | ~numpy._typing._nested_sequence._NestedSequence[bool | int | float | complex | str | bytes], cell_type=CellType.triangle, dtype: ~numpy.dtype[~typing.Any] | None | type[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | str | tuple[~typing.Any, int] | tuple[~typing.Any, ~typing.SupportsIndex | ~collections.abc.Sequence[~typing.SupportsIndex]] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | tuple[~typing.Any, ~typing.Any] = <class 'numpy.float64'>, ghost_mode=GhostMode.shared_facet, partitioner=None, diagonal: ~dolfinx.cpp.mesh.DiagonalType = DiagonalType.right) Mesh [source]
Create a rectangle mesh.
- Parameters:
comm – MPI communicator.
points – Coordinates of the lower - left and upper - right corners of the rectangle.
n – Number of cells in each direction.
cell_type – Mesh cell type.
dtype – Float type for the mesh geometry(
numpy.float32
ornumpy.float64
)ghost_mode – Ghost mode used in the mesh partitioning.
partitioner – Function that computes the parallel distribution of cells across MPI ranks.
diagonal – Direction of diagonal of triangular meshes. The options are
left
,right
,crossed
,left / right
,right / left
.
- Returns:
A mesh of a rectangle.
- dolfinx.mesh.create_submesh(msh: Mesh, dim: int, entities: ndarray[Any, dtype[int32]]) tuple[Mesh, ndarray[Any, dtype[int32]], ndarray[Any, dtype[int32]], ndarray[Any, dtype[int32]]] [source]
Create a mesh with specified entities from another mesh.
- Parameters:
msh – Mesh to create the sub-mesh from.
dim – Topological dimension of the entities in
msh
to include in the sub-mesh.entities – Indices of entities in
msh
to include in the sub-mesh.
- Returns:
The (1) sub mesh, (2) entity map, (3) vertex map, and (4) node map (geometry). Each of the maps a local index of the sub mesh to a local index of
msh
.
- dolfinx.mesh.create_unit_cube(comm: ~mpi4py.MPI.Comm, nx: int, ny: int, nz: int, cell_type=CellType.tetrahedron, dtype: ~numpy.dtype[~typing.Any] | None | type[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | str | tuple[~typing.Any, int] | tuple[~typing.Any, ~typing.SupportsIndex | ~collections.abc.Sequence[~typing.SupportsIndex]] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | tuple[~typing.Any, ~typing.Any] = <class 'numpy.float64'>, ghost_mode=GhostMode.shared_facet, partitioner=None) Mesh [source]
Create a mesh of a unit cube.
- Parameters:
comm – MPI communicator.
nx – Number of cells in “x” direction.
ny – Number of cells in “y” direction.
nz – Number of cells in “z” direction.
cell_type – Mesh cell type
dtype – Float type for the mesh geometry(
numpy.float32
ornumpy.float64
).ghost_mode – Ghost mode used in the mesh partitioning.
partitioner – Function that computes the parallel distribution of cells across MPI ranks.
- Returns:
A mesh of an axis-aligned unit cube with corners at
(0, 0, 0)
and(1, 1, 1)
.
- dolfinx.mesh.create_unit_interval(comm: ~mpi4py.MPI.Comm, nx: int, dtype: ~numpy.dtype[~typing.Any] | None | type[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | str | tuple[~typing.Any, int] | tuple[~typing.Any, ~typing.SupportsIndex | ~collections.abc.Sequence[~typing.SupportsIndex]] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | tuple[~typing.Any, ~typing.Any] = <class 'numpy.float64'>, ghost_mode=GhostMode.shared_facet, partitioner=None) Mesh [source]
Create a mesh on the unit interval.
- Parameters:
comm – MPI communicator.
nx – Number of cells.
points – Coordinates of the end points.
dtype – Float type for the mesh geometry(
numpy.float32
ornumpy.float64
).ghost_mode – Ghost mode used in the mesh partitioning. Options are
GhostMode.none
andGhostMode.shared_facet
.partitioner – Partitioning function to use for determining the parallel distribution of cells across MPI ranks.
- Returns:
A unit interval mesh with end points at 0 and 1.
- dolfinx.mesh.create_unit_square(comm: ~mpi4py.MPI.Comm, nx: int, ny: int, cell_type=CellType.triangle, dtype: ~numpy.dtype[~typing.Any] | None | type[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | str | tuple[~typing.Any, int] | tuple[~typing.Any, ~typing.SupportsIndex | ~collections.abc.Sequence[~typing.SupportsIndex]] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | tuple[~typing.Any, ~typing.Any] = <class 'numpy.float64'>, ghost_mode=GhostMode.shared_facet, partitioner=None, diagonal: ~dolfinx.cpp.mesh.DiagonalType = DiagonalType.right) Mesh [source]
Create a mesh of a unit square.
- Parameters:
comm – MPI communicator.
nx – Number of cells in the “x” direction.
ny – Number of cells in the “y” direction.
cell_type – Mesh cell type.
dtype – Float type for the mesh geometry(
numpy.float32
ornumpy.float64
).ghost_mode – Ghost mode used in the mesh partitioning.
partitioner – Function that computes the parallel distribution of cells across MPI ranks.
diagonal – Direction of diagonal.
- Returns:
A mesh of a square with corners at (0, 0) and (1, 1).
- dolfinx.mesh.entities_to_geometry(msh: Mesh, dim: int, entities: ndarray[Any, dtype[int32]], permute=False) ndarray[Any, dtype[int32]] [source]
Compute the geometric DOFs associated with the closure of the given mesh entities.
- Parameters:
msh – The mesh.
dim – Topological dimension of the entities of interest.
entities – Entity indices (local to the process).
permute – Permute the DOFs such that they are consistent with the orientation of dim-dimensional mesh entities. This requires create_entity_permutations to be called first.
- Returns:
The geometric DOFs associated with the closure of the entities in entities.
- dolfinx.mesh.exterior_facet_indices(topology: Topology) ndarray[Any, dtype[int32]] [source]
Compute the indices of all exterior facets that are owned by the caller.
An exterior facet (co-dimension 1) is one that is connected globally to only one cell of co-dimension 0).
Note
This is a collective operation that should be called on all processes.
- Parameters:
topology – The topology
- Returns:
Sorted list of owned facet indices that are exterior facets of the mesh.
- dolfinx.mesh.locate_entities(msh: Mesh, dim: int, marker: Callable) ndarray [source]
Compute mesh entities satisfying a geometric marking function.
- Parameters:
msh – Mesh to locate entities on.
dim – Topological dimension of the mesh entities to consider.
marker – A function that takes an array of points
x
with shape(gdim, num_points)
and returns an array of booleans of lengthnum_points
, evaluating to True for entities to be located.
- Returns:
Indices (local to the process) of marked mesh entities.
- dolfinx.mesh.locate_entities_boundary(msh: Mesh, dim: int, marker: Callable) ndarray [source]
Compute mesh entities that are connected to an owned boundary facet and satisfy a geometric marking function.
For vertices and edges, in parallel this function will not necessarily mark all entities that are on the exterior boundary. For example, it is possible for a process to have a vertex that lies on the boundary without any of the attached facets being a boundary facet. When used to find degrees-of-freedom, e.g. using
dolfinx.fem.locate_dofs_topological()
, the function that uses the data returned by this function must typically perform some parallel communication.- Parameters:
msh – Mesh to locate boundary entities on.
dim – Topological dimension of the mesh entities to consider
marker – Function that takes an array of points
x
with shape(gdim, num_points)
and returns an array of booleans of lengthnum_points
, evaluating toTrue
for entities to be located.
- Returns:
Indices (local to the process) of marked mesh entities.
- dolfinx.mesh.meshtags(msh: Mesh, dim: int, entities: ndarray[Any, dtype[int32]], values: ndarray | int | float) MeshTags [source]
Create a MeshTags object that associates data with a subset of mesh entities.
- Parameters:
msh – The mesh.
dim – Topological dimension of the mesh entity.
entities – Indices(local to process) of entities to associate values with . The array must be sorted and must not contain duplicates.
values – The corresponding value for each entity.
- Returns:
A mesh tags object.
Note
The type of the returned MeshTags is inferred from the type of
values
.
- dolfinx.mesh.meshtags_from_entities(msh: Mesh, dim: int, entities: AdjacencyList_int32, values: ndarray[Any, dtype[Any]])[source]
Create a :class:dolfinx.mesh.MeshTags` object that associates data with a subset of mesh entities, where the entities are defined by their vertices.
- Parameters:
msh – The mesh.
dim – Topological dimension of the mesh entity.
entities – Entities to associated values with, with entities defined by their vertices.
values – The corresponding value for each entity.
- Returns:
A mesh tags object.
Note
The type of the returned MeshTags is inferred from the type of
values
.
- dolfinx.mesh.refine(msh: ~dolfinx.mesh.Mesh, edges: ~numpy.ndarray | None = None, partitioner: ~typing.Callable | None = <nanobind.nb_func object>, option: ~dolfinx.cpp.refinement.RefinementOption = RefinementOption.none) tuple[Mesh, ndarray[Any, dtype[int32]], ndarray[Any, dtype[int8]]] [source]
Refine a mesh.
Passing
None
forpartitioner
, refined cells will be on the same process as the parent cell.Note
Using the default partitioner for the refined mesh, the refined mesh will not include ghosts cells (cells connected by facet to an owned cells) even if the parent mesh is ghosted.
Warning
Passing
None
forpartitioner
, the refined mesh will not have ghosts cells even if the parent mesh has ghost cells. The possibility to not re-partition the refined mesh and include ghost cells in the refined mesh will be added in a future release.- Parameters:
msh – Mesh from which to create the refined mesh.
edges – Indices of edges to split during refinement. If
None
, mesh refinement is uniform.partitioner – Partitioner to distribute the refined mesh. If
None
no redistribution is performed, i.e. refined cells remain on the same process as the parent cell.option – Controls whether parent cells and/or parent facets are computed.
- Returns:
Refined mesh, (optional) parent cells, (optional) parent facets
- dolfinx.mesh.transfer_meshtag(meshtag: MeshTags, msh1: Mesh, parent_cell: ndarray[Any, dtype[int32]], parent_facet: ndarray[Any, dtype[int8]] | None = None) MeshTags [source]
Generate cell mesh tags on a refined mesh from the mesh tags on the coarse parent mesh.
- Parameters:
meshtag – Mesh tags on the coarse, parent mesh.
msh1 – The refined mesh.
parent_cell – Index of the parent cell for each cell in the refined mesh.
parent_facet – Index of the local parent facet for each cell in the refined mesh. Only required for transfer tags on facets.
- Returns:
Mesh tags on the refined mesh.