dolfinx.mesh

Creation, refining and marking of meshes

Functions

MeshTags(mesh, dim, indices, values)

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

Create a mesh from topology and geometry data

locate_entities(mesh, dim, marker)

Compute list of mesh entities satisfying a geometric marking function.

locate_entities_boundary(mesh, dim, marker)

Compute list of mesh entities that are attached to an owned boundary facet and satisfy a geometric marking function.

refine(mesh[, cell_markers, redistribute])

Refine a mesh

dolfinx.mesh.MeshTags(mesh, dim, indices, values)[source]
dolfinx.mesh.create_mesh(comm, cells, x, domain, ghost_mode=<GhostMode.shared_facet: 1>, partitioner=<built-in method partition_cells_graph of PyCapsule object>)[source]

Create a mesh from topology and geometry data

dolfinx.mesh.create_meshtags(*args, **kwargs)

Overloaded function.

  1. create_meshtags(arg0: dolfinx.cpp.mesh.Mesh, arg1: int, arg2: dolfinx::graph::AdjacencyList<int>, arg3: numpy.ndarray[numpy.int8]) -> dolfinx.cpp.mesh.MeshTags_int8

  2. create_meshtags(arg0: dolfinx.cpp.mesh.Mesh, arg1: int, arg2: dolfinx::graph::AdjacencyList<int>, arg3: numpy.ndarray[numpy.int32]) -> dolfinx.cpp.mesh.MeshTags_int32

  3. create_meshtags(arg0: dolfinx.cpp.mesh.Mesh, arg1: int, arg2: dolfinx::graph::AdjacencyList<int>, arg3: numpy.ndarray[numpy.float64]) -> dolfinx.cpp.mesh.MeshTags_double

  4. create_meshtags(arg0: dolfinx.cpp.mesh.Mesh, arg1: int, arg2: dolfinx::graph::AdjacencyList<int>, arg3: numpy.ndarray[numpy.int64]) -> dolfinx.cpp.mesh.MeshTags_int64

dolfinx.mesh.locate_entities(mesh: dolfinx.cpp.mesh.Mesh, dim: int, marker: function)[source]

Compute list of mesh entities satisfying a geometric marking function.

Parameters
  • mesh – The mesh

  • dim – The 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.

Return type

numpy.ndarray

dolfinx.mesh.locate_entities_boundary(mesh: dolfinx.cpp.mesh.Mesh, dim: int, marker: function)[source]

Compute list of mesh entities that are attached 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 – The mesh

  • dim – The 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.

Return type

numpy.ndarray

dolfinx.mesh.refine(mesh, cell_markers=None, redistribute=True)[source]

Refine a mesh