Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/dd/d7d/namespacedolfinx_1_1mesh.html
DOLFINx 0.7.3
DOLFINx C++ interface
Loading...
Searching...
No Matches
Classes | Concepts | Typedefs | Enumerations | Functions
dolfinx::mesh Namespace Reference

Mesh data structures and algorithms on meshes. More...

Classes

class  Geometry
 Geometry stores the geometry imposed on a mesh. More...
 
class  Mesh
 A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data. More...
 
class  MeshTags
 MeshTags associate values with mesh topology entities. More...
 
class  Topology
 Topology stores the topology of a mesh, consisting of mesh entities and connectivity (incidence relations for the mesh entities). More...
 

Concepts

concept  MarkerFn
 Requirements on function for geometry marking.
 

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.
 

Enumerations

enum class  CellType : int {
  point = 1 , interval = 2 , triangle = 3 , tetrahedron = 4 ,
  quadrilateral = -4 , pyramid = -5 , prism = -6 , hexahedron = -8
}
 Cell type identifier.
 
enum class  DiagonalType {
  left , right , crossed , shared_facet ,
  left_right , right_left
}
 Enum for different diagonal types.
 
enum class  GhostMode : int { none , shared_facet , shared_vertex }
 Enum for different partitioning ghost modes.
 

Functions

std::string to_string (CellType type)
 Get the cell string type for a cell type.
 
CellType to_type (const std::string &cell)
 Get the cell type from a cell string.
 
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.
 
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.
 
bool is_simplex (CellType type)
 Check if cell is a simplex.
 
int num_cell_vertices (CellType type)
 Number vertices for a cell type.
 
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.
 
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.
 
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.
 
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].
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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).
 
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.
 
std::vector< std::int32_t > exterior_facet_indices (const Topology &topology)
 Compute the indices of all exterior facets that are owned by the caller.
 
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.
 
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).
 
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)
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 

Detailed Description

Mesh data structures and algorithms on meshes.

Representations of meshes and support for operations on meshes.

Typedef Documentation

◆ CellPartitionFunction

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.

Parameters
[in]commMPI Communicator
[in]npartsNumber of partitions
[in]tdimTopological dimension
[in]cellsCells 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.
[in]ghost_modeHow to overlap the cell partitioning: none, shared_facet or shared_vertex
Returns
Destination ranks for each cell on this process

Function Documentation

◆ build_dual_graph()

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
[in]commThe MPI communicator
[in]cellsCollection of cells, defined by the cell vertices from which to build the dual graph
[in]tdimThe topological dimension of the cells
Returns
The dual graph

◆ build_local_dual_graph()

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.

Parameters
[in]cellsArray for cell vertices adjacency list
[in]offsetsAdjacency list offsets, i.e. index of the first entry of cell i in cell_vertices is offsets[i]
[in]tdimThe 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).

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.

◆ cell_facet_type()

mesh::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
[in]typeThe cell type
[in]indexThe facet index
Returns
The type of facet for this cell at this index

◆ cell_normals()

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.

◆ cell_num_entities()

int cell_num_entities ( CellType  type,
int  dim 
)

Number of entities of dimension dim.

Parameters
[in]dimEntity dimension
[in]typeCell type
Returns
Number of entities in cell

◆ compute_connectivity()

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
[in]topologyThe topology
[in]d0The dimension of the nodes in the adjacency list
[in]d1The 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.

◆ compute_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
[in]commMPI Communicator
[in]topologyMesh topology
[in]dimThe 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.

◆ compute_entity_permutations()

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

◆ compute_incident_entities()

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
[in]topologyThe topology
[in]entitiesList of indices of topological dimension d0
[in]d0Topological dimension
[in]d1Topological dimension
Returns
List of entities of topological dimension d1 that are incident to entities in entities (topological dimension d0)

◆ compute_midpoints()

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.

◆ create_box()

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 be 6*n[0]*n[1]*n[2] cells. For hexahedra the number of cells will be n[0]*n[1]*n[2].

Parameters
[in]commMPI communicator to build mesh on.
[in]pCorner of the box.
[in]nNumber of cells in each direction.
[in]celltypeCell shape.
[in]partitionerPartitioning function for distributing cells across MPI ranks.
Returns
Mesh

◆ create_cell_partitioner()

mesh::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

◆ create_geometry()

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
[in]commThe MPI communicator to build the Geometry on
[in]topologyThe mesh topology
[in]elementsThe elements that defines the geometry map for each cell
[in]cell_nodesThe mesh cells, including higher-order geometry 'nodes'
[in]xThe 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.
[in]dimThe geometric dimension (1, 2, or 3)
[in]reorder_fnFunction for re-ordering the degree-of-freedom map associated with the geometry data

◆ create_interval()

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 be n and the total number of vertices will be n + 1.

Parameters
[in]commMPI communicator to build the mesh on
[in]nxThe number of cells
[in]xThe end points of the interval
[in]partitionerPartitioning function for distributing cells across MPI ranks.
Returns
A mesh

