Mesh (dolfinx::mesh)

Under development

namespace dolfinx::mesh

Mesh data structures and algorithms on meshes.

Representations of meshes and support for operations on meshes.

Typedefs

using CellPartitionFunction = std::function<dolfinx::graph::AdjacencyList<std::int32_t>(MPI_Comm comm, int nparts, int tdim, const dolfinx::graph::AdjacencyList<std::int64_t> &cells, dolfinx::mesh::GhostMode ghost_mode)>

Signature for the cell partitioning function. The function should compute the destination rank for cells currently on this rank.

Param comm

[in] MPI Communicator

Param nparts

[in] Number of partitions

Param tdim

[in] Topological dimension

Param cells

[in] Cells on this process. The ith entry in list contains 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.

Param ghost_mode

[in] How to overlap the cell partitioning: none, shared_facet or shared_vertex

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.

Mesh create_box(MPI_Comm comm, const std::array<std::array<double, 3>, 2> &p, std::array<std::size_t, 3> n, CellType celltype, GhostMode ghost_mode, const mesh::CellPartitionFunction &partitioner = create_cell_partitioner())

Create a uniform mesh::Mesh over the 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 build mesh on

  • p[in] Points of box

  • n[in] Number of cells in each direction

  • celltype[in] Cell shape

  • ghost_mode[in] Ghost mode

  • partitioner[in] Partitioning function to use for determining the parallel distribution of cells across MPI ranks

Returns

Mesh

Mesh create_interval(MPI_Comm comm, std::size_t n, std::array<double, 2> x, GhostMode ghost_mode, const CellPartitionFunction &partitioner = create_cell_partitioner())

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] The number of cells

  • x[in] The end points of the interval

  • ghost_mode[in] Ghosting mode

  • partitioner[in] Partitioning function to use for determining the parallel distribution of cells across MPI ranks

Returns

A mesh

