55 std::vector<Element*>
const& input_elements)
60 std::unordered_map<std::size_t, Node*> id_node_hash_map;
61 id_node_hash_map.reserve(
64 for (
auto& e : elements)
67 unsigned const n_nodes = e->getNumberOfNodes();
68 for (
unsigned i = 0; i < n_nodes; ++i)
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())
74 auto new_node_in_map = id_node_hash_map[n->
getID()] =
76 e->setNode(i, new_node_in_map);
80 e->setNode(i, it->second);
85 std::map<std::size_t, Node*> nodes_map;
86 for (
auto const& n : id_node_hash_map)
88 nodes_map[n.first] = n.second;
93 nodes_map | ranges::views::values | ranges::to<std::vector>;
95 return std::make_pair(element_nodes, elements);
121 auto const& subdomain_nodes = subdomain_mesh->
getNodes();
122 auto const& local_bulk_node_ids_for_subdomain =
136 std::size_t
const number_of_regular_nodes =
140 MPI_Comm mpi_comm = MPI_COMM_WORLD;
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,
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));
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(),
165 if (isRegularNode(*bulk_mesh, local_bulk_node_ids_for_subdomain,
id))
167 subdomain_global_node_ids.emplace_back(partition_offset +
id);
170 return subdomain_global_node_ids;
175 Mesh const* subdomain_mesh,
176 std::vector<std::size_t>
const& global_regular_base_node_ids)
179 MPI_Comm
const mpi_comm = MPI_COMM_WORLD;
183 auto const local_bulk_node_ids_for_subdomain =
187 std::vector<std::size_t> subdomain_node_id_to_bulk_node_id;
190 if (!isRegularNode(*bulk_mesh, local_bulk_node_ids_for_subdomain,
id))
192 auto const bulk_node_id = local_bulk_node_ids_for_subdomain[id];
195 subdomain_node_id_to_bulk_node_id.push_back(
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);
208 DBUG(
"[{}] global_number_of_subdomain_node_id_to_bulk_node_id: '{}' ",
210 global_number_of_subdomain_node_id_to_bulk_node_id);
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,
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(),
224 ghost_node_ids_of_all_ranks.data(),
225 numbers_of_ids_at_ranks.data(),
231 std::map<std::size_t, std::size_t> global_to_local_bulk_node_ids;
241 std::map<std::size_t, std::size_t> local_subdomain_node_ids;
244 local_subdomain_node_ids[local_bulk_node_ids_for_subdomain[id]] = id;
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());
251 for (
int rank = 0; rank < mpi_comm_size; ++rank)
253 if (rank == mpi_comm_rank)
257 for (
int i = offsets[rank]; i < offsets[rank + 1]; ++i)
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())
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];
270 "[{}] found global subdomain node id: '{}' for global bulk "
273 local_subdomain_node_ids_of_all_ranks[i],
274 ghost_node_ids_of_all_ranks[i]);
279 std::vector<std::size_t> computed_global_ids_for_subdomain_ghost_nodes(
280 global_number_of_subdomain_node_id_to_bulk_node_id);
282 local_subdomain_node_ids_of_all_ranks.data(),
283 computed_global_ids_for_subdomain_ghost_nodes
285 global_number_of_subdomain_node_id_to_bulk_node_id,
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;
299 Mesh const* subdomain_mesh)
301 auto const number_of_regular_base_nodes =
304 MPI_Comm
const mpi_comm = MPI_COMM_WORLD;
307 std::vector<std::size_t> gathered_number_of_regular_base_nodes(
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);
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));
319 return numbers_of_regular_base_nodes_at_rank;
324 Mesh const* subdomain_mesh)
327 MPI_Comm
const mpi_comm = MPI_COMM_WORLD;
330 auto const number_of_regular_base_nodes =
333 std::vector<std::size_t> gathered_number_of_regular_higher_order_nodes(
335 auto const number_of_regular_higher_order_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);
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);
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));
348 return numbers_of_regular_higher_order_nodes_at_rank;
399 Mesh const*
const subdomain_mesh)
401 DBUG(
"Creating NodePartitionedMesh from '{}'", subdomain_mesh->
getName());
403 auto const subdomain_global_node_ids =
408 unsigned long number_of_global_base_nodes = 0;
410 unsigned long const number_of_global_nodes =
412 auto numbers_of_regular_base_nodes_at_rank =
414 auto numbers_of_regular_higher_order_nodes_at_rank =
417 auto const number_of_regular_nodes =
419 DBUG(
"[{}] number_of_regular_nodes: {}", subdomain_mesh->
getName(),
420 number_of_regular_nodes);
421 auto const [nodes, elements] =
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));