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"
26#include "MeshLib/Mesh.h"
27#include "MeshLib/Node.h"
30
31namespace
32{
34 MeshLib::NodePartitionedMesh const& bulk_mesh,
35 std::vector<std::size_t> const& local_bulk_node_ids_for_subdomain,
36 std::size_t const subdomain_node_id)
37{
38 return (!bulk_mesh.isGhostNode(
39 local_bulk_node_ids_for_subdomain[subdomain_node_id]));
40};
41} // namespace
42
43namespace MeshLib
44{
45std::pair<int, int> getMPIRankAndSize(MPI_Comm const& mpi_comm)
46{
47 int mpi_comm_size;
48 MPI_Comm_size(mpi_comm, &mpi_comm_size);
49 int mpi_comm_rank;
50 MPI_Comm_rank(mpi_comm, &mpi_comm_rank);
51 return {mpi_comm_rank, mpi_comm_size};
52}
53
54std::pair<std::vector<Node*>, std::vector<Element*>> copyNodesAndElements(
55 std::vector<Element*> const& input_elements)
56{
57 auto elements = cloneElements(input_elements);
58
59 // original node ids to newly created nodes.
60 std::unordered_map<std::size_t, Node*> id_node_hash_map;
61 id_node_hash_map.reserve(
62 elements.size()); // There will be at least one node per element.
63
64 for (auto& e : elements)
65 {
66 // For each node find a cloned node in map or create if there is none.
67 unsigned const n_nodes = e->getNumberOfNodes();
68 for (unsigned i = 0; i < n_nodes; ++i)
69 {
70 Node const* n = e->getNode(i);
71 auto const it = id_node_hash_map.find(n->getID());
72 if (it == id_node_hash_map.end())
73 {
74 auto new_node_in_map = id_node_hash_map[n->getID()] =
75 new Node(*n);
76 e->setNode(i, new_node_in_map);
77 }
78 else
79 {
80 e->setNode(i, it->second);
81 }
82 }
83 }
84
85 std::map<std::size_t, Node*> nodes_map;
86 for (auto const& n : id_node_hash_map)
87 {
88 nodes_map[n.first] = n.second;
89 }
90
91 // Copy the unique nodes pointers.
92 auto element_nodes =
93 nodes_map | ranges::views::values | ranges::to<std::vector>;
94
95 return std::make_pair(element_nodes, elements);
96}
97
98unsigned long computeNumberOfRegularNodes(NodePartitionedMesh const* bulk_mesh,
99 Mesh const* subdomain_mesh)
100{
101 auto const& subdomain_nodes = subdomain_mesh->getNodes();
102 auto const& local_bulk_node_ids_for_subdomain =
103 *bulkNodeIDs(*subdomain_mesh);
104 unsigned long const number_of_regular_nodes = ranges::count_if(
105 subdomain_nodes | MeshLib::views::ids,
106 [&](std::size_t const id) {
107 return isRegularNode(*bulk_mesh, local_bulk_node_ids_for_subdomain,
108 id);
109 });
110
111 DBUG("[{}] number of regular nodes: {}", subdomain_mesh->getName(),
112 number_of_regular_nodes);
113 return number_of_regular_nodes;
114}
115
116std::vector<std::size_t>
118 NodePartitionedMesh const* bulk_mesh, Mesh const* subdomain_mesh)
119{
120 // create/fill global nodes information of the sub-domain partition
121 auto const& subdomain_nodes = subdomain_mesh->getNodes();
122 auto const& local_bulk_node_ids_for_subdomain =
123 *bulkNodeIDs(*subdomain_mesh);
124
125 // global node ids for the sub-domain:
126 // | regular base node ids | ghost base node ids | regular higher
127 // order node ids | ghost higher order node ids |
128 // | partition offset + local ids | regular node ids from other
129 // partitions |
130
131 // In order to compute the global node ids for the subdomain mesh nodes the
132 // sum of the number of regular nodes of the previous partitions is required
133 // first.
134 // Count the local regular base nodes to compute offset for other subsequent
135 // partitions
136 std::size_t const number_of_regular_nodes =
137 computeNumberOfRegularNodes(bulk_mesh, subdomain_mesh);
138
139 // in the following information exchange with other ranks is required
140 MPI_Comm mpi_comm = MPI_COMM_WORLD;
141 auto const [mpi_comm_rank, mpi_comm_size] = getMPIRankAndSize(mpi_comm);
142
143 // send own number of regular nodes to all others
144 std::vector<std::size_t> gathered_number_of_regular_nodes(mpi_comm_size);
145 MPI_Allgather(&number_of_regular_nodes, 1, MPI_UNSIGNED_LONG,
146 gathered_number_of_regular_nodes.data(), 1, MPI_UNSIGNED_LONG,
147 mpi_comm);
148 // compute the 'offset' in the global_node_ids
149 std::vector<std::size_t> numbers_of_regular_nodes_at_rank;
150 numbers_of_regular_nodes_at_rank.push_back(0);
151 std::partial_sum(begin(gathered_number_of_regular_nodes),
152 end(gathered_number_of_regular_nodes),
153 back_inserter(numbers_of_regular_nodes_at_rank));
154
155 // add the offset to the partitioned-owned subdomain
156 std::vector<std::size_t> subdomain_global_node_ids;
157 subdomain_global_node_ids.reserve(subdomain_nodes.size());
158 auto const partition_offset =
159 numbers_of_regular_nodes_at_rank[mpi_comm_rank];
160 DBUG("[{}] partition offset: {}", subdomain_mesh->getName(),
161 partition_offset);
162 // set the global id for the regular base nodes
163 for (auto const id : subdomain_nodes | MeshLib::views::ids)
164 {
165 if (isRegularNode(*bulk_mesh, local_bulk_node_ids_for_subdomain, id))
166 {
167 subdomain_global_node_ids.emplace_back(partition_offset + id);
168 }
169 }
170 return subdomain_global_node_ids;
171}
172
174 NodePartitionedMesh const* bulk_mesh,
175 Mesh const* subdomain_mesh,
176 std::vector<std::size_t> const& global_regular_base_node_ids)
177{
178 // in the following information exchange with other ranks is required
179 MPI_Comm const mpi_comm = MPI_COMM_WORLD;
180 auto const [mpi_comm_rank, mpi_comm_size] = getMPIRankAndSize(mpi_comm);
181
182 // count regular nodes that is the offset in the mapping
183 auto const local_bulk_node_ids_for_subdomain =
184 *bulkNodeIDs(*subdomain_mesh);
185
186 // compute mapping subdomain ghost node id <-> bulk ghost node id
187 std::vector<std::size_t> subdomain_node_id_to_bulk_node_id;
188 for (auto const id : subdomain_mesh->getNodes() | MeshLib::views::ids)
189 {
190 if (!isRegularNode(*bulk_mesh, local_bulk_node_ids_for_subdomain, id))
191 {
192 auto const bulk_node_id = local_bulk_node_ids_for_subdomain[id];
193 // this puts the global bulk node id at pos:
194 // subdomain_node->getID() - number_of_regular_nodes
195 subdomain_node_id_to_bulk_node_id.push_back(
196 bulk_mesh->getGlobalNodeID(bulk_node_id));
197 }
198 }
199
200 // send ids of bulk ghost nodes belonging to the subdomain mesh to all
201 // other ranks and at the same time receive from all other ranks
202 // first send the sizes to all other to are able to allocate buffer
203 auto const size = subdomain_node_id_to_bulk_node_id.size();
204 std::size_t global_number_of_subdomain_node_id_to_bulk_node_id = 0;
205 MPI_Allreduce(&size, &global_number_of_subdomain_node_id_to_bulk_node_id, 1,
206 MPI_UNSIGNED_LONG, MPI_SUM, mpi_comm);
207
208 DBUG("[{}] global_number_of_subdomain_node_id_to_bulk_node_id: '{}' ",
209 subdomain_mesh->getName(),
210 global_number_of_subdomain_node_id_to_bulk_node_id);
211
212 std::vector<int> numbers_of_ids_at_ranks(mpi_comm_size);
213 MPI_Allgather(&size, 1, MPI_INT, numbers_of_ids_at_ranks.data(), 1, MPI_INT,
214 mpi_comm);
215 std::vector<int> offsets;
216 offsets.push_back(0);
217 std::partial_sum(begin(numbers_of_ids_at_ranks),
218 end(numbers_of_ids_at_ranks), back_inserter(offsets));
219 std::vector<std::size_t> ghost_node_ids_of_all_ranks(
220 global_number_of_subdomain_node_id_to_bulk_node_id);
221 MPI_Allgatherv(subdomain_node_id_to_bulk_node_id.data(), /* sendbuf */
222 size, /* sendcount */
223 MPI_UNSIGNED_LONG, /* sendtype */
224 ghost_node_ids_of_all_ranks.data(), /* recvbuf (out) */
225 numbers_of_ids_at_ranks.data(), /* recvcounts */
226 offsets.data(), /* displs */
227 MPI_UNSIGNED_LONG, /* recvtype */
228 mpi_comm);
229
230 // construct a map for fast search of local bulk node ids
231 std::map<std::size_t, std::size_t> global_to_local_bulk_node_ids;
232 for (auto const id : bulk_mesh->getNodes() | MeshLib::views::ids)
233 {
234 if (!bulk_mesh->isGhostNode(id))
235 {
236 global_to_local_bulk_node_ids[bulk_mesh->getGlobalNodeID(id)] = id;
237 }
238 }
239
240 // construct a map for fast search of local sub-domain node ids
241 std::map<std::size_t, std::size_t> local_subdomain_node_ids;
242 for (auto const id : subdomain_mesh->getNodes() | MeshLib::views::ids)
243 {
244 local_subdomain_node_ids[local_bulk_node_ids_for_subdomain[id]] = id;
245 }
246
247 std::vector<std::size_t> local_subdomain_node_ids_of_all_ranks(
248 global_number_of_subdomain_node_id_to_bulk_node_id,
249 std::numeric_limits<std::size_t>::max());
250 // search in all ranks within the bulk ids for the corresponding id
251 for (int rank = 0; rank < mpi_comm_size; ++rank)
252 {
253 if (rank == mpi_comm_rank)
254 {
255 continue;
256 }
257 for (int i = offsets[rank]; i < offsets[rank + 1]; ++i)
258 {
259 auto it = global_to_local_bulk_node_ids.find(
260 ghost_node_ids_of_all_ranks[i]);
261 if (it == global_to_local_bulk_node_ids.end())
262 {
263 continue;
264 }
265 auto const local_bulk_node_id = it->second;
266 local_subdomain_node_ids_of_all_ranks[i] =
267 global_regular_base_node_ids
268 [local_subdomain_node_ids.find(local_bulk_node_id)->second];
269 DBUG(
270 "[{}] found global subdomain node id: '{}' for global bulk "
271 "node id {} ",
272 subdomain_mesh->getName(),
273 local_subdomain_node_ids_of_all_ranks[i],
274 ghost_node_ids_of_all_ranks[i]);
275 }
276 }
277
278 // send the computed ids back
279 std::vector<std::size_t> computed_global_ids_for_subdomain_ghost_nodes(
280 global_number_of_subdomain_node_id_to_bulk_node_id);
281 MPI_Allreduce(
282 local_subdomain_node_ids_of_all_ranks.data(), /* sendbuf */
283 computed_global_ids_for_subdomain_ghost_nodes
284 .data(), /* recvbuf (out) */
285 global_number_of_subdomain_node_id_to_bulk_node_id, /* sendcount */
286 MPI_UNSIGNED_LONG, /* sendtype */
287 MPI_MAX, /* operation */
288 mpi_comm);
289
290 std::vector<std::size_t> global_ids_for_subdomain_ghost_nodes(
291 computed_global_ids_for_subdomain_ghost_nodes.begin() +
292 offsets[mpi_comm_rank],
293 computed_global_ids_for_subdomain_ghost_nodes.begin() +
294 offsets[mpi_comm_rank + 1]);
295 return global_ids_for_subdomain_ghost_nodes;
296}
297
299 Mesh const* subdomain_mesh)
300{
301 auto const number_of_regular_base_nodes =
302 subdomain_mesh->computeNumberOfBaseNodes();
303
304 MPI_Comm const mpi_comm = MPI_COMM_WORLD;
305 auto const [mpi_comm_rank, mpi_comm_size] = getMPIRankAndSize(mpi_comm);
306
307 std::vector<std::size_t> gathered_number_of_regular_base_nodes(
308 mpi_comm_size);
309 MPI_Allgather(&number_of_regular_base_nodes, 1, MPI_UNSIGNED_LONG,
310 gathered_number_of_regular_base_nodes.data(), 1,
311 MPI_UNSIGNED_LONG, mpi_comm);
312
313 std::vector<std::size_t> numbers_of_regular_base_nodes_at_rank;
314 numbers_of_regular_base_nodes_at_rank.push_back(0);
315 std::partial_sum(begin(gathered_number_of_regular_base_nodes),
316 end(gathered_number_of_regular_base_nodes),
317 back_inserter(numbers_of_regular_base_nodes_at_rank));
318
319 return numbers_of_regular_base_nodes_at_rank;
320}
321
322// similar to the above only with regular higher order nodes
324 Mesh const* subdomain_mesh)
325{
326 // in the following information exchange with other ranks is required
327 MPI_Comm const mpi_comm = MPI_COMM_WORLD;
328 auto [mpi_comm_rank, mpi_comm_size] = getMPIRankAndSize(mpi_comm);
329
330 auto const number_of_regular_base_nodes =
331 subdomain_mesh->computeNumberOfBaseNodes();
332
333 std::vector<std::size_t> gathered_number_of_regular_higher_order_nodes(
334 mpi_comm_size);
335 auto const number_of_regular_higher_order_nodes =
336 subdomain_mesh->getNumberOfNodes() - number_of_regular_base_nodes;
337 MPI_Allgather(&number_of_regular_higher_order_nodes, 1, MPI_UNSIGNED_LONG,
338 gathered_number_of_regular_higher_order_nodes.data(), 1,
339 MPI_UNSIGNED_LONG, mpi_comm);
340
341 std::vector<std::size_t> numbers_of_regular_higher_order_nodes_at_rank;
342 numbers_of_regular_higher_order_nodes_at_rank.push_back(0);
343 std::partial_sum(
344 begin(gathered_number_of_regular_higher_order_nodes),
345 end(gathered_number_of_regular_higher_order_nodes),
346 back_inserter(numbers_of_regular_higher_order_nodes_at_rank));
347
348 return numbers_of_regular_higher_order_nodes_at_rank;
349}
350
352 NodePartitionedMesh const* bulk_mesh, Mesh const* subdomain_mesh)
353{
354 auto global_regular_base_node_ids =
356 subdomain_mesh);
357 auto const local_ids_for_subdomain_ghost_nodes =
359 bulk_mesh, subdomain_mesh, global_regular_base_node_ids);
360
361 // At the moment higher order bulk meshes aren't handled correctly.
362 // If it should become necessary in the future, the necessary things must be
363 // implemented at this point.
364 /*
365 auto const regular_higher_order_node_global_node_ids =
366 computeRegularHigherOrderNodeGlobalNodeIDsofSubDomainPartition(
367 bulk_mesh, subdomain_mesh);
368 auto const ghost_higher_order_node_global_node_ids =
369 computeGhostHigherOrderNodeGlobalNodeIDsofSubDomainPartition(
370 bulk_mesh, subdomain_mesh);
371 */
372
373 std::vector<std::size_t> global_node_ids(
374 std::move(global_regular_base_node_ids));
375 global_node_ids.insert(global_node_ids.end(),
376 local_ids_for_subdomain_ghost_nodes.begin(),
377 local_ids_for_subdomain_ghost_nodes.end());
378
379 return global_node_ids;
380}
381
382unsigned long getNumberOfGlobalNodes(Mesh const* subdomain_mesh)
383{
384 // sum all nodes over all partitions in number_of_global_nodes
385 unsigned long number_of_local_nodes = subdomain_mesh->getNodes().size();
386 unsigned long number_of_global_nodes = 0;
387
388 MPI_Comm mpi_comm = MPI_COMM_WORLD;
389
390 MPI_Allreduce(&number_of_local_nodes, &number_of_global_nodes, 1,
391 MPI_UNSIGNED_LONG, MPI_SUM, mpi_comm);
392 DBUG("[{}] number_of_global_nodes: {}'", subdomain_mesh->getName(),
393 number_of_global_nodes);
394 return number_of_global_nodes;
395}
396
397std::unique_ptr<NodePartitionedMesh> transformMeshToNodePartitionedMesh(
398 NodePartitionedMesh const* const bulk_mesh,
399 Mesh const* const subdomain_mesh)
400{
401 DBUG("Creating NodePartitionedMesh from '{}'", subdomain_mesh->getName());
402
403 auto const subdomain_global_node_ids =
404 computeGlobalNodeIDsOfSubDomainPartition(bulk_mesh, subdomain_mesh);
405
406 // according to comment in struct PartitionedMeshInfo this value
407 // is unused
408 unsigned long number_of_global_base_nodes = 0;
409
410 unsigned long const number_of_global_nodes =
411 getNumberOfGlobalNodes(subdomain_mesh);
412 auto numbers_of_regular_base_nodes_at_rank =
414 auto numbers_of_regular_higher_order_nodes_at_rank =
416
417 auto const number_of_regular_nodes =
418 computeNumberOfRegularNodes(bulk_mesh, subdomain_mesh);
419 DBUG("[{}] number_of_regular_nodes: {}", subdomain_mesh->getName(),
420 number_of_regular_nodes);
421 auto const [nodes, elements] =
422 copyNodesAndElements(subdomain_mesh->getElements());
423
424 return std::make_unique<NodePartitionedMesh>(
425 subdomain_mesh->getName(), nodes, subdomain_global_node_ids, elements,
426 subdomain_mesh->getProperties(), number_of_global_base_nodes,
427 number_of_global_nodes, number_of_regular_nodes,
428 std::move(numbers_of_regular_base_nodes_at_rank),
429 std::move(numbers_of_regular_higher_order_nodes_at_rank));
430}
431
432} // 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.
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)
std::pair< int, int > getMPIRankAndSize(MPI_Comm const &mpi_comm)
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)