Mesh (dolfinx::mesh)

Under development

namespace mesh

Mesh data structures and algorithms on meshes.

Representations of meshes and support for operations on meshes.

Typedefs

using CellPartitionFunction = std::function<graph::AdjacencyList<std::int32_t>(MPI_Comm comm, int nparts, const std::vector<CellType> &cell_types, const std::vector<std::span<const std::int64_t>> &cells)>

Signature for the cell partitioning function. Function that implement this interface compute the destination rank for cells currently on this rank.

Note

Cells can have multiple destination ranks, when ghosted.

Param comm:

[in] MPI Communicator.

Param nparts:

[in] Number of partitions.

Param cell_types:

[in] Cell types in the mesh.

Param cells:

[in] Lists of cells of each cell type. cells[i] is a flattened row major 2D array of shape (num_cells, num_cell_vertices) for cell_types[i] on this process, containing the global indices for the cell vertices. Each cell can appear only once across all processes. The cell vertex indices are not necessarily contiguous globally, i.e. the maximum index across all processes can be greater than the number of vertices. High-order ‘nodes’, e.g. mid-side points, should not be included.

Return:

Destination ranks for each cell on this process.

Enums

enum class CellType : int

Cell type identifier.

Values:

enumerator point
enumerator interval
enumerator triangle
enumerator tetrahedron
enumerator quadrilateral
enumerator pyramid
enumerator prism
enumerator hexahedron
enum class DiagonalType

Enum for different diagonal types.

Values:

enumerator left
enumerator right
enumerator crossed
enumerator shared_facet
enumerator left_right
enumerator right_left
enum class GhostMode : int

Enum for different partitioning ghost modes.

Values:

enumerator none
enumerator shared_facet
enumerator shared_vertex

Functions

std::string to_string(CellType type)

Get the cell string type for a cell type

Parameters:

type[in] The cell type

Returns:

The cell type string

CellType to_type(const std::string &cell)

Get the cell type from a cell string

Parameters:

cell[in] Cell shape string

Returns:

The cell type

CellType cell_entity_type(CellType type, int d, int index)

Return type of cell for entity of dimension d at given entity index.

CellType cell_facet_type(CellType type, int index)

Return facet type of cell For simplex and hypercube cell types, this is independent of the facet index, but for prism and pyramid, it can be triangle or quadrilateral.

Parameters:
  • type[in] The cell type

  • index[in] The facet index

Returns:

The type of facet for this cell at this index

graph::AdjacencyList<int> get_entity_vertices(CellType type, int dim)

Return list of entities, where entities(e, k) is the local vertex index for the kth vertex of entity e of dimension dim

graph::AdjacencyList<int> get_sub_entities(CellType type, int dim0, int dim1)

Get entities of dimension dim1 and that make up entities of dimension dim0

int cell_dim(CellType type)

Return topological dimension of cell type.

int cell_num_entities(CellType type, int dim)

Number of entities of dimension dim

Parameters:
  • dim[in] Entity dimension

  • type[in] Cell type

Returns:

Number of entities in cell

bool is_simplex(CellType type)

Check if cell is a simplex

Parameters:

type[in] Cell type

Returns:

True is the cell type is a simplex

int num_cell_vertices(CellType type)

Number vertices for a cell type

Parameters:

type[in] Cell type

Returns:

The number of cell vertices

std::map<std::array<int, 2>, std::vector<std::set<int>>> cell_entity_closure(CellType cell_type)

Closure entities for a cell, i.e., all lower-dimensional entities attached to a cell entity. Map from entity {dim_e, entity_e} to closure{sub_dim, (sub_entities)}

basix::cell::type cell_type_to_basix_type(CellType celltype)

Convert a cell type to a Basix cell type.

CellType cell_type_from_basix_type(basix::cell::type celltype)

Get a cell type from a Basix cell type.

template<std::floating_point T = double>
Mesh<T> create_box(MPI_Comm comm, MPI_Comm subcomm, std::array<std::array<T, 3>, 2> p, std::array<std::int64_t, 3> n, CellType celltype, CellPartitionFunction partitioner = nullptr)

Create a uniform mesh::Mesh over rectangular prism spanned by the two points p.

The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + / 1)*(n[1] + 1)*(n[2] + 1). For tetrahedra there will be will be 6*n[0]*n[1]*n[2] cells. For hexahedra the number of cells will be n[0]*n[1]*n[2].

Parameters:
  • comm[in] MPI communicator to distribute the mesh on.

  • subcomm[in] MPI communicator to construct and partition the mesh topology on. If the process should not be involved in the topology creation and partitioning then this communicator should be MPI_COMM_NULL.

  • p[in] Corner of the box.

  • n[in] Number of cells in each direction.

  • celltype[in] Cell shape.

  • partitioner[in] Partitioning function for distributing cells across MPI ranks.

Returns:

Mesh

template<std::floating_point T = double>
Mesh<T> create_box(MPI_Comm comm, std::array<std::array<T, 3>, 2> p, std::array<std::int64_t, 3> n, CellType celltype, const CellPartitionFunction &partitioner = nullptr)

Create a uniform mesh::Mesh over rectangular prism spanned by the two points p.

The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + / 1)*(n[1] + 1)*(n[2] + 1). For tetrahedra there will be will be 6*n[0]*n[1]*n[2] cells. For hexahedra the number of cells will be n[0]*n[1]*n[2].

Parameters:
  • comm[in] MPI communicator to distribute the mesh on.

  • p[in] Corner of the box.

  • n[in] Number of cells in each direction.

  • celltype[in] Cell shape.

  • partitioner[in] Partitioning function for distributing cells across MPI ranks.

Returns:

Mesh

