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, int tdim, const graph::AdjacencyList<std::int64_t> &cells)>
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
-
enumerator point
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_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, const std::array<std::array<double, 3>, 2> &p, std::array<std::size_t, 3> n, CellType celltype, mesh::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 be6*n[0]*n[1]*n[2]
cells. For hexahedra the number of cells will ben[0]*n[1]*n[2]
.- Parameters:
comm – [in] MPI communicator to build 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:
-
template<std::floating_point T = double>
Mesh<T> create_rectangle(MPI_Comm comm, const std::array<std::array<double, 2>, 2> &p, std::array<std::size_t, 2> n, CellType celltype, 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 be2*n[0]*n[1]
cells. For quadrilaterals the number of cells will ben[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:
-
template<std::floating_point T = double>
Mesh<T> create_rectangle(MPI_Comm comm, const std::array<std::array<double, 2>, 2> &p, std::array<std::size_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 be2*n[0]*n[1]
cells. For quadrilaterals the number of cells will ben[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:
-
template<std::floating_point T = double>
Mesh<T> create_interval(MPI_Comm comm, std::size_t nx, std::array<double, 2> x, CellPartitionFunction partitioner = nullptr) Interval mesh of the 1D line
[a, b]
.Given
n
cells in the axial direction, the total number of intervals will ben
and the total number of vertices will ben + 1
.- Parameters:
comm – [in] MPI communicator to build the mesh on
nx – [in] The number of cells
x – [in] The end points of the interval
partitioner – [in] Partitioning function for distributing cells across MPI ranks.
- Returns:
A mesh
-
template<typename U>
mesh::Geometry<typename std::remove_reference_t<typename U::value_type>> create_geometry(MPI_Comm comm, const Topology &topology, const std::vector<fem::CoordinateElement<std::remove_reference_t<typename U::value_type>>> &elements, const graph::AdjacencyList<std::int64_t> &cell_nodes, 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. 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
elements – [in] The elements that defines the geometry map for each cell
cell_nodes – [in] The mesh cells, including higher-order geometry ‘nodes’
x – [in] The node coordinates (row-major, with shape
(num_nodes, dim)
. The global index of each node isi + rank_offset
, wherei
is the local row index inx
andrank_offset
is the sum ofx
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
-
template<std::floating_point T>
std::pair<mesh::Geometry<T>, std::vector<int32_t>> create_subgeometry(const Topology &topology, const Geometry<T> &geometry, int dim, std::span<const std::int32_t> subentity_to_entity) Create a sub-geometry for a subset of entities.
- Parameters:
topology – Full mesh topology
geometry – Full mesh geometry
dim – Topological dimension of the sub-topology
subentity_to_entity – 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
.
-
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 std::int64_t> cells, std::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]
, wherev_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
incell_vertices
is offsets[i]
tdim – [in] The topological dimension of the cells
- Returns:
Local dual graph
Facets, defined by their vertices, that are shared by only one cell on this rank. The logically 2D is array flattened (row-major).
The number of columns for the facet data array (2).
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
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 ofvalues
must be equal to number of rows inentities
.
-
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.
Get the permutation numbers to apply to facets. The permutations are numbered so that:
n % 2
gives the number of reflections to applyn // 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.
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
-
Topology create_topology(MPI_Comm comm, const graph::AdjacencyList<std::int64_t> &cells, std::span<const std::int64_t> original_cell_index, std::span<const int> ghost_owners, const std::vector<CellType> &cell_type, const std::vector<std::int32_t> &cell_group_offsets, std::span<const std::int64_t> boundary_vertices)
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] A vector with cell shapes
cell_group_offsets – [in] vector with each group offset, including ghosts.
boundary_vertices – [in] List of vertices on the exterior of the local mesh which may be shared with 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 – Original topology.
dim – Topological dimension of the entities in the new topology.
entities – The indices of the entities in
topology
to include in the new topology.
- Returns:
New topology of dimension
dim
with all entities inentities
, map from entities of dimensiondim
in new sub-topology to entities intopology
, and map from vertices in new sub-topology to vertices intopology
.
-
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>, std::vector<std::int32_t>> 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, list of interprocess entities). Interprocess entities lie on the “true” boundary between owned cells of each process. If the entities already exist, then {nullptr, nullptr, nullptr, std::vector()} 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.
-
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.
-
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
.
-
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 entityentities[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’, andfalse
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’, andfalse
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 orient) Determine the indices in the geometry data for each vertex of the given mesh entities.
Warning
This function should not be used unless there is no alternative. It may be removed in the future.
- Parameters:
mesh – [in] The mesh
dim – [in] Topological dimension of the entities of interest
entities – [in] Entity indices (local) to compute the vertex geometry indices for
orient – [in] If true, in 3D, reorients facets to have consistent normal direction
- Returns:
Indices in the geometry array for the entity vertices. The shape is
(num_entities, num_vertices_per_entity)
and the storage is row-major. The indexindices[i, j]
is the position in the geometry array of thej
-th vertex of theentity[i]
.
-
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 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 Topology &topology, std::span<const std::int32_t> entities, int d0, int d1)
Compute incident indices.
- 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 inentities
(topological dimensiond0
)
-
template<typename U>
Mesh<typename std::remove_reference_t<typename U::value_type>> create_mesh(MPI_Comm comm, const graph::AdjacencyList<std::int64_t> &cells, const std::vector<fem::CoordinateElement<typename std::remove_reference_t<typename U::value_type>>> &elements, const U &x, std::array<std::size_t, 2> xshape, CellPartitionFunction partitioner) Create a mesh using a provided mesh partitioning function.
If the partitioning function is not callable, i.e. it does not store a callable function, no re-distribution of cells is done.
-
template<typename U>
Mesh<typename std::remove_reference_t<typename U::value_type>> create_mesh(MPI_Comm comm, const graph::AdjacencyList<std::int64_t> &cells, const std::vector<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 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 on 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.elements – [in] The coordinate elements that describe the geometric mapping for cells.
x – [in] The coordinates of mesh nodes.
xshape – [in] The shape of
x
. It should be(num_points, gdim)
.ghost_mode – [in] The requested type of cell ghosting/overlap
- Returns:
A distributed Mesh.
-
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] Entity dimension
entities – [in] List of entity indices 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<std::floating_point T>
class Geometry - #include <Geometry.h>
Geometry stores the geometry imposed on a mesh.
Public Functions
- template<typename U, typename V, typename W> inline and std::is_convertible_v< std::remove_cvref_t< V >, std::vector< T > > and std::is_convertible_v< std::remove_cvref_t< W >, std::vector< std::int64_t > > Geometry (std::shared_ptr< const common::IndexMap > index_map, U &&dofmap, 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.
- 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
elements – [in] The elements that describes the cell geometry maps
x – [in] The 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] The ‘global’ input index of each point, commonly from a mesh input file.
-
~Geometry() = default
Destructor.
-
inline int dim() const
Return Euclidean dimension of coordinate system.
-
inline MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<const std::int32_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>> dofmap() const
DOF map.
-
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 maps.
- Returns:
The coordinate/geometry elements
-
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:
The – floating point type for representing the geometry.
Public Functions
Create a mesh.
-
~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
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 and std::is_convertible_v< std::remove_cvref_t< V >, std::vector< T > > 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() = default
Destructor.
-
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 int dim() const
Return topological dimension of tagged entities.
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, there is no clear caching policy implemented and no way of discarding cached data.
Public Functions
-
~Topology() = default
Destructor.
-
int dim() const noexcept
Return the topological dimension of the mesh.
-
void set_entity_group_offsets(int dim, const std::vector<std::int32_t> &offsets)
Set the offsets for each group of entities of a particular dimension. See
entity_group_offsets
.- Parameters:
dim – Dimension of the entities
offsets – The offsets
-
const std::vector<std::int32_t> &entity_group_offsets(int dim) const
Get the offsets for each group of entities of a particular dimension.
The topology may consist of more than one cell type or facet type. In that case, the cells of the same types are grouped together in blocks, firstly for regular cells, then repeated for ghost cells. For example, a mesh with two triangles, three quads and no ghosts would have offsets: 0, 2, 5, 5, 5. A mesh with twenty tetrahedra and two ghost tetrahedra has offsets: 0, 20, 22.
- Parameters:
dim – Dimension of the entities
- Returns:
The offsets
Set the IndexMap for dimension dim.
- Todo:
Merge with set_connectivity
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
. Returnsnullptr
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.
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 applyn // 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
-
std::vector<CellType> cell_types() const noexcept
Cell type.
- Returns:
Cell types 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.
-
const std::vector<std::int32_t> &interprocess_facets() const
List of inter-process facets, if facet topology has been computed.
Public Members
-
std::vector<std::int64_t> original_cell_index
Original cell index.
-
namespace impl
Functions
-
template<std::floating_point T>
Mesh<T> build_tri(MPI_Comm comm, const std::array<std::array<double, 2>, 2> &p, std::array<std::size_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<double, 2>, 2> p, std::array<std::size_t, 2> n, const CellPartitionFunction &partitioner)
-
template<std::floating_point T>
std::vector<T> create_geom(MPI_Comm comm, const std::array<std::array<double, 3>, 2> &p, std::array<std::size_t, 3> n)
-
template<std::floating_point T>
Mesh<T> build_tet(MPI_Comm comm, const std::array<std::array<double, 3>, 2> &p, std::array<std::size_t, 3> n, const CellPartitionFunction &partitioner)
-
template<std::floating_point T>
Mesh<T> build_hex(MPI_Comm comm, const std::array<std::array<double, 3>, 2> &p, std::array<std::size_t, 3> n, const CellPartitionFunction &partitioner)
-
template<std::floating_point T>
Mesh<T> build_prism(MPI_Comm comm, const std::array<std::array<double, 3>, 2> &p, std::array<std::size_t, 3> n, const CellPartitionFunction &partitioner)
-
template<std::floating_point T>
mesh::Mesh<T> build_hex(MPI_Comm comm, const std::array<std::array<double, 3>, 2> &p, std::array<std::size_t, 3> n, const CellPartitionFunction &partitioner)
-
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) The coordinates of ‘vertices’ for for entities of a give 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 on the meh boundary
- Pre:
The provided facets must be on the boundary of the mesh.
- Returns:
(0) Entities attached to the boundary facets, (1) vertex coordinates (shape is
(3, num_vertices)
) and (2) map from vertex in the full mesh to the position (column) 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.
-
template<std::floating_point T>
-
using CellPartitionFunction = std::function<graph::AdjacencyList<std::int32_t>(MPI_Comm comm, int nparts, int tdim, const graph::AdjacencyList<std::int64_t> &cells)>