◆ create_mesh() [1/2]

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
[in]commThe MPI communicator to build the mesh on.
[in]cellsThe 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.
[in]elementsThe coordinate elements that describe the geometric mapping for cells.
[in]xThe coordinates of mesh nodes.
[in]xshapeThe shape of x. It should be (num_points, gdim).
[in]ghost_modeThe requested type of cell ghosting/overlap
Returns
A distributed Mesh.

◆ create_mesh() [2/2]

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.

◆ create_meshtags()

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.

Parameters
[in]topologyMesh topology that the tags are associated with
[in]dimTopological dimension of tagged entities
[in]entitiesLocal vertex indices for tagged entities.
[in]valuesTag values for each entity in entities. The length of values must be equal to number of rows in entities.
Note
Entities that do not exist on this rank are ignored.
Warning
entities must not contain duplicate entities.

◆ create_rectangle() [1/2]

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 be 2*n[0]*n[1] cells. For quadrilaterals the number of cells will be n[0]*n[1].

Parameters
[in]commMPI communicator to build the mesh on.
[in]pBottom-left and top-right corners of the rectangle.
[in]nNumber of cells in each direction.
[in]celltypeCell shape.
[in]partitionerPartitioning function for distributing cells across MPI ranks.
[in]diagonalDirection of diagonals
Returns
Mesh

◆ create_rectangle() [2/2]

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 be 2*n[0]*n[1] cells. For quadrilaterals the number of cells will be n[0]*n[1].

Parameters
[in]commMPI communicator to build the mesh on
[in]pTwo corner points
[in]nNumber of cells in each direction
[in]celltypeCell shape
[in]diagonalDirection of diagonals
Returns
Mesh

◆ create_subgeometry()

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
topologyFull mesh topology
geometryFull mesh geometry
dimTopological dimension of the sub-topology
subentity_to_entityMap 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.

◆ create_submesh()

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
[in]meshThe mesh
[in]dimEntity dimension
[in]entitiesList 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.

◆ create_subtopology()

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
topologyOriginal topology.
dimTopological dimension of the entities in the new topology.
entitiesThe indices of the 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.

◆ create_topology()

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
[in]commMPI communicator across which the topology is distributed
[in]cellsThe 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.
[in]original_cell_indexThe original global index associated with each cell
[in]ghost_ownersThe owning rank of each ghost cell (ghost cells are always at the end of the list of cells)
[in]cell_typeA vector with cell shapes
[in]cell_group_offsetsvector with each group offset, including ghosts.
[in]boundary_verticesList of vertices on the exterior of the local mesh which may be shared with other processes.
Returns
A distributed mesh topology

◆ entities_to_geometry()

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
[in]meshThe mesh
[in]dimTopological dimension of the entities of interest
[in]entitiesEntity indices (local) to compute the vertex geometry indices for
[in]orientIf 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 index indices[i, j] is the position in the geometry array of the j-th vertex of the entity[i].

◆ entities_to_index()

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.

Warning
This function may be removed in the future.
Parameters
[in]topologyThe mesh topology
[in]dimTopological dimension of the entities
[in]entitiesThe mesh entities defined by their vertices
Returns
The index of the ith entity in entities
Note
If an entity cannot be found on this rank, -1 is returned as the index.

◆ exterior_facet_indices()

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
[in]topologyMesh topology
Returns
Sorted list of owned facet indices that are exterior facets of the mesh.

◆ extract_topology()

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
[in]cell_typeThe cell shape
[in]layoutThe layout of geometry 'degrees-of-freedom' on the reference cell
[in]cellsList 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.

◆ h()

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
[in]meshMesh that the entities belong to.
[in]entitiesIndices (local to process) of entities to compute h for.
[in]dimTopological dimension of the entities.
Returns
Greatest distance between any two vertices, h[i] corresponds to the entity entities[i].

◆ is_simplex()

bool is_simplex ( CellType  type)

Check if cell is a simplex.

Parameters
[in]typeCell type
Returns
True is the cell type is a simplex

◆ locate_entities()

template<std::floating_point T, MarkerFn< T > U>
std::vector< std::int32_t > locate_entities ( const Mesh< T > &  mesh,
int  dim,
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
[in]meshMesh to mark entities on.
[in]dimTopological dimension of the entities to be considered.
[in]markerMarking 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)

◆ locate_entities_boundary()

template<std::floating_point T, MarkerFn< T > U>
std::vector< std::int32_t > locate_entities_boundary ( const Mesh< T > &  mesh,
int  dim,
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
[in]meshMesh to mark entities on.
[in]dimTopological dimension of the entities to be considered. Must be less than the topological dimension of the mesh.
[in]markerMarking function, returns true for a point that is 'marked', and false otherwise.
Returns
List of marked entity indices (indices local to the process)

◆ num_cell_vertices()

int num_cell_vertices ( CellType  type)

Number vertices for a cell type.

Parameters
[in]typeCell type
Returns
The number of cell vertices

◆ to_string()

std::string to_string ( CellType  type)

Get the cell string type for a cell type.

Parameters
[in]typeThe cell type
Returns
The cell type string

◆ to_type()

mesh::CellType to_type ( const std::string &  cell)

Get the cell type from a cell string.

Parameters
[in]cellCell shape string
Returns
The cell type