dolfinx.mesh

Creation, refining and marking of meshes

Functions

create_box(comm, points, n[, cell_type, ...])

Create box mesh

create_interval(comm, nx, points[, ...])

Create an interval mesh

create_mesh(comm, cells, x, domain[, ...])

Create a mesh from topology and geometry arrays

create_rectangle(comm, points, n[, ...])

Create rectangle mesh

create_submesh(mesh, dim, entities)

create_unit_cube(comm, nx, ny, nz[, ...])

Create a mesh of a unit cube

create_unit_interval(comm, nx[, ghost_mode, ...])

Create a mesh on the unit interval

create_unit_square(comm, nx, ny[, ...])

Create a mesh of a unit square

locate_entities(mesh, dim, marker)

Compute mesh entities satisfying a geometric marking function

locate_entities_boundary(mesh, dim, marker)

Compute mesh entities that are connected to an owned boundary facet and satisfy a geometric marking function

meshtags(mesh, dim, indices, values)

Create a MeshTags object that associates data with a subset of mesh entities.

meshtags_from_entities(mesh, dim, entities, ...)

Create a MeshTags object that associates data with a subset of mesh entities, where the entities are defined by their vertices.

refine(mesh[, edges, redistribute])

Refine a mesh

Classes

Mesh(comm, topology, geometry, domain)

A class for representing meshes

MeshTagsMetaClass(mesh, dim, indices, values)

A distributed sparse matrix that uses compressed sparse row storage.

class dolfinx.mesh.CellType(self: dolfinx.cpp.mesh.CellType, value: int) None

Bases: pybind11_builtins.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) None

Bases: pybind11_builtins.pybind11_object

Members:

none

shared_facet

shared_vertex

property name
none = <GhostMode.none: 0>
shared_facet = <GhostMode.shared_facet: 1>
shared_vertex = <GhostMode.shared_vertex: 2>
property value
class dolfinx.mesh.Mesh(comm: mpi4py.MPI.Comm, topology: dolfinx.cpp.mesh.Topology, geometry: dolfinx.cpp.mesh.Geometry, domain: ufl.domain.Mesh)[source]

Bases: dolfinx.cpp.mesh.Mesh

A class for representing meshes

Parameters
  • comm – The MPI communicator

  • topology – The mesh topology

  • geometry – The mesh geometry

  • domain – The MPI communicator

Note

Mesh objects are not generally created using this class directly.

classmethod from_cpp(obj: dolfinx.cpp.mesh.Mesh, domain: ufl.domain.Mesh) dolfinx.mesh.Mesh[source]

Create Mesh object from a C++ Mesh object

ufl_cell() ufl.cell.Cell[source]

Return the UFL cell type

ufl_domain() ufl.domain.Mesh[source]

Return the ufl domain corresponding to the mesh.

class dolfinx.mesh.MeshTagsMetaClass(mesh: dolfinx.mesh.Mesh, dim: int, indices: numpy.ndarray, values: numpy.ndarray)[source]

Bases: object

A distributed sparse matrix that uses compressed sparse row storage.

Parameters
  • mesh – The mesh

  • dim – Topological dimension of the mesh entity

  • indices – Entity indices (local to process)

  • values – The corresponding value for each entity

Note

Objects of this type should be created using meshtags() and not created using this initialiser directly.

ufl_id() int[source]

Object identifier.

Notes

This method is used by UFL.

Returns

The id of the object

dolfinx.mesh.build_dual_graph(arg0: MPICommWrapper, arg1: dolfinx::graph::AdjacencyList<long>, arg2: int) dolfinx::graph::AdjacencyList<long>

Build dual graph for cells

dolfinx.mesh.cell_dim(arg0: dolfinx.cpp.mesh.CellType) int
dolfinx.mesh.compute_boundary_facets(arg0: dolfinx::mesh::Topology) List[int]
dolfinx.mesh.compute_incident_entities(arg0: dolfinx.cpp.mesh.Mesh, arg1: numpy.ndarray[numpy.int32], arg2: int, arg3: int) numpy.ndarray[numpy.int32]
dolfinx.mesh.compute_midpoints(arg0: dolfinx::mesh::Mesh, arg1: int, arg2: numpy.ndarray[numpy.int32]) numpy.ndarray[numpy.float64]
dolfinx.mesh.create_box(comm: mpi4py.MPI.Comm, points: typing.List[numpy.array], n: list, cell_type=<CellType.tetrahedron: 4>, ghost_mode=<GhostMode.shared_facet: 1>, partitioner=<built-in method  of PyCapsule object>) dolfinx.mesh.Mesh[source]

