DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Topology.h
1// Copyright (C) 2006-2022 Anders Logg and Garth N. Wells
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include <array>
10#include <cstdint>
11#include <dolfinx/common/MPI.h>
12#include <memory>
13#include <span>
14#include <tuple>
15#include <vector>
16
17namespace dolfinx::common
18{
19class IndexMap;
20}
21
22namespace dolfinx::graph
23{
24template <typename T>
25class AdjacencyList;
26}
27
28namespace dolfinx::mesh
29{
30enum class CellType;
31
44{
45public:
50
55 Topology(MPI_Comm comm, const std::vector<CellType>& cell_type);
56
58 Topology(const Topology& topology) = default;
59
61 Topology(Topology&& topology) = default;
62
64 ~Topology() = default;
65
67 Topology& operator=(const Topology& topology) = delete;
68
70 Topology& operator=(Topology&& topology) = default;
71
73 int dim() const noexcept;
74
79 void set_index_map(int dim, std::shared_ptr<const common::IndexMap> map);
80
89 void set_index_map(std::int8_t dim, std::int8_t i,
90 std::shared_ptr<const common::IndexMap> map);
91
98 std::shared_ptr<const common::IndexMap> index_map(int dim) const;
99
103 std::vector<std::shared_ptr<const common::IndexMap>>
104 index_maps(std::int8_t dim) const;
105
114 std::shared_ptr<const graph::AdjacencyList<std::int32_t>>
115 connectivity(int d0, int d1) const;
116
128 std::shared_ptr<const graph::AdjacencyList<std::int32_t>>
129 connectivity(std::pair<std::int8_t, std::int8_t> d0,
130 std::pair<std::int8_t, std::int8_t> d1) const;
131
134 void set_connectivity(std::shared_ptr<graph::AdjacencyList<std::int32_t>> c,
135 int d0, int d1);
136
147 void set_connectivity(std::shared_ptr<graph::AdjacencyList<std::int32_t>> c,
148 std::pair<std::int8_t, std::int8_t> d0,
149 std::pair<std::int8_t, std::int8_t> d1);
150
152 const std::vector<std::uint32_t>& get_cell_permutation_info() const;
153
168 const std::vector<std::uint8_t>& get_facet_permutations() const;
169
172 CellType cell_type() const;
173
177 std::vector<CellType> entity_types(std::int8_t dim) const;
178
183 std::int32_t create_entities(int dim);
184
189 void create_connectivity(int d0, int d1);
190
193
202 const std::vector<std::int32_t>& interprocess_facets() const;
203
208 const std::vector<std::int32_t>& interprocess_facets(std::int8_t index) const;
209
211 std::vector<std::vector<std::int64_t>> original_cell_index;
212
215 MPI_Comm comm() const;
216
217private:
218 // MPI communicator
219 dolfinx::MPI::Comm _comm;
220
221 // Cell types for entites in Topology, as follows:
222 // [CellType::point, edge_types..., facet_types..., cell_types...]
223 // Only one type is expected for vertices, (and usually edges), but facets
224 // and cells can be a list of multiple types, e.g. [quadrilateral, triangle]
225 // for facets.
226 // Offsets are position in the list for each entity dimension, in
227 // AdjacencyList style.
228 std::vector<CellType> _entity_types;
229 std::vector<std::int8_t> _entity_type_offsets;
230
231 // Parallel layout of entities for each dimension and cell type
232 // flattened in the same layout as _entity_types above.
233 std::vector<std::shared_ptr<const common::IndexMap>> _index_map;
234
235 // Connectivity between entity dimensions and cell types, arranged as
236 // a 2D array. The indexing follows the order in _entity_types, i.e.
237 // increasing in topological dimension. There may be multiple types in each
238 // dimension, e.g. triangle and quadrilateral facets.
239 // Connectivity between different entity types of same dimension will always
240 // be nullptr.
241 std::vector<std::vector<std::shared_ptr<graph::AdjacencyList<std::int32_t>>>>
242 _connectivity;
243
244 // The facet permutations (local facet, cell))
245 // [cell0_0, cell0_1, ,cell0_2, cell1_0, cell1_1, ,cell1_2, ...,
246 // celln_0, celln_1, ,celln_2,]
247 std::vector<std::uint8_t> _facet_permutations;
248
249 // Cell permutation info. See the documentation for
250 // get_cell_permutation_info for documentation of how this is encoded.
251 std::vector<std::uint32_t> _cell_permutations;
252
253 // List of facets that are on the inter-process boundary for each facet type
254 std::vector<std::vector<std::int32_t>> _interprocess_facets;
255};
256
280Topology create_topology(MPI_Comm comm, std::span<const std::int64_t> cells,
281 std::span<const std::int64_t> original_cell_index,
282 std::span<const int> ghost_owners, CellType cell_type,
283 std::span<const std::int64_t> boundary_vertices);
284
297create_topology(MPI_Comm comm, const std::vector<CellType>& cell_type,
298 std::vector<std::span<const std::int64_t>> cells,
299 std::vector<std::span<const std::int64_t>> original_cell_index,
300 std::vector<std::span<const int>> ghost_owners,
301 std::span<const std::int64_t> boundary_vertices);
302
314std::tuple<Topology, std::vector<int32_t>, std::vector<int32_t>>
315create_subtopology(const Topology& topology, int dim,
316 std::span<const std::int32_t> entities);
317
328std::vector<std::int32_t>
329entities_to_index(const Topology& topology, int dim,
330 std::span<const std::int32_t> entities);
331} // namespace dolfinx::mesh
Topology stores the topology of a mesh, consisting of mesh entities and connectivity (incidence relat...
Definition Topology.h:44
std::shared_ptr< const common::IndexMap > index_map(int dim) const
Get the IndexMap that described the parallel distribution of the mesh entities.
Definition Topology.cpp:815
void create_connectivity(int d0, int d1)
Create connectivity between given pair of dimensions, d0 / -> d1.
Definition Topology.cpp:870
std::shared_ptr< const graph::AdjacencyList< std::int32_t > > connectivity(int d0, int d1) const
Return connectivity from entities of dimension d0 to entities of dimension d1. Assumes only one entit...
Definition Topology.cpp:931
void set_connectivity(std::shared_ptr< graph::AdjacencyList< std::int32_t > > c, int d0, int d1)
Set connectivity for given pair of topological dimensions.
Definition Topology.cpp:955
std::vector< std::vector< std::int64_t > > original_cell_index
Original cell index for each cell type.
Definition Topology.h:211
void create_entity_permutations()
Compute entity permutations and reflections.
Definition Topology.cpp:910
std::int32_t create_entities(int dim)
Create entities of given topological dimension.
Definition Topology.cpp:831
const std::vector< std::uint32_t > & get_cell_permutation_info() const
Returns the permutation information.
Definition Topology.cpp:981
Topology(MPI_Comm comm, CellType cell_type)
Empty Topology constructor.
Definition Topology.cpp:716
const std::vector< std::int32_t > & interprocess_facets() const
List of inter-process facets.
Definition Topology.cpp:1010
const std::vector< std::uint8_t > & get_facet_permutations() const
Get the numbers that encode the number of permutations to apply to facets.
Definition Topology.cpp:996
Topology(const Topology &topology)=default
Copy constructor.
~Topology()=default
Destructor.
std::vector< CellType > entity_types(std::int8_t dim) const
Get the entity types in the topology for a given dimension.
Definition Topology.cpp:1023
int dim() const noexcept
Return the topological dimension of the mesh.
Definition Topology.cpp:794
Topology(Topology &&topology)=default
Move constructor.
void set_index_map(int dim, std::shared_ptr< const common::IndexMap > map)
Set the IndexMap for dimension dim.
Definition Topology.cpp:796
Topology & operator=(const Topology &topology)=delete
Assignment.
MPI_Comm comm() const
Definition Topology.cpp:1031
std::vector< std::shared_ptr< const common::IndexMap > > index_maps(std::int8_t dim) const
Definition Topology.cpp:822
CellType cell_type() const
Cell type.
Definition Topology.cpp:1021
Topology & operator=(Topology &&topology)=default
Assignment.
Miscellaneous classes, functions and types.
Definition dolfinx_common.h:8
Graph data structures and algorithms.
Definition dofmapbuilder.h:26
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32
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.
Definition Topology.cpp:1319
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.
Definition Topology.cpp:1331
CellType
Cell type identifier.
Definition cell_types.h:22
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.
Definition Topology.cpp:1425
Top-level namespace.
Definition defines.h:12