DOLFINx
0.5.1
DOLFINx C++ interface
|
Mesh data structures and algorithms on meshes. More...
Classes | |
class | MeshTags |
MeshTags associate values with mesh entities. More... | |
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 | Topology |
Topology stores the topology of a mesh, consisting of mesh entities and connectivity (incidence relations for the mesh entities). More... | |
Typedefs | |
using | CellPartitionFunction = std::function< dolfinx::graph::AdjacencyList< std::int32_t >(MPI_Comm comm, int nparts, int tdim, const dolfinx::graph::AdjacencyList< std::int64_t > &cells, dolfinx::mesh::GhostMode ghost_mode)> |
Signature for the cell partitioning function. The function should compute the destination rank for cells currently on this rank. More... | |
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. More... | |
CellType | to_type (const std::string &cell) |
Get the cell type from a cell string. More... | |
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. More... | |
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. More... | |
bool | is_simplex (CellType type) |
Check if cell is a simplex. More... | |
int | num_cell_vertices (CellType type) |
Number vertices for a cell type. More... | |
std::map< std::array< int, 2 >, std::vector< std::set< int > > > | cell_entity_closure (CellType cell_type) |
Closure entities for a cell, i.e., all lower-dimensional entities attached to a cell entity. Map from entity {dim_e, entity_e} to closure{sub_dim, (sub_entities)}. | |
basix::cell::type | cell_type_to_basix_type (CellType celltype) |
Convert a cell type to a Basix cell type. | |
CellType | cell_type_from_basix_type (basix::cell::type celltype) |
Get a cell type from a Basix cell type. | |
Mesh | create_box (MPI_Comm comm, const std::array< std::array< double, 3 >, 2 > &p, std::array< std::size_t, 3 > n, CellType celltype, GhostMode ghost_mode, const mesh::CellPartitionFunction &partitioner=create_cell_partitioner()) |
Create a uniform mesh::Mesh over the rectangular prism spanned by the two points p . The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + 1)*(n[1] + 1)*(n[2] + 1) . For tetrahedra there will be will be 6*n[0]*n[1]*n[2] cells. For hexahedra the number of cells will be n[0]*n[1]*n[2] . More... | |
Mesh | create_interval (MPI_Comm comm, std::size_t n, std::array< double, 2 > x, GhostMode ghost_mode, const CellPartitionFunction &partitioner=create_cell_partitioner()) |
Interval mesh of the 1D line [a, b] . Given n cells in the axial direction, the total number of intervals will be n and the total number of vertices will be n + 1 . More... | |
Mesh | create_rectangle (MPI_Comm comm, const std::array< std::array< double, 2 >, 2 > &p, std::array< std::size_t, 2 > n, CellType celltype, GhostMode ghost_mode, DiagonalType diagonal=DiagonalType::right) |
Create a uniform mesh::Mesh over the rectangle spanned by the two points p . The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + 1)*(n[1] + 1) . For triangles there will be will be 2*n[0]*n[1] cells. For quadrilaterals the number of cells will be n[0]*n[1] . More... | |
Mesh | create_rectangle (MPI_Comm comm, const std::array< std::array< double, 2 >, 2 > &p, std::array< std::size_t, 2 > n, CellType celltype, GhostMode ghost_mode, const CellPartitionFunction &partitioner, DiagonalType diagonal=DiagonalType::right) |
Create a uniform mesh::Mesh over the rectangle spanned by the two points p . The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + 1)*(n[1] + 1) . For triangles there will be will be 2*n[0]*n[1] cells. For quadrilaterals the number of cells will be n[0]*n[1] . More... | |
mesh::Geometry | create_geometry (MPI_Comm comm, const Topology &topology, const fem::CoordinateElement &element, const graph::AdjacencyList< std::int64_t > &cells, const std::span< const double > &x, int dim, const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &reorder_fn=nullptr) |
Build Geometry from input data. More... | |
std::tuple< graph::AdjacencyList< std::int32_t >, std::vector< std::int64_t >, std::size_t, std::vector< std::int32_t > > | build_local_dual_graph (const std::span< const std::int64_t > &cells, const 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. More... | |
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. More... | |
Mesh | create_mesh (MPI_Comm comm, const graph::AdjacencyList< std::int64_t > &cells, const fem::CoordinateElement &element, std::span< const double > x, std::array< std::size_t, 2 > xshape, GhostMode ghost_mode) |
Create a mesh using the default partitioner. More... | |
Mesh | create_mesh (MPI_Comm comm, const graph::AdjacencyList< std::int64_t > &cells, const fem::CoordinateElement &element, std::span< const double > x, std::array< std::size_t, 2 > xshape, GhostMode ghost_mode, const CellPartitionFunction &cell_partitioner) |
Create a mesh using a provided mesh partitioning function. | |
std::tuple< Mesh, std::vector< std::int32_t >, std::vector< std::int32_t >, std::vector< std::int32_t > > | create_submesh (const Mesh &mesh, int dim, const std::span< const std::int32_t > &entities) |
Create a new mesh consisting of a subset of entities in a mesh. More... | |
template<typename T > | |
MeshTags< T > | create_meshtags (const std::shared_ptr< const Mesh > &mesh, int dim, const graph::AdjacencyList< std::int32_t > &entities, const std::span< const T > &values) |
Create MeshTags from arrays. More... | |
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. More... | |
Topology | create_topology (MPI_Comm comm, const graph::AdjacencyList< std::int64_t > &cells, const std::span< const std::int64_t > &original_cell_index, const std::span< const int > &ghost_owners, const CellType &cell_type, GhostMode ghost_mode) |
Create a distributed mesh topology. More... | |
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. More... | |
std::tuple< std::shared_ptr< graph::AdjacencyList< std::int32_t > >, std::shared_ptr< graph::AdjacencyList< std::int32_t > >, std::shared_ptr< common::IndexMap > > | compute_entities (MPI_Comm comm, const Topology &topology, int dim) |
Compute mesh entities of given topological dimension by computing entity-to-vertex connectivity (dim, 0), and cell-to-entity connectivity (tdim, dim) More... | |
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. More... | |
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. More... | |
std::vector< double > | h (const Mesh &mesh, const std::span< const std::int32_t > &entities, int dim) |
Compute greatest distance between any two vertices of the mesh entities (h ). More... | |
std::vector< double > | cell_normals (const Mesh &mesh, int dim, const std::span< const std::int32_t > &entities) |
Compute normal to given cell (viewed as embedded in 3D) More... | |
std::vector< double > | compute_midpoints (const Mesh &mesh, int dim, const std::span< const std::int32_t > &entities) |
Compute the midpoints for mesh entities of a given dimension. More... | |
std::vector< std::int32_t > | locate_entities (const Mesh &mesh, int dim, const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> &marker) |
Compute indices of all mesh entities that evaluate to true for the provided geometric marking function. An entity is considered marked if the marker function evaluates true for all of its vertices. More... | |
std::vector< std::int32_t > | locate_entities_boundary (const Mesh &mesh, int dim, const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> &marker) |
Compute indices of all mesh entities that are attached to an owned boundary facet and evaluate to true for the provided geometric marking function. An entity is considered marked if the marker function evaluates true for all of its vertices. More... | |
std::vector< std::int32_t > | entities_to_geometry (const Mesh &mesh, int dim, const std::span< const std::int32_t > &entities, bool orient) |
Determine the indices in the geometry data for each vertex of the given mesh entities. More... | |
std::vector< std::int32_t > | exterior_facet_indices (const Topology &topology) |
Compute the indices of all exterior facets that are owned by the caller. More... | |
CellPartitionFunction | create_cell_partitioner (const graph::partition_fn &partfn=&graph::partition_graph) |
Create a function that computes destination rank for mesh cells in this rank by applying the default graph partitioner to the dual graph of the mesh. More... | |
std::vector< std::int32_t > | compute_incident_entities (const Mesh &mesh, const std::span< const std::int32_t > &entities, int d0, int d1) |
Compute incident indices. More... | |
Mesh data structures and algorithms on meshes.
Representations of meshes and support for operations on meshes.
using CellPartitionFunction = std::function<dolfinx::graph::AdjacencyList<std::int32_t>( MPI_Comm comm, int nparts, int tdim, const dolfinx::graph::AdjacencyList<std::int64_t>& cells, dolfinx::mesh::GhostMode ghost_mode)> |
Signature for the cell partitioning function. The function should compute the destination rank for cells currently on this rank.
[in] | comm | MPI Communicator |
[in] | nparts | Number of partitions |
[in] | tdim | Topological dimension |
[in] | cells | 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. |
[in] | ghost_mode | How to overlap the cell partitioning: none, shared_facet or shared_vertex |
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.
[in] | comm | The MPI communicator |
[in] | cells | Collection of cells, defined by the cell vertices from which to build the dual graph |
[in] | tdim | The topological dimension of the cells |
std::tuple< graph::AdjacencyList< std::int32_t >, std::vector< std::int64_t >, std::size_t, std::vector< std::int32_t > > build_local_dual_graph | ( | const std::span< const std::int64_t > & | cells, |
const 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.
[in] | cells | Array for cell vertices adjacency list |
[in] | offsets | Adjacency list offsets, i.e. index of the first entry of cell i in cell_vertices is offsets[i] |
[in] | tdim | The topological dimension of the cells |
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).
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< double > cell_normals | ( | const Mesh & | mesh, |
int | dim, | ||
const 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, |
int | d0, | ||
int | d1 | ||
) |
Compute connectivity (d0 -> d1) for given pair of topological dimensions.
[in] | topology | The topology |
[in] | d0 | The dimension of the nodes in the adjacency list |
[in] | d1 | The dimension of the edges in the adjacency list |
std::tuple< std::shared_ptr< graph::AdjacencyList< std::int32_t > >, std::shared_ptr< graph::AdjacencyList< std::int32_t > >, std::shared_ptr< common::IndexMap > > compute_entities | ( | MPI_Comm | comm, |
const Topology & | topology, | ||
int | dim | ||
) |
Compute mesh entities of given topological dimension by computing entity-to-vertex connectivity (dim, 0), and cell-to-entity connectivity (tdim, dim)
[in] | comm | MPI Communicator |
[in] | topology | Mesh topology |
[in] | dim | The dimension of the entities to create |
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 applyEach 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):
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 Mesh & | mesh, |
const std::span< const std::int32_t > & | entities, | ||
int | d0, | ||
int | d1 | ||
) |
Compute incident indices.
[in] | mesh | The mesh |
[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< double > compute_midpoints | ( | const Mesh & | mesh, |
int | dim, | ||
const 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::Mesh create_box | ( | MPI_Comm | comm, |
const std::array< std::array< double, 3 >, 2 > & | p, | ||
std::array< std::size_t, 3 > | n, | ||
CellType | celltype, | ||
GhostMode | ghost_mode, | ||
const mesh::CellPartitionFunction & | partitioner = create_cell_partitioner() |
||
) |
Create a uniform mesh::Mesh over the rectangular prism spanned by the two points p
. The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + 1)*(n[1] + 1)*(n[2] + 1)
. For tetrahedra there will be will be 6*n[0]*n[1]*n[2]
cells. For hexahedra the number of cells will be n[0]*n[1]*n[2]
.
[in] | comm | MPI communicator to build mesh on |
[in] | p | Points of box |
[in] | n | Number of cells in each direction |
[in] | celltype | Cell shape |
[in] | ghost_mode | Ghost mode |
[in] | partitioner | Partitioning function to use for determining the parallel distribution of cells across MPI ranks |
mesh::CellPartitionFunction create_cell_partitioner | ( | const graph::partition_fn & | partfn = &graph::partition_graph | ) |
Create a function that computes destination rank for mesh cells in this rank by applying the default graph partitioner to the dual graph of the mesh.
mesh::Geometry create_geometry | ( | MPI_Comm | comm, |
const Topology & | topology, | ||
const fem::CoordinateElement & | element, | ||
const graph::AdjacencyList< std::int64_t > & | cells, | ||
const std::span< const double > & | x, | ||
int | dim, | ||
const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> & | reorder_fn = nullptr |
||
) |
Build Geometry from input data.
This function should be called after the mesh topology is built. It distributes the 'node' coordinate data to the required MPI process and then creates a mesh::Geometry object.
[in] | comm | The MPI communicator to build the Geometry on |
[in] | topology | The mesh topology |
[in] | element | The element that defines the geometry map for each cell |
[in] | cells | The mesh cells, including higher-order geometry 'nodes' |
[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 | The geometric dimension (1, 2, or 3) |
[in] | reorder_fn | Function for re-ordering the degree-of-freedom map associated with the geometry data |
mesh::Mesh create_interval | ( | MPI_Comm | comm, |
std::size_t | n, | ||
std::array< double, 2 > | x, | ||
GhostMode | ghost_mode, | ||
const CellPartitionFunction & | partitioner = create_cell_partitioner() |
||
) |
Interval mesh of the 1D line [a, b]
. Given n
cells in the axial direction, the total number of intervals will be n
and the total number of vertices will be n + 1
.
[in] | comm | MPI communicator to build the mesh on |
[in] | n | The number of cells |
[in] | x | The end points of the interval |
[in] | ghost_mode | Ghosting mode |
[in] | partitioner | Partitioning function to use for determining the parallel distribution of cells across MPI ranks |
Mesh create_mesh | ( | MPI_Comm | comm, |
const graph::AdjacencyList< std::int64_t > & | cells, | ||
const fem::CoordinateElement & | element, | ||
std::span< const double > | x, | ||
std::array< std::size_t, 2 > | xshape, | ||
mesh::GhostMode | ghost_mode | ||
) |
Create a mesh using the default partitioner.
This function takes mesh input data that is distributed across processes and creates a mesh::Mesh, with the mesh cell distribution determined by the default cell partitioner. The default partitioner is based a graph partitioning.
[in] | comm | The MPI communicator to build the mesh on |
[in] | cells | 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. |
[in] | element | The coordinate element that describes the geometric mapping for cells |
[in] | x | The coordinates of mesh nodes |
[in] | xshape | The shape of x |
[in] | ghost_mode | The requested type of cell ghosting/overlap |
MeshTags<T> dolfinx::mesh::create_meshtags | ( | const std::shared_ptr< const Mesh > & | mesh, |
int | dim, | ||
const graph::AdjacencyList< std::int32_t > & | entities, | ||
const std::span< const T > & | values | ||
) |
Create MeshTags from arrays.
[in] | mesh | The Mesh 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 conatined duplicate entities. mesh::Mesh create_rectangle | ( | MPI_Comm | comm, |
const std::array< std::array< double, 2 >, 2 > & | p, | ||
std::array< std::size_t, 2 > | n, | ||
CellType | celltype, | ||
GhostMode | ghost_mode, | ||
const CellPartitionFunction & | partitioner, | ||
DiagonalType | diagonal = DiagonalType::right |
||
) |
Create a uniform mesh::Mesh over the rectangle spanned by the two points p
. The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + 1)*(n[1] + 1)
. For triangles there will be will be 2*n[0]*n[1]
cells. For quadrilaterals the number of cells will be n[0]*n[1]
.
[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] | ghost_mode | Mesh ghosting mode |
[in] | partitioner | Partitioning function to use for determining the parallel distribution of cells across MPI ranks |
[in] | diagonal | Direction of diagonals |
mesh::Mesh create_rectangle | ( | MPI_Comm | comm, |
const std::array< std::array< double, 2 >, 2 > & | p, | ||
std::array< std::size_t, 2 > | n, | ||
CellType | celltype, | ||
GhostMode | ghost_mode, | ||
DiagonalType | diagonal = DiagonalType::right |
||
) |
Create a uniform mesh::Mesh over the rectangle spanned by the two points p
. The order of the two points is not important in terms of minimum and maximum coordinates. The total number of vertices will be (n[0] + 1)*(n[1] + 1)
. For triangles there will be will be 2*n[0]*n[1]
cells. For quadrilaterals the number of cells will be n[0]*n[1]
.
[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] | ghost_mode | Mesh ghosting mode |
[in] | diagonal | Direction of diagonals |
std::tuple< Mesh, std::vector< std::int32_t >, std::vector< std::int32_t >, std::vector< std::int32_t > > create_submesh | ( | const Mesh & | mesh, |
int | dim, | ||
const 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 indicies in mesh to include in the new mesh |
Topology create_topology | ( | MPI_Comm | comm, |
const graph::AdjacencyList< std::int64_t > & | cells, | ||
const std::span< const std::int64_t > & | original_cell_index, | ||
const std::span< const int > & | ghost_owners, | ||
const CellType & | cell_type, | ||
GhostMode | ghost_mode | ||
) |
Create a distributed mesh topology.
[in] | comm | MPI communicator across which the topology is distributed |
[in] | cells | 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. |
[in] | original_cell_index | The original global index associated with each cell |
[in] | ghost_owners | The owning rank of each ghost cell (ghost cells are always at the end of the list of cells ) |
[in] | cell_type | The cell shape |
[in] | ghost_mode | Type of cell ghosting: none, shared_facet or shared_vertex |
std::vector< std::int32_t > entities_to_geometry | ( | const Mesh & | mesh, |
int | dim, | ||
const std::span< const std::int32_t > & | entities, | ||
bool | orient | ||
) |
Determine the indices in the geometry data for each vertex of the given mesh entities.
[in] | mesh | The mesh |
[in] | dim | Topological dimension of the entities of interest |
[in] | entities | Entity indices (local) to compute the vertex geometry indices for |
[in] | orient | If true, in 3D, reorients facets to have consistent normal direction |
(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]
. 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.
[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 | The mesh 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.
[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
. std::vector< double > h | ( | const Mesh & | mesh, |
const std::span< const std::int32_t > & | entities, | ||
int | dim | ||
) |
Compute greatest distance between any two vertices of the mesh entities (h
).
[in] | mesh | The 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 & | mesh, |
int | dim, | ||
const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> & | marker | ||
) |
Compute indices of all mesh entities that evaluate to true for the provided geometric marking function. An entity is considered marked if the marker function evaluates true for all of its vertices.
[in] | mesh | The mesh |
[in] | dim | The topological dimension of the entities to be considered |
[in] | marker | The marking function |
std::vector< std::int32_t > locate_entities_boundary | ( | const Mesh & | mesh, |
int | dim, | ||
const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> & | marker | ||
) |
Compute indices of all mesh entities that are attached to an owned boundary facet and evaluate to true for the provided geometric marking function. An entity is considered marked if the marker function evaluates true for all of its vertices.
[in] | mesh | The mesh |
[in] | dim | The topological dimension of the entities to be considered. Must be less than the topological dimension of the mesh. |
[in] | marker | The marking function |
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 |