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>
36 std::vector<std::size_t>
const& local_bulk_node_ids_for_subdomain,
37 std::size_t
const subdomain_node_id)
40 local_bulk_node_ids_for_subdomain[subdomain_node_id]));
47 std::vector<Element*>
const& input_elements)
52 std::unordered_map<std::size_t, Node*> id_node_hash_map;
53 id_node_hash_map.reserve(
56 for (
auto& e : elements)
59 unsigned const n_nodes = e->getNumberOfNodes();
60 for (
unsigned i = 0; i < n_nodes; ++i)
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())
66 auto new_node_in_map = id_node_hash_map[n->
getID()] =
68 e->setNode(i, new_node_in_map);
72 e->setNode(i, it->second);
77 std::map<std::size_t, Node*> nodes_map;
78 for (
auto const& n : id_node_hash_map)
80 nodes_map[n.first] = n.second;
85 nodes_map | ranges::views::values | ranges::to<std::vector>;
87 return std::make_pair(element_nodes, elements);
91 Mesh const* subdomain_mesh)
93 auto const& subdomain_nodes = subdomain_mesh->
getNodes();
94 auto const& local_bulk_node_ids_for_subdomain =
96 unsigned long const number_of_regular_nodes = ranges::count_if(
98 [&](std::size_t
const id) {
99 return isRegularNode(*bulk_mesh, local_bulk_node_ids_for_subdomain,
103 DBUG(
"[{}] number of regular nodes: {}", subdomain_mesh->
getName(),
104 number_of_regular_nodes);
105 return number_of_regular_nodes;
108std::vector<std::size_t>
113 auto const& subdomain_nodes = subdomain_mesh->
getNodes();
114 auto const& local_bulk_node_ids_for_subdomain =
128 std::size_t
const number_of_regular_nodes =
134 std::vector<std::size_t>
const gathered_number_of_regular_nodes =
138 std::vector<std::size_t>
const numbers_of_regular_nodes_at_rank =
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(),
150 if (isRegularNode(*bulk_mesh, local_bulk_node_ids_for_subdomain,
id))
152 subdomain_global_node_ids.emplace_back(partition_offset +
id);
155 return subdomain_global_node_ids;
160 Mesh const* subdomain_mesh,
161 std::vector<std::size_t>
const& global_regular_base_node_ids)
166 auto const local_bulk_node_ids_for_subdomain =
170 std::vector<std::size_t> subdomain_node_id_to_bulk_node_id;
173 if (!isRegularNode(*bulk_mesh, local_bulk_node_ids_for_subdomain,
id))
175 auto const bulk_node_id = local_bulk_node_ids_for_subdomain[id];
178 subdomain_node_id_to_bulk_node_id.push_back(
185 std::vector<std::size_t> ghost_node_ids_of_all_ranks;
186 std::vector<int>
const offsets =
188 ghost_node_ids_of_all_ranks, mpi);
190 std::size_t
const global_number_of_subdomain_node_id_to_bulk_node_id =
192 DBUG(
"[{}] global_number_of_subdomain_node_id_to_bulk_node_id: '{}' ",
194 global_number_of_subdomain_node_id_to_bulk_node_id);
197 std::map<std::size_t, std::size_t> global_to_local_bulk_node_ids;
207 std::map<std::size_t, std::size_t> local_subdomain_node_ids;
210 local_subdomain_node_ids[local_bulk_node_ids_for_subdomain[id]] = id;
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());
217 for (
int rank = 0; rank < mpi.
size; ++rank)
219 if (rank == mpi.
rank)
223 for (
int i = offsets[rank]; i < offsets[rank + 1]; ++i)
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())
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];
235 "[{}] found global subdomain node id: '{}' for global bulk "
238 subdomain_node_ids_of_all_ranks[i],
239 ghost_node_ids_of_all_ranks[i]);
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;
254 Mesh const* subdomain_mesh)
256 auto const number_of_regular_base_nodes =
261 std::vector<std::size_t>
const gathered_number_of_regular_base_nodes =
269 Mesh const* subdomain_mesh)
273 auto const number_of_regular_base_nodes =
276 auto const number_of_regular_higher_order_nodes =
278 std::vector<std::size_t> gathered_number_of_regular_higher_order_nodes =
282 gathered_number_of_regular_higher_order_nodes);
288 auto global_regular_base_node_ids =
291 auto const local_ids_for_subdomain_ghost_nodes =
293 bulk_mesh, subdomain_mesh, global_regular_base_node_ids);
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());
313 return global_node_ids;
321 DBUG(
"[{}] number_of_global_nodes: {}'", subdomain_mesh->
getName(),
322 number_of_global_nodes);
323 return number_of_global_nodes;
328 Mesh const*
const subdomain_mesh)
330 DBUG(
"Creating NodePartitionedMesh from '{}'", subdomain_mesh->
getName());
332 auto const subdomain_global_node_ids =
337 unsigned long number_of_global_base_nodes = 0;
339 unsigned long const number_of_global_nodes =
341 auto numbers_of_regular_base_nodes_at_rank =
343 auto numbers_of_regular_higher_order_nodes_at_rank =
346 auto const number_of_regular_nodes =
348 DBUG(
"[{}] number_of_regular_nodes: {}", subdomain_mesh->
getName(),
349 number_of_regular_nodes);
350 auto const [nodes, elements] =
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));
Definition of Duplicate functions.
Definition of the Element class.
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
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.
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
std::size_t computeNumberOfBaseNodes() const
Get the number of base nodes.
Properties & getProperties()
const std::string getName() const
Get name of the mesh.
std::size_t getNumberOfNodes() const
Get the number of nodes.
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)
static std::vector< int > allgatherv(std::span< T > const send_buffer, std::vector< std::remove_const_t< T > > &receive_buffer, Mpi const &mpi)
static T allreduce(T const &value, MPI_Op const &mpi_op, Mpi const &mpi)
static std::vector< T > allgather(T const &value, Mpi const &mpi)
std::vector< ranges::range_value_t< R > > sizesToOffsets(R const &sizes)
constexpr ranges::views::view_closure ids
For an element of a range view return its id.
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)
std::unique_ptr< NodePartitionedMesh > transformMeshToNodePartitionedMesh(NodePartitionedMesh const *const bulk_mesh, Mesh const *const subdomain_mesh)