dolfinx.mesh
Creation, refining and marking of meshes
Functions
| 
 | |
| 
 | |
| 
 | Create a box mesh. | 
| 
 | Create an interval mesh. | 
| 
 | Create a mesh from topology and geometry arrays. | 
| 
 | Create a rectangle mesh. | 
| 
 | |
| 
 | Create a mesh of a unit cube. | 
| 
 | Create a mesh on the unit interval. | 
| 
 | Create a mesh of a unit square. | 
| 
 | 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. | 
| 
 | Create a MeshTags object that associates data with a subset of mesh entities. | 
| 
 | Create a :class:dolfinx.mesh.MeshTags` object that associates data with a subset of mesh entities, where the entities are defined by their vertices. | 
| 
 | Refine a mesh. | 
| 
 | Refine a mesh. | 
| 
 | Generate cell mesh tags on a refined mesh from the mesh tags on the coarse parent mesh. | 
Classes
| 
 | A class for representing meshes. | 
| 
 | Mesh tags associate data (markers) with a subset of mesh entities of a given dimension. | 
- class dolfinx.mesh.CellType(self: dolfinx.cpp.mesh.CellType, value: int)
- Bases: - pybind11_object- Members: - point - interval - triangle - quadrilateral - tetrahedron - pyramid - prism - hexahedron - hexahedron = <CellType.hexahedron: -8>
 - interval = <CellType.interval: 2>
 - property name
 - point = <CellType.point: 1>
 - prism = <CellType.prism: -6>
 - pyramid = <CellType.pyramid: -5>
 - quadrilateral = <CellType.quadrilateral: -4>
 - tetrahedron = <CellType.tetrahedron: 4>
 - triangle = <CellType.triangle: 3>
 - property value
 
- class dolfinx.mesh.GhostMode(self: dolfinx.cpp.mesh.GhostMode, value: int)
- Bases: - pybind11_object- Members: - none - shared_facet - shared_vertex - property name
 - none = <GhostMode.none: 0>
 - property value
 
- class dolfinx.mesh.Mesh(mesh, domain: Mesh)[source]
- Bases: - object- A class for representing meshes. - Initialize mesh from a C++ mesh. - Parameters:
- mesh – The C++ mesh object. 
- domain – The UFL domain. 
 
 - Note - Mesh objects should not usually be created using this class directly. - property comm
 - property geometry
- Mesh geometry. 
 - 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 - dimto compute size measure of.
 
- Returns:
- Size measure for each requested entity. 
 
 - property name
 - property topology
- Mesh topology. 
 
- 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 + + - meshtagsobject. If mesh is passed,- meshand- meshtagsmust 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. 
 
- dolfinx.mesh.build_dual_graph(comm: MPICommWrapper, cells: dolfinx::graph::AdjacencyList<long>, tdim: int) dolfinx::graph::AdjacencyList<long>
- Build dual graph for cells 
- dolfinx.mesh.cell_dim(type: dolfinx.cpp.mesh.CellType) int
- dolfinx.mesh.compute_incident_entities(topology, entities: ndarray[Any, dtype[int32]], d0: int, d1: int)[source]
- dolfinx.mesh.create_box(comm: ~mpi4py.MPI.Comm, points: ~typing.List[~typing.Union[~numpy._typing._array_like._SupportsArray[~numpy.dtype], ~numpy._typing._nested_sequence._NestedSequence[~numpy._typing._array_like._SupportsArray[~numpy.dtype]], bool, int, float, complex, str, bytes, ~numpy._typing._nested_sequence._NestedSequence[~typing.Union[bool, int, float, complex, str, bytes]]]], n: list, cell_type=<CellType.tetrahedron: 4>, dtype: ~typing.Union[~numpy.dtype[~typing.Any], None, ~typing.Type[~typing.Any], ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]], str, ~typing.Tuple[~typing.Any, int], ~typing.Tuple[~typing.Any, ~typing.Union[~typing.SupportsIndex, ~typing.Sequence[~typing.SupportsIndex]]], ~typing.List[~typing.Any], ~numpy._typing._dtype_like._DTypeDict, ~typing.Tuple[~typing.Any, ~typing.Any]] = <class 'numpy.float64'>, ghost_mode=<GhostMode.shared_facet: 1>, 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.float32or- numpy.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_cell_partitioner(*args, **kwargs)
- Overloaded function. - create_cell_partitioner(arg0: dolfinx.cpp.mesh.GhostMode) -> Callable[[MPICommWrapper, int, int, dolfinx::graph::AdjacencyList<long>], dolfinx::graph::AdjacencyList<int>] 
 - Create default cell partitioner. - create_cell_partitioner(part: Callable[[MPICommWrapper, int, dolfinx::graph::AdjacencyList<long>, bool], dolfinx::graph::AdjacencyList<int>], ghost_mode: dolfinx.cpp.mesh.GhostMode = <GhostMode.none: 0>) -> Callable[[MPICommWrapper, int, int, dolfinx::graph::AdjacencyList<long>], dolfinx::graph::AdjacencyList<int>] 
 - Create a cell partitioner from a graph partitioning function. 
- dolfinx.mesh.create_interval(comm: ~mpi4py.MPI.Comm, nx: int, points: ~typing.Union[~numpy._typing._array_like._SupportsArray[~numpy.dtype], ~numpy._typing._nested_sequence._NestedSequence[~numpy._typing._array_like._SupportsArray[~numpy.dtype]], bool, int, float, complex, str, bytes, ~numpy._typing._nested_sequence._NestedSequence[~typing.Union[bool, int, float, complex, str, bytes]]], dtype: ~typing.Union[~numpy.dtype[~typing.Any], None, ~typing.Type[~typing.Any], ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]], str, ~typing.Tuple[~typing.Any, int], ~typing.Tuple[~typing.Any, ~typing.Union[~typing.SupportsIndex, ~typing.Sequence[~typing.SupportsIndex]]], ~typing.List[~typing.Any], ~numpy._typing._dtype_like._DTypeDict, ~typing.Tuple[~typing.Any, ~typing.Any]] = <class 'numpy.float64'>, ghost_mode=<GhostMode.shared_facet: 1>, 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.float32or- numpy.float64).
- ghost_mode – Ghost mode used in the mesh partitioning. Options are - GhostMode.noneand- GhostMode.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: Union[ndarray, AdjacencyList_int64], x: ndarray, domain: Mesh, partitioner=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]is the ‘nodes’ of cell- i.
- x – Mesh geometry (‘node’ coordinates), with shape - (num_nodes, gdim).
- domain – UFL mesh. 
- partitioner – Function that computes the parallel distribution of cells across MPI ranks. 
 
- Returns:
- A mesh. 
 
- dolfinx.mesh.create_rectangle(comm: ~mpi4py.MPI.Comm, points: ~typing.Union[~numpy._typing._array_like._SupportsArray[~numpy.dtype], ~numpy._typing._nested_sequence._NestedSequence[~numpy._typing._array_like._SupportsArray[~numpy.dtype]], bool, int, float, complex, str, bytes, ~numpy._typing._nested_sequence._NestedSequence[~typing.Union[bool, int, float, complex, str, bytes]]], n: ~typing.Union[~numpy._typing._array_like._SupportsArray[~numpy.dtype], ~numpy._typing._nested_sequence._NestedSequence[~numpy._typing._array_like._SupportsArray[~numpy.dtype]], bool, int, float, complex, str, bytes, ~numpy._typing._nested_sequence._NestedSequence[~typing.Union[bool, int, float, complex, str, bytes]]], cell_type=<CellType.triangle: 3>, dtype: ~typing.Union[~numpy.dtype[~typing.Any], None, ~typing.Type[~typing.Any], ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]], str, ~typing.Tuple[~typing.Any, int], ~typing.Tuple[~typing.Any, ~typing.Union[~typing.SupportsIndex, ~typing.Sequence[~typing.SupportsIndex]]], ~typing.List[~typing.Any], ~numpy._typing._dtype_like._DTypeDict, ~typing.Tuple[~typing.Any, ~typing.Any]] = <class 'numpy.float64'>, ghost_mode=<GhostMode.shared_facet: 1>, partitioner=None, diagonal: ~dolfinx.cpp.mesh.DiagonalType = <DiagonalType.right: 1>) 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.float32or- numpy.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_unit_cube(comm: ~mpi4py.MPI.Comm, nx: int, ny: int, nz: int, cell_type=<CellType.tetrahedron: 4>, dtype: ~typing.Union[~numpy.dtype[~typing.Any], None, ~typing.Type[~typing.Any], ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]], str, ~typing.Tuple[~typing.Any, int], ~typing.Tuple[~typing.Any, ~typing.Union[~typing.SupportsIndex, ~typing.Sequence[~typing.SupportsIndex]]], ~typing.List[~typing.Any], ~numpy._typing._dtype_like._DTypeDict, ~typing.Tuple[~typing.Any, ~typing.Any]] = <class 'numpy.float64'>, ghost_mode=<GhostMode.shared_facet: 1>, 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.float32or- numpy.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).
 
- A mesh of an axis-aligned unit cube with corners at 
 
- dolfinx.mesh.create_unit_interval(comm: ~mpi4py.MPI.Comm, nx: int, dtype: ~typing.Union[~numpy.dtype[~typing.Any], None, ~typing.Type[~typing.Any], ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]], str, ~typing.Tuple[~typing.Any, int], ~typing.Tuple[~typing.Any, ~typing.Union[~typing.SupportsIndex, ~typing.Sequence[~typing.SupportsIndex]]], ~typing.List[~typing.Any], ~numpy._typing._dtype_like._DTypeDict, ~typing.Tuple[~typing.Any, ~typing.Any]] = <class 'numpy.float64'>, ghost_mode=<GhostMode.shared_facet: 1>, 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.float32or- numpy.float64).
- ghost_mode – Ghost mode used in the mesh partitioning. Options are - GhostMode.noneand- GhostMode.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: 3>, dtype: ~typing.Union[~numpy.dtype[~typing.Any], None, ~typing.Type[~typing.Any], ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]], str, ~typing.Tuple[~typing.Any, int], ~typing.Tuple[~typing.Any, ~typing.Union[~typing.SupportsIndex, ~typing.Sequence[~typing.SupportsIndex]]], ~typing.List[~typing.Any], ~numpy._typing._dtype_like._DTypeDict, ~typing.Tuple[~typing.Any, ~typing.Any]] = <class 'numpy.float64'>, ghost_mode=<GhostMode.shared_facet: 1>, partitioner=None, diagonal: ~dolfinx.cpp.mesh.DiagonalType = <DiagonalType.right: 1>) 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.float32or- numpy.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.exterior_facet_indices(topology: dolfinx.cpp.mesh.Topology) numpy.ndarray[numpy.int32]
- dolfinx.mesh.locate_entities(mesh: Mesh, dim: int, marker: Callable) ndarray[source]
- Compute mesh entities satisfying a geometric marking function. - Parameters:
- mesh – Mesh to locate entities on. 
- dim – Topological dimension of the mesh entities to consider. 
- marker – A function that takes an array of points - xwith shape- (gdim, num_points)and returns an array of booleans of length- num_points, evaluating to True for entities to be located.
 
- Returns:
- Indices (local to the process) of marked mesh entities. 
 
- dolfinx.mesh.locate_entities_boundary(mesh: 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:
- mesh – Mesh to locate boundary entities on. 
- dim – Topological dimension of the mesh entities to consider 
- marker – Function that takes an array of points - xwith shape- (gdim, num_points)and returns an array of booleans of length- num_points, evaluating to- Truefor entities to be located.
 
- Returns:
- Indices (local to the process) of marked mesh entities. 
 
- dolfinx.mesh.meshtags(mesh: Mesh, dim: int, entities: ndarray[Any, dtype[int32]], values: Union[ndarray, int, float]) MeshTags[source]
- Create a MeshTags object that associates data with a subset of mesh entities. - Parameters:
- mesh – 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(mesh: 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:
- mesh – 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(mesh: Mesh, edges: Optional[ndarray] = None, redistribute: bool = True) Mesh[source]
- Refine a mesh. - Parameters:
- mesh – Mesh from which to create the refined mesh. 
- edges – Indices of edges to split during refinement. If - None, mesh refinement is uniform.
- redistribute – Refined mesh is re-partitioned if - True.
 
- Returns:
- Refined mesh. 
 
- dolfinx.mesh.refine_plaza(mesh: ~dolfinx.mesh.Mesh, edges: ~typing.Optional[~numpy.ndarray] = None, redistribute: bool = True, option: ~dolfinx.cpp.refinement.RefinementOption = <RefinementOption.none: 0>) tuple[dolfinx.mesh.Mesh, numpy.ndarray[Any, numpy.dtype[numpy.int32]], numpy.ndarray[Any, numpy.dtype[numpy.int32]]][source]
- Refine a mesh. - Parameters:
- mesh – Mesh from which to create the refined mesh. 
- edges – Indices of edges to split during refinement. If - None, mesh refinement is uniform.
- redistribute – Refined mesh is re-partitioned if - True.
- option – Control computation of the parent-refined mesh data. 
 
- Returns:
- Refined mesh, list of parent cell for each refine cell, and list of parent facets. 
 
- dolfinx.mesh.to_string(type: dolfinx.cpp.mesh.CellType) str
- dolfinx.mesh.to_type(cell: str) dolfinx.cpp.mesh.CellType
- dolfinx.mesh.transfer_meshtag(meshtag: MeshTags, mesh1: Mesh, parent_cell: ndarray[Any, dtype[int32]], parent_facet: Optional[ndarray[Any, dtype[int8]]] = 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. 
- mesh1 – 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.