DOLFINx 0.9.0
DOLFINx C++ interface
|
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 |
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) |
CellType | to_type (const std::string &cell) |
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) |
graph::AdjacencyList< int > | get_entity_vertices (CellType type, int dim) |
graph::AdjacencyList< int > | get_sub_entities (CellType type, int dim0, int dim1) |
int | cell_dim (CellType type) |
Return topological dimension of cell type. | |
int | cell_num_entities (CellType type, int dim) |
bool | is_simplex (CellType type) |
int | num_cell_vertices (CellType type) |
std::map< std::array< int, 2 >, std::vector< std::set< int > > > | cell_entity_closure (CellType cell_type) |
basix::cell::type | cell_type_to_basix_type (CellType celltype) |
Convert a cell type to a Basix cell type. | |
CellType | cell_type_from_basix_type (basix::cell::type celltype) |
Get a cell type from a Basix cell type. | |
template<std::floating_point T = double> | |
Mesh< T > | create_box (MPI_Comm comm, MPI_Comm subcomm, std::array< std::array< T, 3 >, 2 > p, std::array< std::int64_t, 3 > n, CellType celltype, CellPartitionFunction partitioner=nullptr) |
Create a uniform mesh::Mesh over rectangular prism spanned by the two points p . | |
template<std::floating_point T = double> | |
Mesh< T > | create_box (MPI_Comm comm, std::array< std::array< T, 3 >, 2 > p, std::array< std::int64_t, 3 > n, CellType celltype, const CellPartitionFunction &partitioner=nullptr) |
Create a uniform mesh::Mesh over rectangular prism spanned by the two points p . | |
template<std::floating_point T = double> | |
Mesh< T > | create_rectangle (MPI_Comm comm, std::array< std::array< T, 2 >, 2 > p, std::array< std::int64_t, 2 > n, CellType celltype, CellPartitionFunction partitioner, DiagonalType diagonal=DiagonalType::right) |
Create a uniform mesh::Mesh over the rectangle spanned by the two points p . | |
template<std::floating_point T = double> | |
Mesh< T > | create_rectangle (MPI_Comm comm, std::array< std::array< T, 2 >, 2 > p, std::array< std::int64_t, 2 > n, CellType celltype, DiagonalType diagonal=DiagonalType::right) |
Create a uniform mesh::Mesh over the rectangle spanned by the two points p . | |
template<std::floating_point T = double> | |
Mesh< T > | create_interval (MPI_Comm comm, std::int64_t n, std::array< T, 2 > p, mesh::GhostMode ghost_mode=mesh::GhostMode::none, CellPartitionFunction partitioner=nullptr) |
Interval mesh of the 1D line [a, b] . | |
template<typename U > | |
Geometry< typename std::remove_reference_t< typename U::value_type > > | create_geometry (const Topology &topology, const std::vector< fem::CoordinateElement< std::remove_reference_t< typename U::value_type > > > &elements, std::span< const std::int64_t > nodes, std::span< const std::int64_t > xdofs, const U &x, int dim, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> reorder_fn=nullptr) |
Build Geometry from input data. | |
template<typename U > | |
Geometry< typename std::remove_reference_t< typename U::value_type > > | create_geometry (const Topology &topology, const fem::CoordinateElement< std::remove_reference_t< typename U::value_type > > &element, std::span< const std::int64_t > nodes, std::span< const std::int64_t > xdofs, const U &x, int dim, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> reorder_fn=nullptr) |
Build Geometry from input data. | |
std::tuple< graph::AdjacencyList< std::int32_t >, std::vector< std::int64_t >, std::size_t, std::vector< std::int32_t > > | build_local_dual_graph (std::span< const CellType > celltypes, const std::vector< std::span< const std::int64_t > > &cells) |
Compute the local part of the dual graph (cell-cell connections via facets) and facets with only one attached cell. | |
graph::AdjacencyList< std::int64_t > | build_dual_graph (MPI_Comm comm, std::span< const CellType > celltypes, const std::vector< std::span< const std::int64_t > > &cells) |
Build distributed mesh dual graph (cell-cell connections via facets) from minimal mesh data. | |
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) |
Topology | create_topology (MPI_Comm comm, std::span< const std::int64_t > cells, std::span< const std::int64_t > original_cell_index, std::span< const int > ghost_owners, CellType cell_type, std::span< const std::int64_t > boundary_vertices) |
Create a mesh topology. | |
Topology | create_topology (MPI_Comm comm, const std::vector< CellType > &cell_type, std::vector< std::span< const std::int64_t > > cells, std::vector< std::span< const std::int64_t > > original_cell_index, std::vector< std::span< const int > > ghost_owners, std::span< const std::int64_t > boundary_vertices) |
Create a topology of mixed cell type. | |
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, std::span< const std::int32_t > entities) |
Get entity indices for entities defined by their vertices. | |
std::tuple< std::vector< std::shared_ptr< graph::AdjacencyList< std::int32_t > > >, std::shared_ptr< graph::AdjacencyList< std::int32_t > >, std::shared_ptr< common::IndexMap >, std::vector< std::int32_t > > | compute_entities (MPI_Comm comm, const Topology &topology, int dim, int index) |
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, std::pair< std::int8_t, std::int8_t > d0, std::pair< std::int8_t, std::int8_t > d1) |
Compute connectivity (d0 -> d1) for given pair of entity types, given by topological dimension and index, as found in Topology::entity_types() | |
std::vector< std::int32_t > | exterior_facet_indices (const Topology &topology) |
Compute the indices of all exterior facets that are owned by the caller. | |
std::vector< std::int64_t > | extract_topology (CellType cell_type, const fem::ElementDofLayout &layout, std::span< const std::int64_t > cells) |
Extract topology from cell data, i.e. extract cell vertices. | |
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 permute=false) |
Compute the geometry degrees of freedom associated with the closure of a given set of cell entities. | |
CellPartitionFunction | create_cell_partitioner (mesh::GhostMode ghost_mode=mesh::GhostMode::none, const graph::partition_fn &partfn=&graph::partition_graph) |
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, MPI_Comm commt, std::span< const std::int64_t > cells, const fem::CoordinateElement< typename std::remove_reference_t< typename U::value_type > > &element, MPI_Comm commg, const U &x, std::array< std::size_t, 2 > xshape, const CellPartitionFunction &partitioner) |
Create a distributed mesh from mesh data using a provided graph partitioning function for determining the parallel distribution of the mesh. | |
template<typename U > | |
Mesh< typename std::remove_reference_t< typename U::value_type > > | create_mesh (MPI_Comm comm, MPI_Comm commt, const std::vector< std::span< const std::int64_t > > &cells, const std::vector< fem::CoordinateElement< typename std::remove_reference_t< typename U::value_type > > > &elements, MPI_Comm commg, const U &x, std::array< std::size_t, 2 > xshape, const CellPartitionFunction &partitioner) |
Create a distributed mixed-topology mesh from mesh data using a provided graph partitioning function for determining the parallel distribution of the mesh. | |
template<typename U > | |
Mesh< typename std::remove_reference_t< typename U::value_type > > | create_mesh (MPI_Comm comm, std::span< const std::int64_t > cells, const fem::CoordinateElement< std::remove_reference_t< typename U::value_type > > &elements, const U &x, std::array< std::size_t, 2 > xshape, GhostMode ghost_mode) |
Create a distributed mesh from mesh data using the default graph partitioner to determine the parallel distribution of the mesh. | |
template<std::floating_point T> | |
std::pair< Geometry< T >, std::vector< int32_t > > | create_subgeometry (const Mesh< T > &mesh, int dim, std::span< const std::int32_t > subentity_to_entity) |
Create a sub-geometry from a mesh and a subset of mesh entities to be included. A sub-geometry is simply a Geometry object containing only the geometric information for the subset of entities. The entities may differ in topological dimension from the original 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. | |
Mesh data structures and algorithms on meshes.
Representations of meshes and support for operations on meshes.
using CellPartitionFunction |
Signature for the cell partitioning function. The function should compute the destination rank for cells currently on this rank.
[in] | comm | MPI Communicator |
[in] | nparts | Number of partitions |
[in] | cell_types | Cell types in the mesh |
[in] | cells | Lists of cells of each cell type. cells[i] is a flattened row major 2D array of shape (num_cells, num_cell_vertices) for cell_types[i] on this process, containing the global indices for the cell vertices. Each cell can appear only once across all processes. The cell vertex indices are not necessarily contiguous globally, i.e. the maximum index across all processes can be greater than the number of vertices. High-order 'nodes', e.g. mid-side points, should not be included. |
graph::AdjacencyList< std::int64_t > build_dual_graph | ( | MPI_Comm | comm, |
std::span< const CellType > | celltypes, | ||
const std::vector< std::span< const std::int64_t > > & | cells ) |
Build distributed mesh dual graph (cell-cell connections via facets) from minimal mesh data.
The computed dual graph is typically passed to a graph partitioner.
[in] | comm | The MPI communicator |
[in] | celltypes | List of cell types |
[in] | cells | Collections of cells, defined by the cell vertices from which to build the dual graph, as flattened arrays for each cell type in celltypes . |
cells
and celltypes
must have the same size. std::tuple< graph::AdjacencyList< std::int32_t >, std::vector< std::int64_t >, std::size_t, std::vector< std::int32_t > > build_local_dual_graph | ( | std::span< const CellType > | celltypes, |
const std::vector< std::span< const std::int64_t > > & | cells ) |
Compute the local part of the dual graph (cell-cell connections via facets) and facets with only one attached cell.
[in] | celltypes | List of cell types. |
[in] | cells | Lists of cell vertices (stored as flattened lists, one for each cell type). |
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).
n
cells of type 0
and m
cells of type 1
, then cells of type 0
are numbered 0..(n-1)
and cells of type 1
are numbered n..(n+m-1)
respectively, in the returned dual graph. 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)}
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.
[in] | type | The cell type |
[in] | index | The facet index |
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)
(entities.size(), 3)
and the storage is row-major. int cell_num_entities | ( | CellType | type, |
int | dim ) |
Number of entities of dimension dim
[in] | dim | Entity dimension |
[in] | type | Cell type |
std::array< std::shared_ptr< graph::AdjacencyList< std::int32_t > >, 2 > compute_connectivity | ( | const Topology & | topology, |
std::pair< std::int8_t, std::int8_t > | d0, | ||
std::pair< std::int8_t, std::int8_t > | d1 ) |
Compute connectivity (d0 -> d1) for given pair of entity types, given by topological dimension and index, as found in Topology::entity_types()
[in] | topology | The topology |
[in] | d0 | The dimension and index of the entities |
[in] | d1 | The dimension and index of the incident entities |
std::tuple< std::vector< std::shared_ptr< graph::AdjacencyList< std::int32_t > > >, std::shared_ptr< graph::AdjacencyList< std::int32_t > >, std::shared_ptr< common::IndexMap >, std::vector< std::int32_t > > compute_entities | ( | MPI_Comm | comm, |
const Topology & | topology, | ||
int | dim, | ||
int | index ) |
Compute mesh entities of given topological dimension by computing entity-to-vertex connectivity (dim, 0)
, and cell-to-entity connectivity (tdim, dim)
.
Computed entities are oriented such that their local (to the process) orientation agrees with their global orientation
[in] | comm | MPI Communicator |
[in] | topology | Mesh topology |
[in] | dim | The dimension of the entities to create |
[in] | index | Index of entity in dimension dim as listed in Topology::entity_types(dim) . |
std::pair< std::vector< std::uint8_t >, std::vector< std::uint32_t > > compute_entity_permutations | ( | const Topology & | topology | ) |
Compute (1) facet rotation and reflection data, and (2) cell permutation data. This information is used in the assembly of (1) facet integrals and (2) most non-Lagrange elements.
The facet rotation and reflection data is encoded so that:
n % 2
gives the number of reflections to applyn // 2
gives the number of rotations to applyThe data is stored in a flattened 2D array, so that data[cell_index * / facets_per_cell + facet_index]
contains the facet with index facet_index
of the cell with index cell_index
. This data passed to FFCx kernels, where it is used to permute the quadrature points on facet integrals when data from the cells on both sides of the facet is used.
The cell permutation data contains information about the entities of each cell, relative to a low-to-high ordering. This data is packed so that a 32-bit int is used for each cell. For 2D cells, one bit is used for each edge, to represent whether or not the edge is reversed: the least significant bit is for edge 0, the next for edge 1, etc. For 3D cells, three bits are used for each face, and for each edge: the least significant bit says whether or not face 0 is reflected, the next 2 bits say how many times face 0 is rotated; the next three bits are for face 1, then three for face 2, etc; after all the faces, there is 1 bit for each edge to say whether or not they are reversed.
For example, if a quadrilateral has cell permutation info ....0111
then (from right to left):
and if a tetrahedron has cell permutation info ....011010010101001000
then (from right to left):
This data is used to correct the direction of vector function on permuted facets.
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.
[in] | topology | The topology |
[in] | entities | List of indices of topological dimension d0 |
[in] | d0 | Topological dimension |
[in] | d1 | Topological dimension |
d1
that are incident to entities in entities
(topological dimension d0
) 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.
(entities.size(), 3)
and the storage is row-major. Mesh< T > create_box | ( | MPI_Comm | comm, |
MPI_Comm | subcomm, | ||
std::array< std::array< T, 3 >, 2 > | p, | ||
std::array< std::int64_t, 3 > | n, | ||
CellType | celltype, | ||
CellPartitionFunction | partitioner = nullptr ) |
Create a uniform mesh::Mesh over rectangular prism spanned by the two points p
.
The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + / 1)*(n[1] + 1)*(n[2] + 1)
. For tetrahedra there will be will be 6*n[0]*n[1]*n[2]
cells. For hexahedra the number of cells will be n[0]*n[1]*n[2]
.
[in] | comm | MPI communicator to distribute the mesh on. |
[in] | subcomm | MPI communicator to construct and partition the mesh topology on. If the process should not be involved in the topology creation and partitioning then this communicator should be MPI_COMM_NULL . |
[in] | p | Corner of the box. |
[in] | n | Number of cells in each direction. |
[in] | celltype | Cell shape. |
[in] | partitioner | Partitioning function for distributing cells across MPI ranks. |
Mesh< T > create_box | ( | MPI_Comm | comm, |
std::array< std::array< T, 3 >, 2 > | p, | ||
std::array< std::int64_t, 3 > | n, | ||
CellType | celltype, | ||
const CellPartitionFunction & | partitioner = nullptr ) |
Create a uniform mesh::Mesh over rectangular prism spanned by the two points p
.
The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + / 1)*(n[1] + 1)*(n[2] + 1)
. For tetrahedra there will be will be 6*n[0]*n[1]*n[2]
cells. For hexahedra the number of cells will be n[0]*n[1]*n[2]
.
[in] | comm | MPI communicator to distribute the mesh on. |
[in] | p | Corner of the box. |
[in] | n | Number of cells in each direction. |
[in] | celltype | Cell shape. |
[in] | partitioner | Partitioning function for distributing cells across MPI ranks. |
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
Geometry< typename std::remove_reference_t< typename U::value_type > > create_geometry | ( | const Topology & | topology, |
const fem::CoordinateElement< std::remove_reference_t< typename U::value_type > > & | element, | ||
std::span< const std::int64_t > | nodes, | ||
std::span< const std::int64_t > | xdofs, | ||
const U & | x, | ||
int | dim, | ||
std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> | reorder_fn = nullptr ) |
Build Geometry from input data.
This function should be called after the mesh topology is built and 'node' coordinate data has been distributed to the processes where it is required.
[in] | topology | Mesh topology. |
[in] | element | Element that defines the geometry map for each cell. |
[in] | nodes | Geometry node global indices for cells on this process. Must be sorted. |
[in] | xdofs | Geometry degree-of-freedom map (using global indices) for cells on this process. nodes is a sorted and unique list of the indices in xdofs . |
[in] | x | 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. |
[in] | dim | Geometric dimension (1, 2, or 3). |
[in] | reorder_fn | Function for re-ordering the degree-of-freedom map associated with the geometry data. |
Geometry< typename std::remove_reference_t< typename U::value_type > > create_geometry | ( | const Topology & | topology, |
const std::vector< fem::CoordinateElement< std::remove_reference_t< typename U::value_type > > > & | elements, | ||
std::span< const std::int64_t > | nodes, | ||
std::span< const std::int64_t > | xdofs, | ||
const U & | x, | ||
int | dim, | ||
std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> | reorder_fn = nullptr ) |
Build Geometry from input data.
This function should be called after the mesh topology is built and 'node' coordinate data has been distributed to the processes where it is required.
[in] | topology | Mesh topology. |
[in] | elements | List of elements that defines the geometry map for each cell type. |
[in] | nodes | Geometry node global indices for cells on this process. |
[in] | xdofs | Geometry degree-of-freedom map (using global indices) for cells on this process. nodes is a sorted and unique list of the indices in xdofs . |
[in] | x | 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. |
[in] | dim | Geometric dimension (1, 2, or 3). |
[in] | reorder_fn | Function for re-ordering the degree-of-freedom map associated with the geometry data. |
Mesh< T > create_interval | ( | MPI_Comm | comm, |
std::int64_t | n, | ||
std::array< T, 2 > | p, | ||
mesh::GhostMode | ghost_mode = mesh::GhostMode::none, | ||
CellPartitionFunction | partitioner = nullptr ) |
Interval mesh of the 1D line [a, b]
.
Given n
cells in the axial direction, the total number of intervals will be n
and the total number of vertices will be n + 1
.
[in] | comm | MPI communicator to build the mesh on. |
[in] | n | Number of cells. |
[in] | p | End points of the interval. |
[in] | ghost_mode | ghost mode of the created mesh, defaults to none |
[in] | partitioner | Partitioning function for distributing cells across MPI ranks. |
Mesh< typename std::remove_reference_t< typename U::value_type > > create_mesh | ( | MPI_Comm | comm, |
MPI_Comm | commt, | ||
const std::vector< std::span< const std::int64_t > > & | cells, | ||
const std::vector< fem::CoordinateElement< typename std::remove_reference_t< typename U::value_type > > > & | elements, | ||
MPI_Comm | commg, | ||
const U & | x, | ||
std::array< std::size_t, 2 > | xshape, | ||
const CellPartitionFunction & | partitioner ) |
Create a distributed mixed-topology mesh from mesh data using a provided graph partitioning function for determining the parallel distribution of the mesh.
From mesh input data that is distributed across processes, a distributed mesh::Mesh is created. If the partitioning function is not callable, i.e. it does not store a callable function, no re-distribution of cells is done.
create_mesh
for mixed topology meshes, and does not include cell reordering.[in] | comm | Communicator to build the mesh on. |
[in] | commt | Communicator that the topology data (cells ) is distributed on. This should be MPI_COMM_NULL for ranks that should not participate in computing the topology partitioning. |
[in] | cells | Cells on the calling process, as a list of lists, one list for each cell type (or an empty list if there are no cells of that type on this process). The cells are defined by their 'nodes' (using global indices) following the Basix ordering, and concatenated to form a flattened list. For lowest order cells this will be just the cell vertices. For higher-order cells, other cells 'nodes' will be included. See dolfinx::io::cells for examples of the Basix ordering. |
[in] | elements | Coordinate elements for the cells, one for each cell type in the mesh. In parallel, these must be the same on all processes. |
[in] | commg | Communicator for geometry |
[in] | x | Geometry data ('node' coordinates). Row-major storage. The global index of the i th node (row) in x is taken as i plus the process offset oncomm , The offset is the sum of x rows on all processed with a lower rank than the caller. |
[in] | xshape | Shape of the x data. |
[in] | partitioner | Graph partitioner that computes the owning rank for each cell. If not callable, cells are not redistributed. |
comm
. Mesh< typename std::remove_reference_t< typename U::value_type > > create_mesh | ( | MPI_Comm | comm, |
MPI_Comm | commt, | ||
std::span< const std::int64_t > | cells, | ||
const fem::CoordinateElement< typename std::remove_reference_t< typename U::value_type > > & | element, | ||
MPI_Comm | commg, | ||
const U & | x, | ||
std::array< std::size_t, 2 > | xshape, | ||
const CellPartitionFunction & | partitioner ) |
Create a distributed mesh from mesh data using a provided graph partitioning function for determining the parallel distribution of the mesh.
From mesh input data that is distributed across processes, a distributed mesh::Mesh is created. If the partitioning function is not callable, i.e. it does not store a callable function, no re-distribution of cells is done.
[in] | comm | Communicator to build the mesh on. |
[in] | commt | Communicator that the topology data (cells ) is distributed on. This should be MPI_COMM_NULL for ranks that should not participate in computing the topology partitioning. |
[in] | cells | Cells on the calling process. Each cell (node in the AdjacencyList ) is defined by its 'nodes' (using global indices) following the Basix ordering. For lowest order cells this will be just the cell vertices. For higher-order cells, other cells 'nodes' will be included. See dolfinx::io::cells for examples of the Basix ordering. |
[in] | element | Coordinate element for the cells. |
[in] | commg | Communicator for geometry |
[in] | x | Geometry data ('node' coordinates). Row-major storage. The global index of the i th node (row) in x is taken as i plus the process offset oncomm , The offset is the sum of x rows on all processed with a lower rank than the caller. |
[in] | xshape | Shape of the x data. |
[in] | partitioner | Graph partitioner that computes the owning rank for each cell. If not callable, cells are not redistributed. |
comm
. Mesh< typename std::remove_reference_t< typename U::value_type > > create_mesh | ( | MPI_Comm | comm, |
std::span< const std::int64_t > | cells, | ||
const fem::CoordinateElement< std::remove_reference_t< typename U::value_type > > & | elements, | ||
const U & | x, | ||
std::array< std::size_t, 2 > | xshape, | ||
GhostMode | ghost_mode ) |
Create a distributed mesh from mesh data using the default graph partitioner to determine the parallel distribution of the mesh.
This function takes mesh input data that is distributed across processes and creates a mesh::Mesh, with the mesh cell distribution determined by the default cell partitioner. The default partitioner is based on graph partitioning.
[in] | comm | MPI communicator to build the mesh on. |
[in] | cells | Cells on the calling process. See create_mesh for a detailed description. |
[in] | elements | Coordinate elements for the cells. |
[in] | x | Geometry data ('node' coordinates). See create_mesh for a detailed description. |
[in] | xshape | The shape of x . It should be (num_points, gdim) . |
[in] | ghost_mode | The requested type of cell ghosting/overlap |
comm
. 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.
[in] | topology | Mesh topology that the tags are associated with. |
[in] | dim | Topological dimension of tagged entities. |
[in] | entities | Local vertex indices for tagged entities. |
[in] | values | Tag values for each entity in entities . The length of values must be equal to number of rows in entities . |
entities
must not contain duplicate entities. Mesh< T > create_rectangle | ( | MPI_Comm | comm, |
std::array< std::array< T, 2 >, 2 > | p, | ||
std::array< std::int64_t, 2 > | n, | ||
CellType | celltype, | ||
CellPartitionFunction | partitioner, | ||
DiagonalType | diagonal = DiagonalType::right ) |
Create a uniform mesh::Mesh over the rectangle spanned by the two points p
.
The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + / 1)*(n[1] + 1)
. For triangles there will be will be 2*n[0]*n[1]
cells. For quadrilaterals the number of cells will be n[0]*n[1]
.
[in] | comm | MPI communicator to build the mesh on. |
[in] | p | Bottom-left and top-right corners of the rectangle. |
[in] | n | Number of cells in each direction. |
[in] | celltype | Cell shape. |
[in] | partitioner | Partitioning function for distributing cells across MPI ranks. |
[in] | diagonal | Direction of diagonals |
Mesh< T > create_rectangle | ( | MPI_Comm | comm, |
std::array< std::array< T, 2 >, 2 > | p, | ||
std::array< std::int64_t, 2 > | n, | ||
CellType | celltype, | ||
DiagonalType | diagonal = DiagonalType::right ) |
Create a uniform mesh::Mesh over the rectangle spanned by the two points p
.
The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + / 1)*(n[1] + 1)
. For triangles there will be will be 2*n[0]*n[1]
cells. For quadrilaterals the number of cells will be n[0]*n[1]
.
[in] | comm | MPI communicator to build the mesh on |
[in] | p | Two corner points |
[in] | n | Number of cells in each direction |
[in] | celltype | Cell shape |
[in] | diagonal | Direction of diagonals |
std::pair< Geometry< T >, std::vector< int32_t > > create_subgeometry | ( | const Mesh< T > & | mesh, |
int | dim, | ||
std::span< const std::int32_t > | subentity_to_entity ) |
Create a sub-geometry from a mesh and a subset of mesh entities to be included. A sub-geometry is simply a Geometry
object containing only the geometric information for the subset of entities. The entities may differ in topological dimension from the original mesh.
[in] | mesh | The full mesh. |
[in] | dim | Topological dimension of the sub-topology. |
[in] | subentity_to_entity | Map from sub-topology entity to the entity in the parent topology. |
geometry
. 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.
[in] | mesh | The mesh |
[in] | dim | Entity dimension |
[in] | entities | List of entity indices in mesh to include in the new mesh |
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.
topology | Original (parent) topology. |
dim | Topological dimension of the entities in the new topology. |
entities | Indices of entities in topology to include in the new topology. |
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
. Topology create_topology | ( | MPI_Comm | comm, |
const std::vector< CellType > & | cell_type, | ||
std::vector< std::span< const std::int64_t > > | cells, | ||
std::vector< std::span< const std::int64_t > > | original_cell_index, | ||
std::vector< std::span< const int > > | ghost_owners, | ||
std::span< const std::int64_t > | boundary_vertices ) |
Create a topology of mixed cell type.
comm | MPI Communicator |
cell_type | List of cell types |
cells | Lists of cells, using vertex indices, flattened, for each cell type. |
original_cell_index | Input cell index for each cell type |
ghost_owners | Owning rank for ghost cells (at end of each list of cells). |
boundary_vertices | Vertices of undetermined ownership on external or inter-process boundary. |
Topology create_topology | ( | MPI_Comm | comm, |
std::span< const std::int64_t > | cells, | ||
std::span< const std::int64_t > | original_cell_index, | ||
std::span< const int > | ghost_owners, | ||
CellType | cell_type, | ||
std::span< const std::int64_t > | boundary_vertices ) |
Create a mesh topology.
This function creates a Topology from cells that have been distributed to the processes that own or ghost the cell.
[in] | comm | Communicator across which the topology is distributed. |
[in] | cells | Cell topology (list of vertices for each cell) using global indices for the vertices. It contains cells that have been distributed to this rank, e.g. via a graph partitioner. It must also contain all ghost cells via facet, i.e. cells that are on a neighboring process and which share a facet with a local cell. Ghost cells are the last n entries in cells , where n is given by the length of ghost_owners . |
[in] | original_cell_index | Original global index associated with each cell. |
[in] | ghost_owners | Owning rank of each ghost cell (ghost cells are always at the end of the list of cells ). |
[in] | cell_type | A vector with cell shapes. |
[in] | boundary_vertices | Vertices on the 'exterior' (boundary) of the local topology. These vertices might appear on other processes. |
std::vector< std::int32_t > entities_to_geometry | ( | const Mesh< T > & | mesh, |
int | dim, | ||
std::span< const std::int32_t > | entities, | ||
bool | permute = false ) |
Compute the geometry degrees of freedom associated with the closure of a given set of cell entities.
[in] | mesh | The mesh. |
[in] | dim | Topological dimension of the entities of interest. |
[in] | entities | Entity indices (local to process). |
[in] | permute | If true , permute the DOFs such that they are consistent with the orientation of dim -dimensional mesh entities. This requires create_entity_permutations to be called first. |
entities
. The shape is (num_entities, num_xdofs_per_entity)
and the storage is row-major. The index indices[i, j]
is the position in the geometry array of the j
-th vertex of the entity[i]
.dim -> mesh.topology().dim()
and mesh.topology().dim() -> dim
must have been computed. Otherwise an exception is thrown. std::vector< std::int32_t > entities_to_index | ( | const Topology & | topology, |
int | dim, | ||
std::span< const std::int32_t > | entities ) |
Get entity indices for entities defined by their vertices.
[in] | topology | The mesh topology |
[in] | dim | Topological dimension of the entities |
[in] | entities | The mesh entities defined by their vertices |
entities
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).
[in] | topology | Mesh topology |
std::vector< std::int64_t > extract_topology | ( | CellType | cell_type, |
const fem::ElementDofLayout & | layout, | ||
std::span< const std::int64_t > | cells ) |
Extract topology from cell data, i.e. extract cell vertices.
[in] | cell_type | The cell shape |
[in] | layout | The layout of geometry 'degrees-of-freedom' on the reference cell |
[in] | cells | List of 'nodes' for each cell using global indices. The layout must be consistent with layout . |
cell
. 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
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
).
[in] | mesh | Mesh that the entities belong to. |
[in] | entities | Indices (local to process) of entities to compute h for. |
[in] | dim | Topological dimension of the entities. |
h[i]
corresponds to the entity entities[i]
. bool is_simplex | ( | CellType | type | ) |
Check if cell is a simplex
[in] | type | Cell type |
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.
[in] | mesh | Mesh to mark entities on. |
[in] | dim | Topological dimension of the entities to be considered. |
[in] | marker | Marking function, returns true for a point that is 'marked', and false otherwise. |
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.
[in] | mesh | Mesh to mark entities on. |
[in] | dim | Topological dimension of the entities to be considered. Must be less than the topological dimension of the mesh. |
[in] | marker | Marking function, returns true for a point that is 'marked', and false otherwise. |
int num_cell_vertices | ( | CellType | type | ) |
Number vertices for a cell type
[in] | type | Cell type |
std::string to_string | ( | CellType | type | ) |
Get the cell string type for a cell type
[in] | type | The cell type |
mesh::CellType to_type | ( | const std::string & | cell | ) |
Get the cell type from a cell string
[in] | cell | Cell shape string |