Create 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

  • 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.

  1. create_cell_partitioner() -> Callable[[MPICommWrapper, int, int, dolfinx::graph::AdjacencyList<long>, dolfinx.cpp.mesh.GhostMode], dolfinx::graph::AdjacencyList<int>]

  2. create_cell_partitioner(arg0: Callable[[MPICommWrapper, int, dolfinx::graph::AdjacencyList<long>, bool], dolfinx::graph::AdjacencyList<int>]) -> Callable[[MPICommWrapper, int, int, dolfinx::graph::AdjacencyList<long>, dolfinx.cpp.mesh.GhostMode], dolfinx::graph::AdjacencyList<int>]

dolfinx.mesh.create_interval(comm: mpi4py.MPI.Comm, nx: int, points: list, ghost_mode=<GhostMode.shared_facet: 1>, partitioner=<built-in method  of PyCapsule object>) dolfinx.mesh.Mesh[source]

Create an interval mesh

Parameters
  • comm – MPI communicator

  • nx – Number of cells

  • points – Coordinates of the end points

  • ghost_mode – Ghost mode used in the mesh partitioning. Options are GhostMode.none’ and `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: mpi4py.MPI.Comm, cells: typing.Union[numpy.ndarray, dolfinx.cpp.graph.AdjacencyList_int64], x: numpy.ndarray, domain: ufl.domain.Mesh, ghost_mode=<GhostMode.shared_facet: 1>, partitioner=<built-in method  of PyCapsule object>) dolfinx.mesh.Mesh[source]

Create a mesh from topology and geometry arrays

Parameters
  • comm – MPI communicator to define the mesh on

  • cells – Cells of the mesh

  • x – Mesh geometry (‘node’ coordinates), with shape (gdim, num_nodes)

  • domain – UFL mesh

  • ghost_mode – The ghost mode used in the mesh partitioning

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

Returns

A new mesh

dolfinx.mesh.create_rectangle(comm: mpi4py.MPI.Comm, points: typing.List[numpy.array], n: list, cell_type=<CellType.triangle: 3>, ghost_mode=<GhostMode.shared_facet: 1>, partitioner=<built-in method  of PyCapsule object>, diagonal: dolfinx.cpp.mesh.DiagonalType = <DiagonalType.right: 1>) dolfinx.mesh.Mesh[source]

Create 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

  • 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>, ghost_mode=<GhostMode.shared_facet: 1>, partitioner=<built-in method  of PyCapsule object>) dolfinx.mesh.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

  • 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, ghost_mode=<GhostMode.shared_facet: 1>, partitioner=<built-in method  of PyCapsule object>) dolfinx.mesh.Mesh[source]

Create a mesh on the unit interval

Parameters
  • comm – MPI communicator

  • nx – Number of cells

  • points – Coordinates of the end points

  • ghost_mode – Ghost mode used in the mesh partitioning. Options are GhostMode.none’ and `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>, ghost_mode=<GhostMode.shared_facet: 1>, partitioner=<built-in method  of PyCapsule object>, diagonal: dolfinx.cpp.mesh.DiagonalType = <DiagonalType.right: 1>) dolfinx.mesh.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

  • 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.locate_entities(mesh: dolfinx.mesh.Mesh, dim: int, marker: function) numpy.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 x with 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: dolfinx.mesh.Mesh, dim: int, marker: function) numpy.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 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 x with 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.meshtags(mesh: dolfinx.mesh.Mesh, dim: int, indices: numpy.ndarray, values: numpy.ndarray) dolfinx.mesh.MeshTagsMetaClass[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

  • indices – Entity indices (local to process)

  • values – The corresponding value for each entity

Returns

A MeshTags object

Note

The type of the returned MeshTags is inferred from the type of values.

dolfinx.mesh.meshtags_from_entities(mesh: dolfinx.mesh.Mesh, dim: int, entities: dolfinx.cpp.graph.AdjacencyList_int32, values: numpy.ndarray)[source]

Create a 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 – Tagged entities, with entities defined by their vertices

  • values – The corresponding value for each entity

Returns

A MeshTags object

Note

The type of the returned MeshTags is inferred from the type of values.

dolfinx.mesh.refine(mesh: dolfinx.mesh.Mesh, edges: Optional[numpy.ndarray] = None, redistribute: bool = True) dolfinx.mesh.Mesh[source]

Refine a mesh

Parameters
  • mesh – The mesh from which to build a refined mesh

  • edges – Optional argument to specify which edges should be refined. If not supplied uniform refinement is applied.

  • redistribute – Optional argument to redistribute the refined mesh if mesh is a distributed mesh.

Returns

A refined mesh

dolfinx.mesh.to_string(arg0: dolfinx.cpp.mesh.CellType) str
dolfinx.mesh.to_type(arg0: str) dolfinx.cpp.mesh.CellType