DOLFINx
0.1.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 |
A MeshTags are used to associate mesh entities with values. The entity index (local to process) identifies the entity. MeshTags is a sparse data storage class; it allows tags to be associated with an arbitrary subset of mesh entities. An entity can have only one associated tag. 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< const 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)> |
Enumerations | |
enum | CellType : int { point = 1, interval = 2, triangle = 3, tetrahedron = 4, quadrilateral = -4, pyramid = -5, prism = -6, hexahedron = -8 } |
Cell type identifier. | |
enum | 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) |
Return type of cell for entity of dimension d. | |
CellType | cell_facet_type (CellType type) |
Return facet type of cell. 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. | |
xt::xtensor< int, 2 > | 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 (mesh::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 (mesh::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::Geometry | create_geometry (MPI_Comm comm, const Topology &topology, const fem::CoordinateElement &coordinate_element, const graph::AdjacencyList< std::int64_t > &cells, const xt::xtensor< double, 2 > &x) |
Build Geometry FIXME: document. | |
std::pair< graph::AdjacencyList< std::int64_t >, std::array< std::int32_t, 2 > > | build_dual_graph (const MPI_Comm comm, const graph::AdjacencyList< std::int64_t > &cell_vertices, int tdim) |
Build distributed dual graph (cell-cell connections) from minimal mesh data, and return (graph, ghost_vertices, [num local edges, num non-local edges]) | |
std::pair< graph::AdjacencyList< std::int32_t >, xt::xtensor< std::int64_t, 2 > > | build_local_dual_graph (const graph::AdjacencyList< std::int64_t > &cell_vertices, int tdim) |
Compute local part of the dual graph, and return (local_graph, facet_cell_map, number of local edges in the graph (undirected) | |
Mesh | create_mesh (MPI_Comm comm, const graph::AdjacencyList< std::int64_t > &cells, const fem::CoordinateElement &element, const xt::xtensor< double, 2 > &x, GhostMode ghost_mode) |
Create a mesh using the default partitioner. This function takes mesh input data that is distributed across processes and creates a Mesh , with the cell distribution determined by the default cell partitioner. The default partitioner is based a graph partitioning. More... | |
Mesh | create_mesh (MPI_Comm comm, const graph::AdjacencyList< std::int64_t > &cells, const fem::CoordinateElement &element, const xt::xtensor< double, 2 > &x, GhostMode ghost_mode, const CellPartitionFunction &cell_partitioner) |
Create a mesh using a provided mesh partitioning function. | |
template<typename T > | |
mesh::MeshTags< T > | create_meshtags (const std::shared_ptr< const mesh::Mesh > &mesh, const int dim, const graph::AdjacencyList< std::int32_t > &entities, const xtl::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... | |
std::vector< bool > | compute_boundary_facets (const Topology &topology) |
Compute marker for owned facets that are on the exterior of the domain, i.e. are connected to only one cell. The function does not require parallel communication. More... | |
Topology | create_topology (MPI_Comm comm, const graph::AdjacencyList< std::int64_t > &cells, const std::vector< std::int64_t > &original_cell_index, const std::vector< int > &ghost_owners, const CellType &cell_type, mesh::GhostMode ghost_mode) |
Create distributed topology. 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 xtl::span< const std::int32_t > &entities, int dim) |
Compute greatest distance between any two vertices. | |
xt::xtensor< double, 2 > | cell_normals (const Mesh &mesh, int dim, const xtl::span< const std::int32_t > &entities) |
Compute normal to given cell (viewed as embedded in 3D) | |
xt::xtensor< double, 2 > | midpoints (const mesh::Mesh &mesh, int dim, const xtl::span< const std::int32_t > &entities) |
Compute midpoints or mesh entities of a given dimension. | |
std::vector< std::int32_t > | locate_entities (const mesh::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 &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... | |
xt::xtensor< std::int32_t, 2 > | entities_to_geometry (const mesh::Mesh &mesh, int dim, const xtl::span< const std::int32_t > &entity_list, bool orient) |
Compute the indices the geometry data for the vertices of the given mesh entities. More... | |
std::vector< std::int32_t > | exterior_facet_indices (const Mesh &mesh) |
Compute the indices (local) of all exterior facets. An exterior facet (co-dimension 1) is one that is connected globally to only one cell of co-dimension 0). More... | |
graph::AdjacencyList< std::int32_t > | partition_cells_graph (MPI_Comm comm, int n, int tdim, const graph::AdjacencyList< std::int64_t > &cells, mesh::GhostMode ghost_mode) |
Compute destination rank for mesh cells in this rank by applying the default graph partitioner to the dual graph of the mesh. More... | |
graph::AdjacencyList< std::int32_t > | partition_cells_graph (MPI_Comm comm, int n, int tdim, const graph::AdjacencyList< std::int64_t > &cells, mesh::GhostMode ghost_mode, const graph::partition_fn &partfn) |
Compute destination rank for mesh cells on this rank by applying the a provided graph partitioner to the dual graph of the mesh. | |
Mesh data structures and algorithms on meshes.
Representations of meshes and support for operations on meshes.
using dolfinx::mesh::CellPartitionFunction = typedef std::function<const 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.
mesh::CellType dolfinx::mesh::cell_facet_type | ( | mesh::CellType | type | ) |
Return facet type of cell.
[in] | type | The cell type |
int dolfinx::mesh::cell_num_entities | ( | mesh::CellType | type, |
int | dim | ||
) |
Number of entities of dimension dim.
[in] | dim | Entity dimension |
[in] | type | Cell type |
std::vector< bool > dolfinx::mesh::compute_boundary_facets | ( | const Topology & | topology | ) |
Compute marker for owned facets that are on the exterior of the domain, i.e. are connected to only one cell. The function does not require parallel communication.
[in] | topology | The topology |
std::array< std::shared_ptr< graph::AdjacencyList< std::int32_t > >, 2 > dolfinx::mesh::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 > > dolfinx::mesh::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 > > dolfinx::mesh::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.
Mesh dolfinx::mesh::create_mesh | ( | MPI_Comm | comm, |
const graph::AdjacencyList< std::int64_t > & | cells, | ||
const fem::CoordinateElement & | element, | ||
const xt::xtensor< double, 2 > & | x, | ||
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
, with the 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] | ghost_mode | The requested type of cell ghosting/overlap |
mesh::MeshTags<T> dolfinx::mesh::create_meshtags | ( | const std::shared_ptr< const mesh::Mesh > & | mesh, |
const int | dim, | ||
const graph::AdjacencyList< std::int32_t > & | entities, | ||
const xtl::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. |
Topology dolfinx::mesh::create_topology | ( | MPI_Comm | comm, |
const graph::AdjacencyList< std::int64_t > & | cells, | ||
const std::vector< std::int64_t > & | original_cell_index, | ||
const std::vector< int > & | ghost_owners, | ||
const CellType & | cell_type, | ||
mesh::GhostMode | ghost_mode | ||
) |
Create distributed topology.
[in] | comm | MPI communicator across which the topology is distributed |
[in] | cells | The cell topology (list of cell vertices) 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 which 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 ownership of the ghost cells (ghost cells are always at the end of the list of cells, above) |
[in] | cell_type | The cell shape |
[in] | ghost_mode | How to partition the cell overlap: none, shared_facet or shared_vertex. |
xt::xtensor< std::int32_t, 2 > dolfinx::mesh::entities_to_geometry | ( | const mesh::Mesh & | mesh, |
int | dim, | ||
const xtl::span< const std::int32_t > & | entity_list, | ||
bool | orient | ||
) |
Compute the indices the geometry data for the vertices of the given mesh entities.
[in] | mesh | Mesh |
[in] | dim | Topological dimension of the entities of interest |
[in] | entity_list | List of entity indices (local) |
[in] | orient | If true, in 3D, reorients facets to have consistent normal direction |
std::vector< std::int32_t > dolfinx::mesh::exterior_facet_indices | ( | const Mesh & | mesh | ) |
Compute the indices (local) of all exterior facets. An exterior facet (co-dimension 1) is one that is connected globally to only one cell of co-dimension 0).
[in] | mesh | Mesh |
graph::AdjacencyList< std::int64_t > dolfinx::mesh::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
. bool dolfinx::mesh::is_simplex | ( | mesh::CellType | type | ) |
Check if cell is a simplex.
[in] | type | Cell type |
std::vector< std::int32_t > dolfinx::mesh::locate_entities | ( | const mesh::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 > dolfinx::mesh::locate_entities_boundary | ( | const mesh::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 dolfinx::mesh::num_cell_vertices | ( | mesh::CellType | type | ) |
Number vertices for a cell type.
[in] | type | Cell type |
graph::AdjacencyList< std::int32_t > dolfinx::mesh::partition_cells_graph | ( | MPI_Comm | comm, |
int | n, | ||
int | tdim, | ||
const graph::AdjacencyList< std::int64_t > & | cells, | ||
mesh::GhostMode | ghost_mode | ||
) |
Compute destination rank for mesh cells in this rank by applying the default graph partitioner to the dual graph of the mesh.
[in] | comm | MPI Communicator |
[in] | n | 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 |
std::string dolfinx::mesh::to_string | ( | mesh::CellType | type | ) |
Get the cell string type for a cell type.
[in] | type | The cell type |
mesh::CellType dolfinx::mesh::to_type | ( | const std::string & | cell | ) |
Get the cell type from a cell string.
[in] | cell | Cell shape string |