template<std::floating_point T = double>
Mesh<T> create_rectangle(MPI_Comm comm, std::array<std::array<T, 2>, 2> p, std::array<std::int64_t, 2> n, CellType celltype, CellPartitionFunction partitioner, DiagonalType diagonal = DiagonalType::right)

Create a uniform mesh::Mesh over the rectangle spanned by the two points p.

The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + / 1)*(n[1] + 1). For triangles there will be will be 2*n[0]*n[1] cells. For quadrilaterals the number of cells will be n[0]*n[1].

Parameters:
  • comm[in] MPI communicator to build the mesh on.

  • p[in] Bottom-left and top-right corners of the rectangle.

  • n[in] Number of cells in each direction.

  • celltype[in] Cell shape.

  • partitioner[in] Partitioning function for distributing cells across MPI ranks.

  • diagonal[in] Direction of diagonals

Returns:

Mesh

template<std::floating_point T = double>
Mesh<T> create_rectangle(MPI_Comm comm, std::array<std::array<T, 2>, 2> p, std::array<std::int64_t, 2> n, CellType celltype, DiagonalType diagonal = DiagonalType::right)

Create a uniform mesh::Mesh over the rectangle spanned by the two points p.

The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + / 1)*(n[1] + 1). For triangles there will be will be 2*n[0]*n[1] cells. For quadrilaterals the number of cells will be n[0]*n[1].

Parameters:
  • comm[in] MPI communicator to build the mesh on

  • p[in] Two corner points

  • n[in] Number of cells in each direction

  • celltype[in] Cell shape

  • diagonal[in] Direction of diagonals

Returns:

Mesh

template<std::floating_point T = double>
Mesh<T> create_interval(MPI_Comm comm, std::int64_t n, std::array<T, 2> p, mesh::GhostMode ghost_mode = mesh::GhostMode::none, CellPartitionFunction partitioner = nullptr)

Interval mesh of the 1D line [a, b].

Given n cells in the axial direction, the total number of intervals will be n and the total number of vertices will be n + 1.

Parameters:
  • comm[in] MPI communicator to build the mesh on.

  • n[in] Number of cells.

  • p[in] End points of the interval.

  • ghost_mode[in] ghost mode of the created mesh, defaults to none

  • partitioner[in] Partitioning function for distributing cells across MPI ranks.

Returns:

A mesh.

template<typename U>
Geometry<typename std::remove_reference_t<typename U::value_type>> create_geometry(const Topology &topology, const std::vector<fem::CoordinateElement<std::remove_reference_t<typename U::value_type>>> &elements, std::span<const std::int64_t> nodes, std::span<const std::int64_t> xdofs, const U &x, int dim, std::function<std::vector<int>(const graph::AdjacencyList<std::int32_t>&)> reorder_fn = nullptr)

Build Geometry from input data.

This function should be called after the mesh topology is built and ‘node’ coordinate data has been distributed to the processes where it is required.

Note

Experimental new interface for multiple cmap/dofmap