Mesh create_rectangle(MPI_Comm comm, const std::array<std::array<double, 2>, 2> &p, std::array<std::size_t, 2> n, CellType celltype, GhostMode ghost_mode, 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

  • ghost_mode[in] Mesh ghosting mode

  • diagonal[in] Direction of diagonals

Returns

Mesh

Mesh create_rectangle(MPI_Comm comm, const std::array<std::array<double, 2>, 2> &p, std::array<std::size_t, 2> n, CellType celltype, GhostMode ghost_mode, const 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] Two corner points

  • n[in] Number of cells in each direction

  • celltype[in] Cell shape

  • ghost_mode[in] Mesh ghosting mode

  • partitioner[in] Partitioning function to use for determining the parallel distribution of cells across MPI ranks

  • diagonal[in] Direction of diagonals

Returns

Mesh

mesh::Geometry create_geometry(MPI_Comm comm, const Topology &topology, const fem::CoordinateElement &element, const graph::AdjacencyList<std::int64_t> &cells, const xtl::span<const double> &x, int dim, const 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. It distributes the ‘node’ coordinate data to the required MPI process and then creates a mesh::Geometry object.

Parameters
  • comm[in] The MPI communicator to build the Geometry on

  • topology[in] The mesh topology

  • element[in] The element that defines the geometry map for each cell

  • cells[in] The mesh cells, including higher-ordder geometry ‘nodes’

  • 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] The geometric dimension (1, 2, or 3)

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

std::tuple<graph::AdjacencyList<std::int32_t>, std::vector<std::int64_t>, std::size_t, std::vector<std::int32_t>> build_local_dual_graph(const xtl::span<const std::int64_t> &cells, const xtl::span<const std::int32_t> &offsets, int tdim)

Compute the local part of the dual graph (cell-cell connections via facets) and facet 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 padding value (all padding values will be equal).

Note

The return data will likely change once we support mixed topology meshes.

Parameters
  • cells[in] Array for cell vertices adjacency list

  • offsets[in] Adjacency list offsets, i.e. index of the first entry of cell i in cell_vertices is offsets[i]

  • tdim[in] The topological dimension of the cells

Returns

  1. Local dual graph

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

  3. The number of columns for the facet data array (2).

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

graph::AdjacencyList<std::int64_t> build_dual_graph(const MPI_Comm comm, const graph::AdjacencyList<std::int64_t> &cells, int tdim)

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

Parameters
  • comm[in] The MPI communicator

  • cells[in] Collection of cells, defined by the cell vertices from which to build the dual graph

  • tdim[in] The topological dimension of the cells

Returns

The dual graph

Mesh create_mesh(MPI_Comm comm, const graph::AdjacencyList<std::int64_t> &cells, const fem::CoordinateElement &element, const xt::xtensor<double, 2> &x, GhostMode ghost_mode)

Create a mesh using the default partitioner.

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 a graph partitioning.

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

  • cells[in] The cells on the this MPI rank. Each cell (node in the AdjacencyList) is defined by its ‘nodes’ (using global indices). For lowest order cells this will be just the cell vertices. For higher-order cells, other cells ‘nodes’ will be included.

  • element[in] The coordinate element that describes the geometric mapping for cells

  • x[in] The coordinates of mesh nodes

  • ghost_mode[in] The requested type of cell ghosting/overlap

Returns

A distributed Mesh.

Mesh create_mesh(MPI_Comm comm, const graph::AdjacencyList<std::int64_t> &cells, const fem::CoordinateElement &element, const xt::xtensor<double, 2> &x, GhostMode ghost_mode, const CellPartitionFunction &cell_partitioner)

Create a mesh using a provided mesh partitioning function.

std::tuple<Mesh, std::vector<std::int32_t>, std::vector<std::int32_t>, std::vector<std::int32_t>> create_submesh(const Mesh &mesh, int dim, const xtl::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] Entity dimension

  • entities[in] List of entity indicies in mesh to include in the new mesh

Returns

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

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

Create MeshTags from arrays.

Note

Entities that do not exist on this rank are ignored.

Warning

entities must not conatined duplicate entities.

Parameters
  • mesh[in] The Mesh 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 assemble of (1) facet inetgrals and (2) vector elements.

  1. Get the permutation numbers to apply to facets. The permutations are numbered so that:

    • n % 2 gives the number of reflections to apply

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

    Each column of the returned array represents a cell, and each row a facet of that cell. This data is used to permute the quadrature point on facet integrals when data from the cells on both sides of the facet is used.

  2. Get the permutation 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

std::vector<std::int8_t> compute_boundary_facets(const Topology &topology)

Compute marker for owned facets that are on the exterior of the domain, i.e. are connected to only one cell.

This function does not require parallel communication.

Parameters

topology[in] The topology

Returns

Vector with length equal to the number of owned facets on this this process. True if the ith facet (local index) is on the exterior of the domain.

Topology create_topology(MPI_Comm comm, const graph::AdjacencyList<std::int64_t> &cells, const xtl::span<const std::int64_t> &original_cell_index, const xtl::span<const int> &ghost_owners, const CellType &cell_type, GhostMode ghost_mode)

Create a distributed mesh topology.

Parameters
  • comm[in] MPI communicator across which the topology is distributed

  • cells[in] The 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 share a facet with a local cell.

  • original_cell_index[in] The original global index associated with each cell

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

  • cell_type[in] The cell shape

  • ghost_mode[in] Type of cell ghosting: none, shared_facet or shared_vertex

Returns

A distributed mesh topology

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

Get entity indices for entities defined by their vertices.

Note

If an entity cannot be found on this rank, -1 is returned as the index.

Warning

This function may be removed in the future.

Parameters
  • topology[in] The mesh topology

  • dim[in] Topological dimension of the entities

  • entities[in] The mesh entities defined by their vertices

Returns

The index of the ith entity in entities

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

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

Parameters
  • comm[in] MPI Communicator

  • topology[in] Mesh topology

  • dim[in] The dimension of the entities to create

Returns

Tuple of (cell-entity connectivity, entity-vertex connectivity, index map). If the entities already exist, then {nullptr, nullptr, nullptr} is returned.

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

Compute connectivity (d0 -> d1) for given pair of topological dimensions.

Parameters
  • topology[in] The topology

  • d0[in] The dimension of the nodes in the adjacency list

  • d1[in] The dimension of the edges in the adjacency list

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.

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

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

Parameters
  • cell_type[in] The cell shape

  • layout[in] The 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.

std::vector<double> h(const Mesh &mesh, const xtl::span<const std::int32_t> &entities, int dim)

Compute greatest distance between any two vertices.

xt::xtensor<double, 2> cell_normals(const Mesh &mesh, int dim, const xtl::span<const std::int32_t> &entities)

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

xt::xtensor<double, 2> compute_midpoints(const Mesh &mesh, int dim, const xtl::span<const std::int32_t> &entities)

Compute the midpoints for mesh entities of a given dimension.

std::vector<std::int32_t> locate_entities(const Mesh &mesh, int dim, const std::function<xt::xtensor<bool, 1>(const xt::xtensor<double, 2>&)> &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 true for all of its vertices.

Parameters
  • mesh[in] The mesh

  • dim[in] The topological dimension of the entities to be considered

  • marker[in] The marking function

Returns

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

std::vector<std::int32_t> locate_entities_boundary(const Mesh &mesh, int dim, const std::function<xt::xtensor<bool, 1>(const xt::xtensor<double, 2>&)> &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 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] The mesh

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

  • marker[in] The marking function

Returns

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

xt::xtensor<std::int32_t, 2> entities_to_geometry(const Mesh &mesh, int dim, const xtl::span<const std::int32_t> &entity_list, bool orient)

Compute the indices the geometry data for the vertices of the given mesh entities.

Parameters
  • mesh[in] Mesh

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

  • entity_list[in] List of entity indices (local)

  • orient[in] If true, in 3D, reorients facets to have consistent normal direction

Returns

Indices in the geometry array for the mesh entity vertices, i.e. indices(i, j) is the position in the geometry array of the j-th vertex of the entity entity_list[i].

std::vector<std::int32_t> exterior_facet_indices(const Mesh &mesh)

Compute the indices (local) of all exterior facets. An exterior facet (co-dimension 1) is one that is connected globally to only one cell of co-dimension 0).

Parameters

mesh[in] Mesh

Returns

List of facet indices of exterior facets of the mesh

CellPartitionFunction create_cell_partitioner(const graph::partition_fn &partfn = &graph::partition_graph)

Create a function that computes destination rank for mesh cells in 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 Mesh &mesh, const xtl::span<const std::int32_t> &entities, int d0, int d1)

