OGS
transformMeshToNodePartitionedMesh.cpp
Go to the documentation of this file.
1
11
12#include <mpi.h>
13
14#include <numeric>
15#include <range/v3/algorithm/count_if.hpp>
16#include <range/v3/numeric.hpp>
17#include <range/v3/range/conversion.hpp>
18#include <range/v3/view/enumerate.hpp>
19#include <range/v3/view/indirect.hpp>
20#include <range/v3/view/map.hpp>
21#include <unordered_map>
22#include <vector>
23
24#include "BaseLib/Logging.h"
25#include "BaseLib/MPI.h"
27#include "MeshLib/Mesh.h"
28#include "MeshLib/Node.h"
31
32namespace
33{
35 MeshLib::NodePartitionedMesh const& bulk_mesh,
36 std::vector<std::size_t> const& local_bulk_node_ids_for_subdomain,
37 std::size_t const subdomain_node_id)
38{
39 return (!bulk_mesh.isGhostNode(
40 local_bulk_node_ids_for_subdomain[subdomain_node_id]));
41};
42} // namespace
43
44namespace MeshLib
45{
46std::pair<std::vector<Node*>, std::vector<Element*>> copyNodesAndElements(
47 std::vector<Element*> const& input_elements)
48{
49 auto elements = cloneElements(input_elements);
50
51 // original node ids to newly created nodes.
52 std::unordered_map<std::size_t, Node*> id_node_hash_map;
53 id_node_hash_map.reserve(
54 elements.size()); // There will be at least one node per element.
55
56 for (auto& e : elements)
57 {
58 // For each node find a cloned node in map or create if there is none.
59 unsigned const n_nodes = e->getNumberOfNodes();
60 for (unsigned i = 0; i < n_nodes; ++i)
61 {
62 Node const* n = e->getNode(i);
63 auto const it = id_node_hash_map.find(n->getID());
64 if (it == id_node_hash_map.end())
65 {
66 auto new_node_in_map = id_node_hash_map[n->getID()] =
67 new Node(*n);
68 e->setNode(i, new_node_in_map);
69 }
70 else
71 {
72 e->setNode(i, it->second);
73 }
74 }
75 }
76
77 std::map<std::size_t, Node*> nodes_map;
78 for (auto const& n : id_node_hash_map)
79 {
80 nodes_map[n.first] = n.second;
81 }
82
83 // Copy the unique nodes pointers.
84 auto element_nodes =
85 nodes_map | ranges::views::values | ranges::to<std::vector>;
86
87 return std::make_pair(element_nodes, elements);
88}
89
90unsigned long computeNumberOfRegularNodes(NodePartitionedMesh const* bulk_mesh,
91 Mesh const* subdomain_mesh)
92{
93 auto const& subdomain_nodes = subdomain_mesh->getNodes();
94 auto const& local_bulk_node_ids_for_subdomain =
95 *bulkNodeIDs(*subdomain_mesh);
96 unsigned long const number_of_regular_nodes = ranges::count_if(
97 subdomain_nodes | MeshLib::views::ids,
98 [&](std::size_t const id) {
99 return isRegularNode(*bulk_mesh, local_bulk_node_ids_for_subdomain,
100 id);
101 });
102
103 DBUG("[{}] number of regular nodes: {}", subdomain_mesh->getName(),
104 number_of_regular_nodes);
105 return number_of_regular_nodes;
106}
107
108std::vector<std::size_t>
110 NodePartitionedMesh const* bulk_mesh, Mesh const* subdomain_mesh)
111{
112 // create/fill global nodes information of the sub-domain partition
113 auto const& subdomain_nodes = subdomain_mesh->getNodes();
114 auto const& local_bulk_node_ids_for_subdomain =
115 *bulkNodeIDs(*subdomain_mesh);
116
117 // global node ids for the sub-domain:
118 // | regular base node ids | ghost base node ids | regular higher
119 // order node ids | ghost higher order node ids |
120 // | partition offset + local ids | regular node ids from other
121 // partitions |
122
123 // In order to compute the global node ids for the subdomain mesh nodes the
124 // sum of the number of regular nodes of the previous partitions is required
125 // first.
126 // Count the local regular base nodes to compute offset for other subsequent
127 // partitions
128 std::size_t const number_of_regular_nodes =
129 computeNumberOfRegularNodes(bulk_mesh, subdomain_mesh);
130
132
133 // send own number of regular nodes to all others
134 std::vector<std::size_t> const gathered_number_of_regular_nodes =
135 BaseLib::MPI::allgather(number_of_regular_nodes, mpi);
136
137 // compute the 'offset' in the global_node_ids
138 std::vector<std::size_t> const numbers_of_regular_nodes_at_rank =
139 BaseLib::sizesToOffsets(gathered_number_of_regular_nodes);
140
141 // add the offset to the partitioned-owned subdomain
142 std::vector<std::size_t> subdomain_global_node_ids;
143 subdomain_global_node_ids.reserve(subdomain_nodes.size());
144 auto const partition_offset = numbers_of_regular_nodes_at_rank[mpi.rank];
145 DBUG("[{}] partition offset: {}", subdomain_mesh->getName(),
146 partition_offset);
147 // set the global id for the regular base nodes
148 for (auto const id : subdomain_nodes | MeshLib::views::ids)
149 {
150 if (isRegularNode(*bulk_mesh, local_bulk_node_ids_for_subdomain, id))
151 {
152 subdomain_global_node_ids.emplace_back(partition_offset + id);
153 }
154 }
155 return subdomain_global_node_ids;
156}
157
159 NodePartitionedMesh const* bulk_mesh,
160 Mesh const* subdomain_mesh,
161 std::vector<std::size_t> const& global_regular_base_node_ids)
162{
164
165 // count regular nodes that is the offset in the mapping
166 auto const local_bulk_node_ids_for_subdomain =
167 *bulkNodeIDs(*subdomain_mesh);
168
169 // compute mapping subdomain ghost node id <-> bulk ghost node id
170 std::vector<std::size_t> subdomain_node_id_to_bulk_node_id;
171 for (auto const id : subdomain_mesh->getNodes() | MeshLib::views::ids)
172 {
173 if (!isRegularNode(*bulk_mesh, local_bulk_node_ids_for_subdomain, id))
174 {
175 auto const bulk_node_id = local_bulk_node_ids_for_subdomain[id];
176 // this puts the global bulk node id at pos:
177 // subdomain_node->getID() - number_of_regular_nodes
178 subdomain_node_id_to_bulk_node_id.push_back(
179 bulk_mesh->getGlobalNodeID(bulk_node_id));
180 }
181 }
182
183 // Send ids of bulk ghost nodes belonging to the subdomain mesh to all
184 // other ranks and at the same time receive from all other ranks.
185 std::vector<std::size_t> ghost_node_ids_of_all_ranks;
186 std::vector<int> const offsets =
187 BaseLib::MPI::allgatherv(std::span{subdomain_node_id_to_bulk_node_id},
188 ghost_node_ids_of_all_ranks, mpi);
189
190 std::size_t const global_number_of_subdomain_node_id_to_bulk_node_id =
191 offsets.back();
192 DBUG("[{}] global_number_of_subdomain_node_id_to_bulk_node_id: '{}' ",
193 subdomain_mesh->getName(),
194 global_number_of_subdomain_node_id_to_bulk_node_id);
195
196 // construct a map for fast search of local bulk node ids
197 std::map<std::size_t, std::size_t> global_to_local_bulk_node_ids;
198 for (auto const id : bulk_mesh->getNodes() | MeshLib::views::ids)
199 {
200 if (!bulk_mesh->isGhostNode(id))
201 {
202 global_to_local_bulk_node_ids[bulk_mesh->getGlobalNodeID(id)] = id;
203 }
204 }
205
206 // construct a map for fast search of local sub-domain node ids
207 std::map<std::size_t, std::size_t> local_subdomain_node_ids;
208 for (auto const id : subdomain_mesh->getNodes() | MeshLib::views::ids)
209 {
210 local_subdomain_node_ids[local_bulk_node_ids_for_subdomain[id]] = id;
211 }
212
213 std::vector<std::size_t> subdomain_node_ids_of_all_ranks(
214 global_number_of_subdomain_node_id_to_bulk_node_id,
215 std::numeric_limits<std::size_t>::max());
216 // search in all ranks within the bulk ids for the corresponding id
217 for (int rank = 0; rank < mpi.size; ++rank)
218 {
219 if (rank == mpi.rank)
220 {
221 continue;
222 }
223 for (int i = offsets[rank]; i < offsets[rank + 1]; ++i)
224 {
225 auto it = global_to_local_bulk_node_ids.find(
226 ghost_node_ids_of_all_ranks[i]);
227 if (it == global_to_local_bulk_node_ids.end())
228 {
229 continue;
230 }
231 auto const local_bulk_node_id = it->second;
232 subdomain_node_ids_of_all_ranks[i] = global_regular_base_node_ids
233 [local_subdomain_node_ids.find(local_bulk_node_id)->second];
234 DBUG(
235 "[{}] found global subdomain node id: '{}' for global bulk "
236 "node id {} ",
237 subdomain_mesh->getName(),
238 subdomain_node_ids_of_all_ranks[i],
239 ghost_node_ids_of_all_ranks[i]);
240 }
241 }
242
243 // find maximum over all ranks
244 BaseLib::MPI::allreduceInplace(subdomain_node_ids_of_all_ranks, MPI_MAX,
245 mpi);
246
247 std::vector<std::size_t> global_ids_for_subdomain_ghost_nodes(
248 subdomain_node_ids_of_all_ranks.begin() + offsets[mpi.rank],
249 subdomain_node_ids_of_all_ranks.begin() + offsets[mpi.rank + 1]);
250 return global_ids_for_subdomain_ghost_nodes;
251}
252
254 Mesh const* subdomain_mesh)
255{
256 auto const number_of_regular_base_nodes =
257 subdomain_mesh->computeNumberOfBaseNodes();
258
260
261 std::vector<std::size_t> const gathered_number_of_regular_base_nodes =
262 BaseLib::MPI::allgather(number_of_regular_base_nodes, mpi);
263
264 return BaseLib::sizesToOffsets(gathered_number_of_regular_base_nodes);
265}
266
267// similar to the above only with regular higher order nodes
269 Mesh const* subdomain_mesh)
270{
272
273 auto const number_of_regular_base_nodes =
274 subdomain_mesh->computeNumberOfBaseNodes();
275
276 auto const number_of_regular_higher_order_nodes =
277 subdomain_mesh->getNumberOfNodes() - number_of_regular_base_nodes;
278 std::vector<std::size_t> gathered_number_of_regular_higher_order_nodes =
279 BaseLib::MPI::allgather(number_of_regular_higher_order_nodes, mpi);
280
282 gathered_number_of_regular_higher_order_nodes);
283}
284
286 NodePartitionedMesh const* bulk_mesh, Mesh const* subdomain_mesh)
287{
288 auto global_regular_base_node_ids =
290 subdomain_mesh);
291 auto const local_ids_for_subdomain_ghost_nodes =
293 bulk_mesh, subdomain_mesh, global_regular_base_node_ids);
294
295 // At the moment higher order bulk meshes aren't handled correctly.
296 // If it should become necessary in the future, the necessary things must be
297 // implemented at this point.
298 /*
299 auto const regular_higher_order_node_global_node_ids =
300 computeRegularHigherOrderNodeGlobalNodeIDsofSubDomainPartition(
301 bulk_mesh, subdomain_mesh);
302 auto const ghost_higher_order_node_global_node_ids =
303 computeGhostHigherOrderNodeGlobalNodeIDsofSubDomainPartition(
304 bulk_mesh, subdomain_mesh);
305 */
306
307 std::vector<std::size_t> global_node_ids(
308 std::move(global_regular_base_node_ids));
309 global_node_ids.insert(global_node_ids.end(),
310 local_ids_for_subdomain_ghost_nodes.begin(),
311 local_ids_for_subdomain_ghost_nodes.end());
312
313 return global_node_ids;
314}
315
316unsigned long getNumberOfGlobalNodes(Mesh const* subdomain_mesh)
317{
318 // sum all nodes over all partitions in number_of_global_nodes
319 unsigned long const number_of_global_nodes = BaseLib::MPI::allreduce(
320 subdomain_mesh->getNodes().size(), MPI_SUM, BaseLib::MPI::Mpi{});
321 DBUG("[{}] number_of_global_nodes: {}'", subdomain_mesh->getName(),
322 number_of_global_nodes);
323 return number_of_global_nodes;
324}
325
326std::unique_ptr<NodePartitionedMesh> transformMeshToNodePartitionedMesh(
327 NodePartitionedMesh const* const bulk_mesh,
328 Mesh const* const subdomain_mesh)
329{
330 DBUG("Creating NodePartitionedMesh from '{}'", subdomain_mesh->getName());
331
332 auto const subdomain_global_node_ids =
333 computeGlobalNodeIDsOfSubDomainPartition(bulk_mesh, subdomain_mesh);
334
335 // according to comment in struct PartitionedMeshInfo this value
336 // is unused
337 unsigned long number_of_global_base_nodes = 0;
338
339 unsigned long const number_of_global_nodes =
340 getNumberOfGlobalNodes(subdomain_mesh);
341 auto numbers_of_regular_base_nodes_at_rank =
343 auto numbers_of_regular_higher_order_nodes_at_rank =
345
346 auto const number_of_regular_nodes =
347 computeNumberOfRegularNodes(bulk_mesh, subdomain_mesh);
348 DBUG("[{}] number_of_regular_nodes: {}", subdomain_mesh->getName(),
349 number_of_regular_nodes);
350 auto const [nodes, elements] =
351 copyNodesAndElements(subdomain_mesh->getElements());
352
353 return std::make_unique<NodePartitionedMesh>(
354 subdomain_mesh->getName(), nodes, subdomain_global_node_ids, elements,
355 subdomain_mesh->getProperties(), number_of_global_base_nodes,
356 number_of_global_nodes, number_of_regular_nodes,
357 std::move(numbers_of_regular_base_nodes_at_rank),
358 std::move(numbers_of_regular_higher_order_nodes_at_rank));
359}
360
361} // namespace MeshLib
Definition of Duplicate functions.
Definition of the Element class.
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Definition of the Mesh class.
Definition of mesh class for partitioned mesh (by node) for parallel computing within the framework o...
Definition of the Node class.
std::size_t getID() const
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
std::size_t computeNumberOfBaseNodes() const
Get the number of base nodes.
Definition Mesh.cpp:238
Properties & getProperties()
Definition Mesh.h:134
const std::string getName() const
Get name of the mesh.
Definition Mesh.h:103
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:100
bool isGhostNode(const std::size_t node_id) const
Check whether a node with ID of node_id is a ghost node.
std::size_t getGlobalNodeID(const std::size_t node_id) const
Get the global node ID of a node with its local ID.
static void allreduceInplace(std::vector< T > &vector, MPI_Op const &mpi_op, Mpi const &mpi)
Definition MPI.h:149
static std::vector< int > allgatherv(std::span< T > const send_buffer, std::vector< std::remove_const_t< T > > &receive_buffer, Mpi const &mpi)
Definition MPI.h:164
static T allreduce(T const &value, MPI_Op const &mpi_op, Mpi const &mpi)
Definition MPI.h:128
static std::vector< T > allgather(T const &value, Mpi const &mpi)
Definition MPI.h:100
std::vector< ranges::range_value_t< R > > sizesToOffsets(R const &sizes)
Definition Algorithm.h:283
constexpr ranges::views::view_closure ids
For an element of a range view return its id.
Definition Mesh.h:225
std::pair< std::vector< Node * >, std::vector< Element * > > copyNodesAndElements(std::vector< Element * > const &input_elements)
unsigned long getNumberOfGlobalNodes(Mesh const *subdomain_mesh)
std::vector< std::size_t > computeNumberOfRegularHigherOrderNodesAtRank(Mesh const *subdomain_mesh)
std::vector< std::size_t > computeGlobalNodeIDsOfSubDomainPartition(NodePartitionedMesh const *bulk_mesh, Mesh const *subdomain_mesh)
std::vector< std::size_t > computeNumberOfRegularBaseNodesAtRank(Mesh const *subdomain_mesh)
std::vector< std::size_t > computeRegularBaseNodeGlobalNodeIDsOfSubDomainPartition(NodePartitionedMesh const *bulk_mesh, Mesh const *subdomain_mesh)
unsigned long computeNumberOfRegularNodes(NodePartitionedMesh const *bulk_mesh, Mesh const *subdomain_mesh)
std::vector< Element * > cloneElements(std::vector< Element * > const &elements)
Clones a vector of elements using the Element::clone() function.
std::vector< std::size_t > computeGhostBaseNodeGlobalNodeIDsOfSubDomainPartition(NodePartitionedMesh const *bulk_mesh, Mesh const *subdomain_mesh, std::vector< std::size_t > const &global_regular_base_node_ids)
PropertyVector< std::size_t > const * bulkNodeIDs(Mesh const &mesh)
Definition Mesh.cpp:292
std::unique_ptr< NodePartitionedMesh > transformMeshToNodePartitionedMesh(NodePartitionedMesh const *const bulk_mesh, Mesh const *const subdomain_mesh)
bool isRegularNode(MeshLib::NodePartitionedMesh const &bulk_mesh, std::vector< std::size_t > const &local_bulk_node_ids_for_subdomain, std::size_t const subdomain_node_id)