Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d9/d93/Topology_8h_source.html
DOLFINx 0.7.3
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:
47 Topology(MPI_Comm comm, std::vector<CellType> type);
48
50 Topology(const Topology& topology) = default;
51
53 Topology(Topology&& topology) = default;
54
56 ~Topology() = default;
57
59 Topology& operator=(const Topology& topology) = delete;
60
62 Topology& operator=(Topology&& topology) = default;
63
65 int dim() const noexcept;
66
72 const std::vector<std::int32_t>& offsets);
73
85 const std::vector<std::int32_t>& entity_group_offsets(int dim) const;
86
91 void set_index_map(int dim, std::shared_ptr<const common::IndexMap> map);
92
99 std::shared_ptr<const common::IndexMap> index_map(int dim) const;
100
109 std::shared_ptr<const graph::AdjacencyList<std::int32_t>>
110 connectivity(int d0, int d1) const;
111
114 void set_connectivity(std::shared_ptr<graph::AdjacencyList<std::int32_t>> c,
115 int d0, int d1);
116
118 const std::vector<std::uint32_t>& get_cell_permutation_info() const;
119
132 const std::vector<std::uint8_t>& get_facet_permutations() const;
133
136 std::vector<CellType> cell_types() const noexcept;
137
142 std::int32_t create_entities(int dim);
143
148 void create_connectivity(int d0, int d1);
149
152
155 const std::vector<std::int32_t>& interprocess_facets() const;
156
158 std::vector<std::int64_t> original_cell_index;
159
162 MPI_Comm comm() const;
163
164private:
165 // MPI communicator
166 dolfinx::MPI::Comm _comm;
167
168 // Cell types
169 std::vector<CellType> _cell_types;
170
171 // Entity group offsets
172 std::array<std::vector<std::int32_t>, 4> _entity_group_offsets;
173
174 // Parallel layout of entities for each dimension
175 std::array<std::shared_ptr<const common::IndexMap>, 4> _index_map;
176
177 // AdjacencyList for pairs [d0][d1] == d0 -> d1 connectivity
178 std::vector<std::vector<std::shared_ptr<graph::AdjacencyList<std::int32_t>>>>
179 _connectivity;
180
181 // The facet permutations (local facet, cell))
182 // [cell0_0, cell0_1, ,cell0_2, cell1_0, cell1_1, ,cell1_2, ...,
183 // celln_0, celln_1, ,celln_2,]
184 std::vector<std::uint8_t> _facet_permutations;
185
186 // Cell permutation info. See the documentation for
187 // get_cell_permutation_info for documentation of how this is encoded.
188 std::vector<std::uint32_t> _cell_permutations;
189
190 // List of facets that are on the inter-process boundary
191 std::vector<std::int32_t> _interprocess_facets;
192};
193
214 const graph::AdjacencyList<std::int64_t>& cells,
215 std::span<const std::int64_t> original_cell_index,
216 std::span<const int> ghost_owners,
217 const std::vector<CellType>& cell_type,
218 const std::vector<std::int32_t>& cell_group_offsets,
219 std::span<const std::int64_t> boundary_vertices);
220
231std::tuple<Topology, std::vector<int32_t>, std::vector<int32_t>>
232create_subtopology(const Topology& topology, int dim,
233 std::span<const std::int32_t> entities);
234
245std::vector<std::int32_t>
246entities_to_index(const Topology& topology, int dim,
247 const graph::AdjacencyList<std::int32_t>& entities);
248} // 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:748
void create_connectivity(int d0, int d1)
Create connectivity between given pair of dimensions, d0 -> d1.
Definition Topology.cpp:787
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.
Definition Topology.cpp:836
std::vector< CellType > cell_types() const noexcept
Cell type.
Definition Topology.cpp:885
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:843
void create_entity_permutations()
Compute entity permutations and reflections.
Definition Topology.cpp:815
std::int32_t create_entities(int dim)
Create entities of given topological dimension.
Definition Topology.cpp:754
const std::vector< std::uint32_t > & get_cell_permutation_info() const
Returns the permutation information.
Definition Topology.cpp:851
void set_entity_group_offsets(int dim, const std::vector< std::int32_t > &offsets)
Set the offsets for each group of entities of a particular dimension. See entity_group_offsets.
Definition Topology.cpp:730
const std::vector< std::int32_t > & interprocess_facets() const
List of inter-process facets, if facet topology has been computed.
Definition Topology.cpp:880
const std::vector< std::uint8_t > & get_facet_permutations() const
Get the permutation number to apply to a facet.
Definition Topology.cpp:866
Topology(const Topology &topology)=default
Copy constructor.
~Topology()=default
Destructor.
int dim() const noexcept
Return the topological dimension of the mesh.
Definition Topology.cpp:728
std::vector< std::int64_t > original_cell_index
Original cell index.
Definition Topology.h:158
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:741
Topology & operator=(const Topology &topology)=delete
Assignment.
MPI_Comm comm() const
Mesh MPI communicator.
Definition Topology.cpp:890
const std::vector< std::int32_t > & entity_group_offsets(int dim) const
Get the offsets for each group of entities of a particular dimension.
Definition Topology.cpp:736
Topology & operator=(Topology &&topology)=default
Assignment.
Miscellaneous classes, functions and types.
Definition dolfinx_common.h:8
Graph data structures and algorithms.
Definition dofmapbuilder.h:25
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32
Topology create_topology(MPI_Comm comm, const graph::AdjacencyList< std::int64_t > &cells, std::span< const std::int64_t > original_cell_index, std::span< const int > ghost_owners, const std::vector< CellType > &cell_type, const std::vector< std::int32_t > &cell_group_offsets, std::span< const std::int64_t > boundary_vertices)
Create a distributed mesh topology.
Definition Topology.cpp:892
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:1157
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.
Definition Topology.cpp:1271
CellType
Cell type identifier.
Definition cell_types.h:22
Top-level namespace.
Definition defines.h:12