Compute incident indices.

Parameters
  • mesh[in] The mesh

  • 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 T>
class MeshTags
#include <MeshTags.h>

MeshTags associate values with mesh 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(const std::shared_ptr<const Mesh> &mesh, int dim, U &&indices, V &&values)

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

Parameters
  • mesh[in] The mesh on which the tags are associated

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

  • indices[in] std::vector<std::int32> of sorted and unique entity indices (indices local to the process)

  • values[in] std::vector<T> of values for each index in indices. The size must be equal to the size of indices.

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

inline const std::vector<std::int32_t> &indices() const

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

inline const std::vector<T> &values() const

Values attached to mesh entities.

inline int dim() const

Return topological dimension of tagged entities.

inline std::shared_ptr<const Mesh> mesh() const

Return mesh.

Public Members

std::string name = "mesh_tags"

Name.

class Geometry
#include <Geometry.h>

Geometry stores the geometry imposed on a mesh.

Public Functions

template<typename AdjacencyList32, typename Array, typename Vector64>
inline Geometry(const std::shared_ptr<const common::IndexMap> &index_map, AdjacencyList32 &&dofmap, const fem::CoordinateElement &element, Array &&x, int dim, Vector64 &&input_global_indices)

Constructor.

Parameters
  • index_map[in] Index map associated with the geometry dofmap

  • dofmap[in] The geometry (point) dofmap. For a cell, it gives the position in the point array of each local geometry node

  • element[in] The element that describes the cell geometry map

  • x[in] The point coordinates. It is a std::vector<double> and uses row-major storage. The shape is (num_points, 3).

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

  • input_global_indices[in] The ‘global’ input index of each point, commonly from a mesh input file. The type is std:vector<std::int64_t>.

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.

int dim() const

Return Euclidean dimension of coordinate system.

const graph::AdjacencyList<std::int32_t> &dofmap() const

DOF map.

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

Index map.

xtl::span<const double> x() const

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

Returns

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

xtl::span<double> x()

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

Returns

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

const fem::CoordinateElement &cmap() const

The element that describes the geometry map.

Returns

The coordinate/geometry element

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

Global user indices.

class Mesh
#include <Mesh.h>

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

Public Functions

template<typename Topology, typename Geometry>
inline Mesh(MPI_Comm comm, Topology &&topology, Geometry &&geometry)

Create a 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

Topology &topology()

Get mesh topology.

Returns

The topology object associated with the mesh.

const Topology &topology() const

Get mesh topology (const version)

Returns

The topology object associated with the mesh.

Topology &topology_mutable() const

Get mesh topology if one really needs the mutable version.

Returns

The topology object associated with the mesh.

Geometry &geometry()

Get mesh geometry.

Returns

The geometry object associated with the mesh

const Geometry &geometry() const

Get mesh geometry (const version)

Returns

The geometry object associated with the mesh

MPI_Comm comm() const

Mesh MPI communicator.

Returns

The communicator on which the mesh is distributed

Public Members

std::string name = "mesh"

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, there is no clear caching policy implemented and no way of discarding cached data.

Public Functions

Topology(MPI_Comm comm, CellType type)

Create empty mesh topology.

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

Return the topological dimension of the mesh.

void set_index_map(int dim, const std::shared_ptr<const common::IndexMap> &map)

Todo:

Merge with set_connectivity

Set the IndexMap for dimension dim

Warning

This is experimental and likely to change

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(int d0, int d1) const

Return connectivity from entities of dimension d0 to entities of dimension d1.

Parameters
  • d0[in]

  • d1[in]

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.

void set_connectivity(std::shared_ptr<graph::AdjacencyList<std::int32_t>> c, int d0, int d1)

Set connectivity for given pair of topological dimensions.

Todo:

Merge with set_index_map

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 permutation number to apply to a facet.

The permutations are numbered so that:

  • n % 2 gives the number of reflections to apply

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

Each column of the returned array represents a cell, and each row a facet of that cell.

Note

An exception is raised if the permutations have not been computed

Returns

The permutation number

CellType cell_type() const noexcept

Cell type.

Returns

Cell type that the topology is for

std::int32_t create_entities(int dim)

Create entities of given topological dimension.

Parameters

dim[in] Topological dimension

Returns

Number of newly created entities, returns -1 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

The communicator on which the topology is distributed

Public Members

std::vector<std::int64_t> original_cell_index

Original cell index.