Parameters:
  • topology[in] Mesh topology.

  • elements[in] List of elements that defines the geometry map for each cell type.

  • nodes[in] Geometry node global indices for cells on this process.

  • xdofs[in] Geometry degree-of-freedom map (using global indices) for cells on this process. nodes is a sorted and unique list of the indices in xdofs.

  • x[in] The node coordinates (row-major, with shape (num_nodes, dim). The global index of each node is i + / rank_offset, where i is the local row index in x and rank_offset is the sum of x rows on all processed with a lower rank than the caller.

  • dim[in] Geometric dimension (1, 2, or 3).

  • reorder_fn[in] Function for re-ordering the degree-of-freedom map associated with the geometry data.

Pre:

Must be sorted.

Returns:

A mesh geometry.

std::tuple<graph::AdjacencyList<std::int32_t>, std::vector<std::int64_t>, std::size_t, std::vector<std::int32_t>> build_local_dual_graph(std::span<const CellType> celltypes, const std::vector<std::span<const std::int64_t>> &cells)

Compute the local part of the dual graph (cell-cell connections via facets) and facets with only one attached cell.

Each row of the returned data (2) contains [v0, ... v_(n-1), x, .., / x], where v_i is a vertex global index, x is a negative value (all padding values will be equal). The vertex global indices are sorted for each facet.

Note

The cells of each cell type are numbered locally consecutively, i.e. if there are n cells of type 0 and m cells of type 1, then cells of type 0 are numbered 0..(n-1) and cells of type 1 are numbered n..(n+m-1) respectively, in the returned dual graph.

Parameters:
  • celltypes[in] List of cell types.

  • cells[in] Lists of cell vertices (stored as flattened lists, one for each cell type).

Returns:

  1. Local dual graph

  2. Facets, defined by their sorted vertices, that are shared by only one cell on this rank. The logically 2D array is flattened (row-major).

  3. Facet data array (2) number of columns

  4. Attached cell (local index) to each returned facet in (2).

graph::AdjacencyList<std::int64_t> build_dual_graph(MPI_Comm comm, std::span<const CellType> celltypes, const std::vector<std::span<const std::int64_t>> &cells)

Build distributed mesh dual graph (cell-cell connections via facets) from minimal mesh data.

The computed dual graph is typically passed to a graph partitioner.

Note

Collective function

Note

cells and celltypes must have the same size.

Parameters:
  • comm[in] The MPI communicator

  • celltypes[in] List of cell types

  • cells[in] Collections of cells, defined by the cell vertices from which to build the dual graph, as flattened arrays for each cell type in celltypes.

Returns:

The dual graph

template<typename T>
MeshTags<T> create_meshtags(std::shared_ptr<const Topology> topology, int dim, const graph::AdjacencyList<std::int32_t> &entities, std::span<const T> values)

Create MeshTags from arrays.

Note

Entities that do not exist on this rank are ignored.

Warning

entities must not contain duplicate entities.

Parameters:
  • topology[in] Mesh topology that the tags are associated with.

  • dim[in] Topological dimension of tagged entities.

  • entities[in] Local vertex indices for tagged entities.

  • values[in] Tag values for each entity in entities. The length of values must be equal to number of rows in entities.

std::pair<std::vector<std::uint8_t>, std::vector<std::uint32_t>> compute_entity_permutations(const Topology &topology)

Compute (1) facet rotation and reflection data, and (2) cell permutation data. This information is used in the assembly of (1) facet integrals and (2) most non-Lagrange elements.

  1. The facet rotation and reflection data is encoded so that:

    • n % 2 gives the number of reflections to apply

    • n // 2 gives the number of rotations to apply

    The data is stored in a flattened 2D array, so that data[cell_index * / facets_per_cell + facet_index] contains the facet with index facet_index of the cell with index cell_index. This data passed to FFCx kernels, where it is used to permute the quadrature points on facet integrals when data from the cells on both sides of the facet is used.

  2. The cell permutation data contains information about the entities of each cell, relative to a low-to-high ordering. This data is packed so that a 32-bit int is used for each cell. For 2D cells, one bit is used for each edge, to represent whether or not the edge is reversed: the least significant bit is for edge 0, the next for edge 1, etc. For 3D cells, three bits are used for each face, and for each edge: the least significant bit says whether or not face 0 is reflected, the next 2 bits say how many times face 0 is rotated; the next three bits are for face 1, then three for face 2, etc; after all the faces, there is 1 bit for each edge to say whether or not they are reversed.

    For example, if a quadrilateral has cell permutation info ....0111 then (from right to left):

    • edge 0 is reflected (1)

    • edge 1 is reflected (1)

    • edge 2 is reflected (1)

    • edge 3 is not permuted (0)

    and if a tetrahedron has cell permutation info ....011010010101001000 then (from right to left):

    • face 0 is not permuted (000)

    • face 1 is reflected (001)

    • face 2 is rotated twice then reflected (101)

    • face 3 is rotated once (010)

    • edge 0 is not permuted (0)

    • edge 1 is reflected (1)

    • edge 2 is not permuted (0)

    • edge 3 is reflected (1)

    • edge 4 is reflected (1)

    • edge 5 is not permuted (0)

    This data is used to correct the direction of vector function on permuted facets.

Returns:

Facet permutation and cells permutations

Topology create_topology(MPI_Comm comm, const std::vector<CellType> &cell_types, std::vector<std::span<const std::int64_t>> cells, std::vector<std::span<const std::int64_t>> original_cell_index, std::vector<std::span<const int>> ghost_owners, std::span<const std::int64_t> boundary_vertices)

Create a mesh topology.

This function creates a Topology from cells that have been already distributed to the processes that own or ghost the cell.

Parameters:
  • comm[in] Communicator across which the topology will be distributed.

  • cell_types[in] List of cell types in the topology.

  • cells[in] Cell topology (list of vertices for each cell) for each cell type using global indices for the vertices. The cell type for cells[i] is cell_types[i], using row-major storage and where the row cells[i][j] is the vertices for cell j of cell type i. Each cells[i] contains cells that have been distributed to this rank, e.g. via a graph partitioner. It must also contain all ghost cells via facet, i.e. cells that are on a neighboring process and which share a facet with a local cell. Ghost cells are the last n entries in cells[i], where n is given by the length of ghost_owners[i].

  • original_cell_index[in] Input cell index for each cell type, e.g. the cell index in an input file. This index remains associated with the cell after any re-ordering and parallel (re)distribution.

  • ghost_owners[in] Owning rank for ghost cells (ghost cells are at end of each list of cells).

  • boundary_vertices[in] Vertices on the ‘exterior’ (boundary) of the local topology. These vertices might appear on other processes.

Returns:

A distributed mesh topology.

Topology create_topology(MPI_Comm comm, std::span<const std::int64_t> cells, std::span<const std::int64_t> original_cell_index, std::span<const int> ghost_owners, CellType cell_type, std::span<const std::int64_t> boundary_vertices)

Create a mesh topology for a single cell type.

This function provides a simplified interface to create_topology for the case that a mesh has one cell type only,

Parameters:
  • comm[in] Communicator across which the topology will be distributed.

  • cells[in] Cell topology (list of vertices for each cell) using global indices for the vertices. It contains cells that have been distributed to this rank, e.g. via a graph partitioner. It must also contain all ghost cells via facet, i.e. cells that are on a neighboring process and which share a facet with a local cell. Ghost cells are the last n entries in cells, where n is given by the length of ghost_owners.

  • original_cell_index[in] Original global index associated with each cell.

  • ghost_owners[in] Owning rank of each ghost cell (ghost cells are always at the end of the list of cells).

  • cell_type[in] A vector with cell shapes.

  • boundary_vertices[in] Vertices on the ‘exterior’ (boundary) of the local topology. These vertices might appear on other processes.

Returns:

A distributed mesh topology.

std::tuple<Topology, std::vector<int32_t>, std::vector<int32_t>> create_subtopology(const Topology &topology, int dim, std::span<const std::int32_t> entities)

Create a topology for a subset of entities of a given topological dimension.

Parameters:
  • topology[in] Original (parent) topology.

  • dim[in] Topological dimension of the entities in the new topology.

  • entities[in] Indices of entities in topology to include in the new topology.

Returns:

New topology of dimension dim with all entities in entities, map from entities of dimension dim in new sub-topology to entities in topology, and map from vertices in new sub-topology to vertices in topology.

std::vector<std::int32_t> entities_to_index(const Topology &topology, int dim, std::span<const std::int32_t> entities)

Get entity indices for entities defined by their vertices.

Warning

This function may be removed in the future.

Parameters:
  • topology[in] Mesh topology.

  • dim[in] Topological dimension of the entities.

  • entities[in] Mesh entities defined by their vertices.

Returns:

Index of the ith entity in entities, If an entity cannot be found on this rank, -1 is returned as the index.

std::tuple<std::vector<std::shared_ptr<graph::AdjacencyList<std::int32_t>>>, std::shared_ptr<graph::AdjacencyList<std::int32_t>>, std::shared_ptr<common::IndexMap>, std::vector<std::int32_t>> compute_entities(const Topology &topology, int dim, CellType entity_type)

Compute mesh entities of given topological dimension by computing cell-to-entity (tdim, i) ->(dim, entity_type)and / entity-to-vertex connectivity(dim, entity_type) -> (0, 0) connectivity.

Computed entities are oriented such that their local (to the process) orientation agrees with their global orientation

Parameters:
  • topology[in] Mesh topology.

  • dim[in] Dimension of the entities to create.

  • entity_type[in] Entity type in dimension dim to create. Entity type must be in the list returned by Topology::entity_types.

Returns:

Tuple of (cell->entity connectivity, entity->vertex connectivity, index map for created entities, list of interprocess entities). Interprocess entities lie on the “true” boundary between owned cells of each process. If entities of type entity_type already exists, then {nullptr, nullptr, nullptr, std::vector()} is returned.

std::array<std::shared_ptr<graph::AdjacencyList<std::int32_t>>, 2> compute_connectivity(const Topology &topology, std::array<int, 2> d0, std::array<int, 2> d1)

Compute connectivity (d0 -> d1) for given pair of entity types, given by topological dimension and index, as found in Topology::entity_types()

Parameters:
  • topology[in] The topology

  • d0[in] Dimension and index of the entities, (dim0, i).

  • d1[in] Dimension and index of the incident entities, (dim1, / j).

Returns:

The connectivities [(d0 -> d1), (d1 -> d0)] if they are computed. If (d0, d1) already exists then a nullptr is returned. If (d0, d1) is computed and the computation of (d1, d0) was required as part of computing (d0, d1), the (d1, d0) is returned as the second entry. The second entry is otherwise nullptr.

std::vector<std::int32_t> exterior_facet_indices(const Topology &topology)

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

Collective.

Parameters:

topology[in] Mesh topology.

Returns:

Sorted list of owned facet indices that are exterior facets of the mesh.

std::vector<std::int64_t> extract_topology(CellType cell_type, const fem::ElementDofLayout &layout, std::span<const std::int64_t> cells)

Extract topology from cell data, i.e. extract cell vertices.

Parameters:
  • cell_type[in] Cell shape.

  • layout[in] Layout of geometry ‘degrees-of-freedom’ on the reference cell.

  • cells[in] List of ‘nodes’ for each cell using global indices. The layout must be consistent with layout.

Returns:

Cell topology. The global indices will, in general, have ‘gaps’ due to mid-side and other higher-order nodes being removed from the input cell.

template<std::floating_point T>
std::vector<T> h(const Mesh<T> &mesh, std::span<const std::int32_t> entities, int dim)

Compute greatest distance between any two vertices of the mesh entities (h).

Parameters:
  • mesh[in] Mesh that the entities belong to.

  • entities[in] Indices (local to process) of entities to compute h for.

  • dim[in] Topological dimension of the entities.

Returns:

Greatest distance between any two vertices, h[i] corresponds to the entity entities[i].

template<std::floating_point T>
std::vector<T> cell_normals(const Mesh<T> &mesh, int dim, std::span<const std::int32_t> entities)

Compute normal to given cell (viewed as embedded in 3D).

Returns:

The entity normals. The shape is (entities.size(), 3) and the storage is row-major.

template<std::floating_point T>
std::vector<T> compute_midpoints(const Mesh<T> &mesh, int dim, std::span<const std::int32_t> entities)

Compute the midpoints for mesh entities of a given dimension.

Returns:

The entity midpoints. The shape is (entities.size(), 3) and the storage is row-major.

template<std::floating_point T, MarkerFn<T> U>
std::vector<std::int32_t> locate_entities(const Mesh<T> &mesh, int dim, U marker)

Compute indices of all mesh entities that evaluate to true for the provided geometric marking function.

An entity is considered marked if the marker function evaluates to true for all of its vertices.

Parameters:
  • mesh[in] Mesh to mark entities on.

  • dim[in] Topological dimension of the entities to be considered.

  • marker[in] Marking function, returns true for a point that is ‘marked’, and false otherwise.

Returns:

List of marked entity indices, including any ghost indices (indices local to the process).

template<std::floating_point T, MarkerFn<T> U>
std::vector<std::int32_t> locate_entities_boundary(const Mesh<T> &mesh, int dim, U marker)

Compute indices of all mesh entities that are attached to an owned boundary facet and evaluate to true for the provided geometric marking function.

An entity is considered marked if the marker function evaluates to true for all of its vertices.

Note

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[in] Mesh to mark entities on.

  • dim[in] Topological dimension of the entities to be considered. Must be less than the topological dimension of the mesh.

  • marker[in] Marking function, returns true for a point that is ‘marked’, and false otherwise.

Returns:

List of marked entity indices (indices local to the process).

template<std::floating_point T>
std::vector<std::int32_t> entities_to_geometry(const Mesh<T> &mesh, int dim, std::span<const std::int32_t> entities, bool permute = false)

Compute the geometry degrees of freedom associated with the closure of a given set of cell entities.

Parameters:
  • mesh[in] The mesh.

  • dim[in] Topological dimension of the entities of interest.

  • entities[in] Entity indices (local to process).

  • permute[in] If true, 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 geometry DOFs associated with the closure of each entity in entities. The shape is (num_entities, num_xdofs_per_entity) and the storage is row-major. The index indices[i, j] is the position in the geometry array of the j-th vertex of the entity[i].

Pre:

The mesh connectivities dim -> mesh.topology().dim() and mesh.topology().dim() -> dim must have been computed. Otherwise an exception is thrown.

CellPartitionFunction create_cell_partitioner(mesh::GhostMode ghost_mode = mesh::GhostMode::none, const graph::partition_fn &partfn = &graph::partition_graph)

Create a function that computes destination rank for mesh cells on this rank by applying the default graph partitioner to the dual graph of the mesh.

Returns:

Function that computes the destination ranks for each cell.

std::vector<std::int32_t> compute_incident_entities(const Topology &topology, std::span<const std::int32_t> entities, int d0, int d1)

Compute incident entities.

Parameters:
  • topology[in] The topology.

  • entities[in] List of indices of topological dimension d0.

  • d0[in] Topological dimension.

  • d1[in] Topological dimension.

Returns:

List of entities of topological dimension d1 that are incident to entities in entities (topological dimension d0).

template<typename U>
Mesh<typename std::remove_reference_t<typename U::value_type>> create_mesh(MPI_Comm comm, MPI_Comm commt, std::vector<std::span<const std::int64_t>> cells, const std::vector<fem::CoordinateElement<typename std::remove_reference_t<typename U::value_type>>> &elements, MPI_Comm commg, const U &x, std::array<std::size_t, 2> xshape, const CellPartitionFunction &partitioner)

Create a distributed mesh::Mesh from mesh data and using the provided graph partitioning function for determining the parallel distribution of the mesh.

The input cells and geometry data can be distributed across the calling ranks, but must be not duplicated across ranks.

The function partitioner computes the parallel distribution, i.e. the destination rank for each cell passed to the constructor. If partitioner is not callable, i.e. it does not store a callable function, no parallel re-distribution of cells is performed.

Note

Collective.

Parameters:
  • comm[in] Communicator to build the mesh on.

  • commt[in] Communicator that the topology data (cells) is distributed on. This should be MPI_COMM_NULL for ranks that should not participate in computing the topology partitioning.

  • cells[in] Cells, grouped by cell type with cells[i] being the cells of the same type. Cells are defined by their ‘nodes’ (using global indices) following the Basix ordering, and for each cell type concatenated to form a flattened list. For lowest-order cells this will be just the cell vertices. For higher-order geometry cells, other cell ‘nodes’ will be included. See io::cells for examples of the Basix ordering.

  • elements[in] Coordinate elements for the cells, where elements[i] is the coordinate element for the cells in cells[i]. The list of elements must be the same on all calling parallel ranks.

  • commg[in] Communicator for geometry.

  • x[in] Geometry data (‘node’ coordinates). Row-major storage. The global index of the ith node (row) in x is taken as i plus the parallel rank offset (on comm), where the offset is the sum of x rows on all lower ranks than the caller.

  • xshape[in] Shape of the x data.

  • partitioner[in] Graph partitioner that computes the owning rank for each cell in cells. If not callable, cells are not redistributed.

Returns:

A mesh distributed on the communicator comm.

template<typename U>
Mesh<typename std::remove_reference_t<typename U::value_type>> create_mesh(MPI_Comm comm, MPI_Comm commt, std::span<const std::int64_t> cells, const fem::CoordinateElement<typename std::remove_reference_t<typename U::value_type>> &element, MPI_Comm commg, const U &x, std::array<std::size_t, 2> xshape, const CellPartitionFunction &partitioner)

Create a distributed mesh with a single cell type from mesh data and using a provided graph partitioning function for determining the parallel distribution of the mesh.

From mesh input data that is distributed across processes, a distributed mesh::Mesh is created. If the partitioning function is not callable, i.e. it does not store a callable function, no re-distribution of cells is done.

This constructor provides a simplified interface to the more general create_mesh constructor, which supports meshes with more than one cell type.

Parameters:
  • comm[in] Communicator to build the mesh on.

  • commt[in] Communicator that the topology data (cells) is distributed on. This should be MPI_COMM_NULL for ranks that should not participate in computing the topology partitioning.

  • cells[in] Cells on the calling process. Each cell (node in the AdjacencyList) is defined by its ‘nodes’ (using global indices) following the Basix ordering. For lowest order cells this will be just the cell vertices. For higher-order cells, other cells ‘nodes’ will be included. See dolfinx::io::cells for examples of the Basix ordering.

  • element[in] Coordinate element for the cells.

  • commg[in] Communicator for geometry.

  • x[in] Geometry data (‘node’ coordinates). Row-major storage. The global index of the ith node (row) in x is taken as i plus the process offset oncomm, The offset is the sum of x rows on all processed with a lower rank than the caller.

  • xshape[in] Shape of the x data.

  • partitioner[in] Graph partitioner that computes the owning rank for each cell. If not callable, cells are not redistributed.

Returns:

A mesh distributed on the communicator comm.

template<typename U>
Mesh<typename std::remove_reference_t<typename U::value_type>> create_mesh(MPI_Comm comm, std::span<const std::int64_t> cells, const fem::CoordinateElement<std::remove_reference_t<typename U::value_type>> &elements, const U &x, std::array<std::size_t, 2> xshape, GhostMode ghost_mode)

Create a distributed mesh from mesh data using the default graph partitioner to determine the parallel distribution of the mesh.

This function takes mesh input data that is distributed across processes and creates a mesh::Mesh, with the mesh cell distribution determined by the default cell partitioner. The default partitioner is based on graph partitioning.

Parameters:
  • comm[in] MPI communicator to build the mesh on.

  • cells[in] Cells on the calling process. See create_mesh for a detailed description.

  • elements[in] Coordinate elements for the cells.

  • x[in] Geometry data (‘node’ coordinates). See create_mesh for a detailed description.

  • xshape[in] Shape of x. It should be (num_points, gdim).

  • ghost_mode[in] Required type of cell ghosting/overlap.

Returns:

A mesh distributed on the communicator comm.

template<std::floating_point T>
std::pair<Geometry<T>, std::vector<int32_t>> create_subgeometry(const Mesh<T> &mesh, int dim, std::span<const std::int32_t> subentity_to_entity)

Create a sub-geometry from a mesh and a subset of mesh entities to be included.

A sub-geometry is simply a mesh::Geometry object containing only the geometric information for the subset of entities. The entities may differ in topological dimension from the original mesh.

Parameters:
  • mesh[in] The full mesh.

  • dim[in] Topological dimension of the sub-topology.

  • subentity_to_entity[in] Map from sub-topology entity to the entity in the parent topology.

Returns:

A sub-geometry and a map from sub-geometry coordinate degree-of-freedom to the coordinate degree-of-freedom in geometry.

template<std::floating_point T>
std::tuple<Mesh<T>, std::vector<std::int32_t>, std::vector<std::int32_t>, std::vector<std::int32_t>> create_submesh(const Mesh<T> &mesh, int dim, std::span<const std::int32_t> entities)

Create a new mesh consisting of a subset of entities in a mesh.

Parameters:
  • mesh[in] The mesh.

  • dim[in] Dimension entities in mesh that will be cells in the sub-mesh.

  • entities[in] Indices of entities in mesh to include in the sub-mesh.

Returns:

A new mesh, and maps from the new mesh entities, vertices, and geometry to the input mesh entities, vertices, and geometry.

template<std::floating_point T>
class Geometry
#include <Geometry.h>

Geometry stores the geometry imposed on a mesh.

Public Types

using value_type = T

Value type.

Public Functions

template<typename U, typename V, typename W>
inline Geometry(std::shared_ptr<const common::IndexMap> index_map, U &&dofmaps, const std::vector<fem::CoordinateElement<typename std::remove_reference_t<typename V::value_type>>> &elements, V &&x, int dim, W &&input_global_indices)

Constructor of object that holds mesh geometry data.

The geometry ‘degree-of-freedom map’ is a logically rectangular array (row-major) where each row corresponds to a cell, and the columns are the indices in the coordinate array of the cell coordinate degree-of-freedom. For each cell type there is a geometry degree-of-freedom map.

Parameters:
  • index_map[in] Index map associated with the geometry degrees-of-freedom.

  • dofmaps[in] The geometry (point) dofmaps for each cell type in the mesh. For a cell of a given type, the dofmap gives the position in the point array of each local geometry node of the cell. Each cell type has its own dofmap. Each dofmap uses row-major storage.

  • elements[in] Elements that describe cell geometry maps. element[i] is the coordinate element associated with dofmaps[i].

  • x[in] Point coordinates. The shape is (num_points, 3) and the storage is row-major.

  • dim[in] The geometric dimension (0 < dim <= 3).

  • input_global_indices[in] ‘Global’ input index of each point, commonly from a mesh input file.

Geometry(const Geometry&) = default

Copy constructor.

Geometry(Geometry&&) = default

Move constructor.

~Geometry() = default

Destructor.

Geometry &operator=(const Geometry&) = delete

Copy assignment.

Geometry &operator=(Geometry&&) = default

Move assignment.

inline int dim() const

Return dimension of the Euclidean coordinate system.

inline MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<const std::int32_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>> dofmap() const

DofMap for the geometry.

Returns:

A 2D array with shape (num_cells, dofs_per_cell).

inline MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<const std::int32_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>> dofmap(std::size_t i) const

Degree-of-freedom map associated with the ith coordinate map element in the geometry.

Parameters:

i[in] Index of the requested degree-of-freedom map. The degree-of-freedom map corresponds to the geometry element cmaps()[i].

Returns:

A dofmap array, with shape (num_cells, dofs_per_cell).

inline std::shared_ptr<const common::IndexMap> index_map() const

Index map for the geometry ‘degrees-of-freedom’.

Returns:

Index map for the geometry dofs.

inline std::span<const value_type> x() const

Access geometry degrees-of-freedom data (const version).

Returns:

The flattened row-major geometry data, where the shape is (num_points, 3).

inline std::span<value_type> x()

Access geometry degrees-of-freedom data (non-const version).

Returns:

The flattened row-major geometry data, where the shape is (num_points, 3).

inline const std::vector<fem::CoordinateElement<value_type>> &cmaps() const

The elements that describes the geometry map.

The coordinate element cmaps()[i] corresponds to the degree-of-freedom map dofmap(i).

Returns:

The coordinate/geometry elements.

inline const fem::CoordinateElement<value_type> &cmap() const

The element that describes the geometry map.

Returns:

The coordinate/geometry element

inline const std::vector<std::int64_t> &input_global_indices() const

Global user indices.

template<std::floating_point T>
class Mesh
#include <Mesh.h>

A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.

Template Parameters:

Floating – point type for representing the geometry.

Public Types

using geometry_type = Geometry<T>

Value type.

Public Functions

template<typename V>
inline Mesh(MPI_Comm comm, std::shared_ptr<Topology> topology, V &&geometry)

Create a Mesh.

Note

This constructor is not normally called by users. User code will normally use create_mesh.

Parameters:
  • comm[in] MPI Communicator.

  • topology[in] Mesh topology.

  • geometry[in] Mesh geometry.

Mesh(const Mesh &mesh) = default

Copy constructor

Parameters:

mesh[in] Mesh to be copied

Mesh(Mesh &&mesh) = default

Move constructor

Parameters:

meshMesh to be moved.

~Mesh() = default

Destructor.

Mesh &operator=(Mesh &&mesh) = default

Assignment move operator

Parameters:

mesh – Another Mesh object

inline std::shared_ptr<Topology> topology()

Get mesh topology.

Returns:

The topology object associated with the mesh.

inline std::shared_ptr<const Topology> topology() const

Get mesh topology (const version).

Returns:

The topology object associated with the mesh.

inline std::shared_ptr<Topology> topology_mutable() const

Get mesh topology if one really needs the mutable version.

Returns:

The topology object associated with the mesh.

inline Geometry<T> &geometry()

Get mesh geometry.

Returns:

The geometry object associated with the mesh.

inline const Geometry<T> &geometry() const

Get mesh geometry (const version).

Returns:

The geometry object associated with the mesh.

inline MPI_Comm comm() const

Mesh MPI communicator.

Returns:

The communicator on which the mesh is distributed.

Public Members

std::string name = "mesh"

Name.

template<typename T>
class MeshTags
#include <MeshTags.h>

MeshTags associate values with mesh topology entities.

The entity index (local to process) identifies the entity. MeshTags is a sparse data storage class; it allows tags to be associated with an arbitrary subset of mesh entities. An entity can have only one associated tag.

Template Parameters:

Type

Public Functions

template<typename U, typename V>
inline MeshTags(std::shared_ptr<const Topology> topology, int dim, U &&indices, V &&values)

Create a MeshTag from entities of given dimension on a mesh.

Parameters:
  • topology[in] Mesh topology on which the tags are associated.

  • dim[in] Topological dimension of mesh entities to tag.

  • indices[in] List of entity indices (indices local to the process).

  • values[in] List of values for each index in indices. The size must be equal to the size of indices.

Pre:

indices must be sorted and unique.

MeshTags(const MeshTags &tags) = default

Copy constructor.

MeshTags(MeshTags &&tags) = default

Move constructor.

~MeshTags() = default

Destructor.

MeshTags &operator=(const MeshTags &tags) = default

Move assignment.

MeshTags &operator=(MeshTags &&tags) = default

Move assignment.

inline std::vector<std::int32_t> find(const T value) const

Find all entities with a given tag value.

Parameters:

value[in] The value

Returns:

Indices of tagged entities. The indices are sorted.

inline std::span<const std::int32_t> indices() const

Indices of tagged topology entities (local-to-process). The indices are sorted.

inline std::span<const T> values() const

Values attached to topology entities.

inline int dim() const

Return topological dimension of tagged entities.

inline std::shared_ptr<const Topology> topology() const

Return topology.

Public Members

std::string name = "mesh_tags"

Name.

class Topology
#include <Topology.h>

Topology stores the topology of a mesh, consisting of mesh entities and connectivity (incidence relations for the mesh entities).

A mesh entity e may be identified globally as a pair e = (dim, i), where dim is the topological dimension and i is the index of the entity within that topological dimension.

Todo:

Rework memory management and associated API. Currently, the caching policy is not clear.

Public Functions

Topology(std::vector<CellType> cell_types, std::shared_ptr<const common::IndexMap> vertex_map, std::vector<std::shared_ptr<const common::IndexMap>> cell_maps, std::vector<std::shared_ptr<graph::AdjacencyList<std::int32_t>>> cells, const std::optional<std::vector<std::vector<std::int64_t>>> &original_cell_index = std::nullopt)

Create a mesh topology.

A Topology represents the connectivity of a mesh. Mesh entities, i.e. vertices, edges, faces and cells, are defined in terms of their vertices. Connectivity represents the relationships between entities, e.g. the cells that are connected to a given edge in the mesh.

Parameters:
  • cell_types[in] Types of cells.

  • vertex_map[in] Index map describing the distribution of mesh vertices.

  • cell_maps[in] Index maps describing the distribution of mesh cells for each cell type in cell_types.

  • cells[in] Cell-to-vertex connectivities for each cell type in cell_types.

  • original_cell_index[in] Original indices for each cell in cells.

Topology(const Topology &topology) = default

Copy constructor.

Topology(Topology &&topology) = default

Move constructor.

~Topology() = default

Destructor.

Topology &operator=(const Topology &topology) = delete

Assignment.

Topology &operator=(Topology &&topology) = default

Assignment.

int dim() const noexcept

Topological dimension of the mesh.

const std::vector<CellType> &entity_types(int dim) const

Entity types in the topology for a given dimension.

Parameters:

dim[in] Topological dimension.

Returns:

Entity types.

CellType cell_type() const

Cell type.

This function is is for topologies with one cell type only.

Returns:

Cell type that the topology is for.

std::vector<std::shared_ptr<const common::IndexMap>> index_maps(int dim) const

Get the index maps that described the parallel distribution of the mesh entities of a given topological dimension.

Parameters:

dim[in] Topological dimension.

Returns:

Index maps, one for each cell type.

std::shared_ptr<const common::IndexMap> index_map(int dim) const

Get the IndexMap that described the parallel distribution of the mesh entities.

Parameters:

dim[in] Topological dimension

Returns:

Index map for the entities of dimension dim. Returns nullptr if index map has not been set.

std::shared_ptr<const graph::AdjacencyList<std::int32_t>> connectivity(std::array<int, 2> d0, std::array<int, 2> d1) const

Get the connectivity from entities of topological dimension d0 to dimension d1.

The entity type and incident entity type are each described by a pair (dim, index). The index within a topological dimension dim, is that of the cell type given in entity_types(dim).

Parameters:
  • d0[in] Pair of (topological dimension of entities, index of “entity type” within topological dimension).

  • d1[in] Pair of (topological dimension of entities, index of incident “entity type” within topological dimension).

Returns:

AdjacencyList of connectivity from entity type in d0 to entity types in d1, or nullptr if not yet computed.

std::shared_ptr<const graph::AdjacencyList<std::int32_t>> connectivity(int d0, int d1) const

Return connectivity from entities of dimension d0 to entities of dimension d1. Assumes only one entity type per dimension.

Parameters:
  • d0[in] Topological dimension.

  • d1[in] Topological dimension.

Returns:

The adjacency list that for each entity of dimension d0 gives the list of incident entities of dimension d1. Returns nullptr if connectivity has not been computed.

const std::vector<std::uint32_t> &get_cell_permutation_info() const

Returns the permutation information.

const std::vector<std::uint8_t> &get_facet_permutations() const

Get the numbers that encode the number of permutations to apply to facets.

The permutations are encoded so that:

  • n % 2 gives the number of reflections to apply

  • n // 2 gives the number of rotations to apply

The data is stored in a flattened 2D array, so that data[cell_index * / facets_per_cell + facet_index] contains the facet with index facet_index of the cell with index cell_index.

Note

An exception is raised if the permutations have not been computed

Returns:

The encoded permutation info

const std::vector<std::int32_t> &interprocess_facets(int index) const

List of inter-process facets of a given type.

“Inter-process” facets are facets that are connected (1) to a cell that is owned by the calling process (rank) and (2) to a cell that is owned by another process.

Facets must have been computed for inter-process facet data to be available.

Parameters:

index[in] Index of facet type, following the order given by entity_types.

Returns:

Indices of the inter-process facets.

const std::vector<std::int32_t> &interprocess_facets() const

List of inter-process facets.

“Inter-process” facets are facets that are connected (1) to a cell that is owned by the calling process (rank) and (2) to a cell that is owned by another process.

Pre:

Inter-process facets are available only if facet topology has been computed.

bool create_entities(int dim)

Create entities of given topological dimension.

Parameters:

dim[in] Topological dimension of entities to compute.

Returns:

True if entities are created, false if entities already existed.

void create_connectivity(int d0, int d1)

Create connectivity between given pair of dimensions, d0 / -> d1.

Parameters:
  • d0[in] Topological dimension.

  • d1[in] Topological dimension.

void create_entity_permutations()

Compute entity permutations and reflections.

MPI_Comm comm() const

Mesh MPI communicator.

Returns:

Communicator on which the topology is distributed.

Public Members

std::vector<std::vector<std::int64_t>> original_cell_index

Original cell index for each cell type.

namespace impl

Functions

template<std::floating_point T>
std::tuple<std::vector<T>, std::vector<std::int64_t>> create_interval_cells(std::array<T, 2> p, std::int64_t n)
template<std::floating_point T>
Mesh<T> build_tri(MPI_Comm comm, std::array<std::array<T, 2>, 2> p, std::array<std::int64_t, 2> n, const CellPartitionFunction &partitioner, DiagonalType diagonal)
template<std::floating_point T>
Mesh<T> build_quad(MPI_Comm comm, const std::array<std::array<T, 2>, 2> p, std::array<std::int64_t, 2> n, const CellPartitionFunction &partitioner)
template<std::floating_point T>
std::vector<T> create_geom(MPI_Comm comm, std::array<std::array<T, 3>, 2> p, std::array<std::int64_t, 3> n)
template<std::floating_point T>
Mesh<T> build_tet(MPI_Comm comm, MPI_Comm subcomm, std::array<std::array<T, 3>, 2> p, std::array<std::int64_t, 3> n, const CellPartitionFunction &partitioner)
template<std::floating_point T>
Mesh<T> build_hex(MPI_Comm comm, MPI_Comm subcomm, std::array<std::array<T, 3>, 2> p, std::array<std::int64_t, 3> n, const CellPartitionFunction &partitioner)
template<std::floating_point T>
Mesh<T> build_prism(MPI_Comm comm, MPI_Comm subcomm, std::array<std::array<T, 3>, 2> p, std::array<std::int64_t, 3> n, const CellPartitionFunction &partitioner)
template<std::floating_point T>
mesh::Mesh<T> build_hex(MPI_Comm comm, MPI_Comm subcomm, std::array<std::array<T, 3>, 2> p, std::array<std::int64_t, 3> n, const CellPartitionFunction &partitioner)
template<typename T>
void reorder_list(std::span<T> list, std::span<const std::int32_t> nodemap)

Re-order the nodes of a fixed-degree adjacency list.

Parameters:
  • list[inout] Fixed-degree adjacency list stored row-major. Degree is equal to list.size() / nodemap.size().

  • nodemap[in] Map from old to new index, i.e. for an old index i the new index is nodemap[i].

template<std::floating_point T>
std::tuple<std::vector<std::int32_t>, std::vector<T>, std::vector<std::int32_t>> compute_vertex_coords_boundary(const mesh::Mesh<T> &mesh, int dim, std::span<const std::int32_t> facets)

Compute the coordinates of ‘vertices’ for entities of a given dimension that are attached to specified facets.

Parameters:
  • mesh[in] Mesh to compute the vertex coordinates for.

  • dim[in] Topological dimension of the entities.

  • facets[in] List of facets (must be on the mesh boundary).

Pre:

The provided facets must be on the boundary of the mesh.

Returns:

(0) Entities attached to the boundary facets (sorted), (1) vertex coordinates (shape is (3, num_vertices)) and (2) map from vertex in the full mesh to the position in the vertex coordinates array (set to -1 if vertex in full mesh is not in the coordinate array).

template<std::floating_point T>
std::pair<std::vector<T>, std::array<std::size_t, 2>> compute_vertex_coords(const mesh::Mesh<T> &mesh)

The coordinates for all ‘vertices’ in the mesh.

Parameters:

mesh[in] Mesh to compute the vertex coordinates for.

Returns:

The vertex coordinates. The shape is (3, num_vertices) and the jth column hold the coordinates of vertex j.