openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
SF_mesh_utils_emi.h
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 // openCARP is an open cardiac electrophysiology simulator.
3 //
4 // Copyright (C) 2020 openCARP project
5 //
6 // This program is licensed under the openCARP Academic Public License (APL)
7 // v1.0: You can use and redistribute it and/or modify it in non-commercial
8 // academic environments under the terms of APL as published by the openCARP
9 // project v1.0, or (at your option) any later version. Commercial use requires
10 // a commercial license (info@opencarp.org).
11 //
12 // This program is distributed without any warranty; see the openCARP APL for
13 // more details.
14 //
15 // You should have received a copy of the openCARP APL along with this program
16 // and can find it online: http://www.opencarp.org/license
17 // ----------------------------------------------------------------------------
18 
26 #if WITH_EMI_MODEL
27 
28 #ifndef _SF_MESH_UTILS_EMI_H
29 #define _SF_MESH_UTILS_EMI_H
30 
31 #include <algorithm>
32 #include <cstdio>
33 #include <cstdlib>
34 
35 #include "SF_mesh_utils.h"
36 #include "SF_mesh_io_emi.h"
37 
38 namespace SF {
39 
40 
41 const int MAX_INTERSECTIONS = 20;
42 // Adaptive merge-rank settings for EMI helper exchanges.
43 // For communicator sizes up to the node-aware cap, we preserve the original
44 // all-ranks behavior. Once size > merge_rank_cap, we reduce the number of
45 // merge ranks based on the estimated global communication volume, bounded by
46 // that cap, and spread those merge ranks across the communicator. The lower
47 // bound avoids concentrating all merge work on too few ranks. These values are
48 // heuristics and may need retuning on different machines.
49 const int EMI_MIN_MERGE_RANKS = 8;
50 const int EMI_MAX_MERGE_RANKS_PER_NODE = 2;
51 const unsigned long long EMI_TARGET_BYTES_PER_MERGE_RANK = 8ull * 1024ull * 1024ull;
52 
53 // Return the maximum number of merge ranks allowed for this communicator.
54 // We estimate the number of participating nodes from MPI shared-memory subcommunicators
55 // and allow only a small number of merge ranks per node. This keeps merge ranks
56 // distributed across the machine instead of concentrating them on a few low ranks.
57 inline int emi_node_aware_merge_rank_cap(MPI_Comm comm)
58 {
59  int size = 0;
60  MPI_Comm_size(comm, &size);
61 
62  // Split the communicator into shared-memory subcommunicators, one per node.
63  MPI_Comm node_comm = MPI_COMM_NULL;
64  MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, &node_comm);
65 
66  // Rank 0 in each node subcommunicator acts as the node representative.
67  int node_rank = 0;
68  MPI_Comm_rank(node_comm, &node_rank);
69 
70  int local_leader = (node_rank == 0) ? 1 : 0;
71  int num_nodes = 0;
72  // Count how many node representatives exist globally, which gives the number of nodes.
73  MPI_Allreduce(&local_leader, &num_nodes, 1, MPI_INT, MPI_SUM, comm);
74 
75  MPI_Comm_free(&node_comm);
76 
77  // Allow only a small number of merge ranks per node, bounded by the communicator size.
78  return std::min(size, std::max(EMI_MIN_MERGE_RANKS, num_nodes * EMI_MAX_MERGE_RANKS_PER_NODE));
79 }
80 
81 // Build the destination ranks for the merge step used by EMI MPI_Exchange()-based helpers.
82 // The only requirement of the algorithm is that every identical key is sent to the
83 // same merge rank, where the partial data are combined, and then sent back. For
84 // communicator sizes up to the node-aware cap, we preserve the original all-ranks
85 // behavior. Once size > merge_rank_cap, we reduce the number of merge ranks based
86 // on the estimated global communication volume, bounded by that cap, and spread
87 // those merge ranks across the communicator.
88 template<class K, class TokenFn>
89 inline int emi_select_merge_destinations(const vector<K>& key_vec,
90  size_t dsize,
91  size_t bytes_per_entry,
92  vector<int>& dest,
93  MPI_Comm comm,
94  const char* label,
95  TokenFn token_fn)
96 {
97  int size = 0, rank = 0;
98  MPI_Comm_size(comm, &size);
99  MPI_Comm_rank(comm, &rank);
100 
101  // Estimate how many bytes this helper exchanges locally and globally.
102  const unsigned long long local_bytes =
103  static_cast<unsigned long long>(dsize) * static_cast<unsigned long long>(bytes_per_entry);
104  unsigned long long global_bytes = 0;
105  MPI_Allreduce(&local_bytes, &global_bytes, 1, MPI_UNSIGNED_LONG_LONG, MPI_SUM, comm);
106 
107  // Compute the node-aware upper bound on how many merge ranks we allow.
108  const int merge_rank_cap = emi_node_aware_merge_rank_cap(comm);
109  int nmerge = size;
110  if (size > merge_rank_cap) {
111  // For large runs, estimate how many merge ranks are needed from the total data volume.
112  const int nmerge_floor = EMI_MIN_MERGE_RANKS;
113  const int nmerge_cap = merge_rank_cap;
114  const unsigned long long nmerge_est_ull =
115  std::max(1ull, (global_bytes + EMI_TARGET_BYTES_PER_MERGE_RANK - 1ull) /
116  EMI_TARGET_BYTES_PER_MERGE_RANK);
117  const int nmerge_est = static_cast<int>(std::min<unsigned long long>(nmerge_cap, nmerge_est_ull));
118  // Clamp the estimate to the allowed range for this communicator.
119  nmerge = std::max(nmerge_floor, nmerge_est);
120  }
121 
122  dest.resize(dsize);
123  for (size_t i = 0; i < dsize; ++i) {
124  if (nmerge == size) {
125  // Preserve the original all-ranks merge behavior for small and medium runs.
126  dest[i] = static_cast<int>(token_fn(key_vec[i]) % static_cast<size_t>(size));
127  } else {
128  // First pick one logical merge bucket in [0, nmerge).
129  const size_t bucket = token_fn(key_vec[i]) % static_cast<size_t>(nmerge);
130  // Then spread those buckets over the full communicator instead of packing them
131  // into ranks [0, nmerge), which reduces concentration on low ranks.
132  dest[i] = static_cast<int>((bucket * static_cast<size_t>(size)) / static_cast<size_t>(nmerge));
133  }
134  }
135 
136  if (rank == 0 && std::getenv("OPENCARP_EMI_LOG_MERGE_RANKS") != NULL) {
137  const char* mode = (nmerge == size) ? "all-ranks" : "limited";
138  std::fprintf(stdout,
139  "[EMI merge ranks] %s: size=%d dsize=%zu global_bytes=%llu cap=%d nmerge=%d mode=%s\n",
140  label, size, dsize, global_bytes, merge_rank_cap, nmerge, mode);
141  std::fflush(stdout);
142  }
143 
144  return nmerge;
145 }
146 
155 struct intersection_tags {
163  int tags[MAX_INTERSECTIONS];
164  intersection_tags() {
165  std::fill_n(tags, MAX_INTERSECTIONS, -1);
166  }
167 
168 };
169 
170 struct intersection_data {
171  int tags[MAX_INTERSECTIONS];
172  int data[MAX_INTERSECTIONS];
173  int ranks[MAX_INTERSECTIONS];
174 
175  intersection_data() {
176  std::fill_n(tags, MAX_INTERSECTIONS, -1);
177  std::fill_n(data, MAX_INTERSECTIONS, -1);
178  std::fill_n(ranks, MAX_INTERSECTIONS, -1);
179  }
180 };
181 
182 struct intersection_indices {
183  int indices[MAX_INTERSECTIONS];
184  intersection_indices() {
185  std::fill_n(indices, MAX_INTERSECTIONS, -1);
186  }
187 };
196 template<class K, class V> inline
197 void assign_counter_face(hashmap::unordered_map<K, V> & map, const MPI_Comm comm)
198 {
199  size_t dsize = map.size();
200  vector<K> key_vec (dsize);
201  vector<V> value_vec (dsize);
202 
203  // make key and value vector of elements without counterparts
204  // which are from different ranks
205  size_t idx=0;
206  for(const auto & v : map) {
207  if (v.second.second.eidx == -1){
208  key_vec[idx] = v.first;
209  value_vec[idx] = v.second;
210  idx++;
211  }
212  }
213  dsize = idx;
214  key_vec.resize(dsize);
215  value_vec.resize(dsize);
216 
217  vector<int> perm, dest;
218  emi_select_merge_destinations(key_vec, dsize, sizeof(K) + sizeof(V), dest, comm,
219  "assign_counter_face",
220  [](const K& key) { return hashmap::hash_ops<K>::hash(key); });
221 
222  commgraph<size_t> grph;
223  grph.configure(dest, comm);
224  size_t nrecv = sum(grph.rcnt);
225 
226  interval(perm, 0, dsize);
227  binary_sort_copy(dest, perm);
228 
229  // fill send buffer and communicate
230  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
231  vector<V> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
232  for(size_t i=0; i<dsize; i++) {
233  sbuff_key[i] = key_vec[perm[i]];
234  sbuff_value[i] = value_vec[perm[i]];
235  }
236  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
237  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
238 
239  // Merge values for identical keys on the receiving side.
240  // Each merged entry should contain both sides and preserve owner index_unique.
242 
243  // Merge all received entries for the same key.
244  for (size_t i = 0; i < nrecv; i++) {
245  auto it = rmap.find(rbuff_key[i]);
246  if(it != rmap.end()) {
247  // Add counter faces only if the tag numbers of the two faces are different.
248  // That is, the counter face must have a different tag number from the current face.
249  // However, in some special cases, both faces may belong to the extracellular domain,
250  // which are removed later.
251  if (it->second.first.tag != rbuff_value[i].first.tag) {
252  it->second.second = rbuff_value[i].first;
253  }
254  } else {
255  rmap.insert({rbuff_key[i], rbuff_value[i]});
256  }
257  }
258 
259  for (size_t i = 0; i < nrecv; i++ ) {
260  auto it = rmap.find(rbuff_key[i]);
261  if(it != rmap.end()) rbuff_value[i] = it->second;
262  }
263 
264  // Send merged values back to the originating ranks.
265  grph.transpose();
266  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
267 
268  for (size_t i = 0; i < dsize; i++ ) {
269  auto it = map.find(sbuff_key[i]);
270  if (it != map.end()) it->second = ibuff_value[i];
271  }
272 }
273 
282 template<class K, class V> inline
283 void assign_unique_first_face(hashmap::unordered_map<K, V> & map, const MPI_Comm comm)
284 {
285  // Exchange and merge interface-face pairs across ranks.
286  // Each key represents one physical interface; after this call, both sides
287  // (first/second) carry rank/index data, and the owner side has index_unique set.
288  //
289  // Steps:
290  // 1) Collect only partial entries (missing counter-face) that need MPI exchange.
291  // 2) Hash-partition these keys to a merge owner rank.
292  // 3) Exchange keys/values to owners and merge both sides for each key.
293  // 4) Send merged entries back so every rank has complete data for its keys.
294  size_t dsize = map.size();
295  vector<K> key_vec (dsize);
296  vector<V> value_vec (dsize);
297 
298  // Collect only entries without counterparts (cross-rank interfaces).
299  size_t idx=0;
300  for(const auto & v : map) {
301  if (v.second.second.index_unique == -1 && v.second.second.index_one == -1){
302  key_vec[idx] = v.first;
303  value_vec[idx] = v.second;
304  idx++;
305  }
306  }
307  dsize = idx;
308  key_vec.resize(dsize);
309  value_vec.resize(dsize);
310 
311  vector<int> perm, dest;
312  emi_select_merge_destinations(key_vec, dsize, sizeof(K) + sizeof(V), dest, comm,
313  "assign_unique_first_face",
314  [](const K& key) { return hashmap::hash_ops<K>::hash(key); });
315 
316  commgraph<size_t> grph;
317  grph.configure(dest, comm);
318  size_t nrecv = sum(grph.rcnt);
319 
320  interval(perm, 0, dsize);
321  binary_sort_copy(dest, perm);
322 
323  // Send keys and partial values to the owner rank for each key.
324  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
325  vector<V> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
326  for(size_t i=0; i<dsize; i++) {
327  sbuff_key[i] = key_vec[perm[i]];
328  sbuff_value[i] = value_vec[perm[i]];
329  }
330  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
331  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
332 
333  // rmap accumulates merged entries per key on the receiving rank.
334  // Multiple ranks may send the same key; we combine them into one pair.
336 
337  // Merge all received entries for the same key.
338  for (size_t i = 0; i < nrecv; i++) {
339  auto it = rmap.find(rbuff_key[i]);
340  if(it != rmap.end()) {
341  // Merge entries without losing existing owner info.
342  auto& existing = it->second;
343  using Face = std::decay_t<decltype(existing.first)>;
344 
345  auto update_slot = [&](Face& slot, const Face& src) {
346  // Merge strategy for one side:
347  // - index_unique is only set by the owner and must never be overwritten by a non-owner.
348  // - rank/index_one are filled lazily when missing.
349  if (src.index_unique >= 0) {
350  slot.index_unique = src.index_unique;
351  slot.rank = src.rank;
352  } else {
353  if (slot.rank < 0) slot.rank = src.rank;
354  }
355  if (slot.index_one < 0 && src.index_one >= 0) slot.index_one = src.index_one;
356  };
357 
358  auto same_face = [&](const Face& a, const Face& b) {
359  // Use (rank,index_one) to identify a side uniquely (index_one is local).
360  return (a.rank >= 0 && b.rank >= 0 &&
361  a.index_one >= 0 && b.index_one >= 0 &&
362  a.rank == b.rank && a.index_one == b.index_one);
363  };
364 
365  auto ingest = [&](const Face& incoming) {
366  // Ignore empty placeholders (no useful info).
367  if (incoming.index_one < 0 && incoming.index_unique < 0 && incoming.rank < 0) return;
368  // Try to match incoming with an existing side (same rank or same index_one).
369  if (same_face(existing.first, incoming) || existing.first.rank == incoming.rank) {
370  update_slot(existing.first, incoming);
371  } else if (same_face(existing.second, incoming) || existing.second.rank == incoming.rank) {
372  update_slot(existing.second, incoming);
373  } else if (existing.second.index_one == -1 && existing.second.index_unique == -1) {
374  // If the second slot is empty, place incoming there (new side).
375  existing.second = incoming;
376  } else if (existing.first.index_one == -1 && existing.first.index_unique == -1) {
377  // If the first slot is empty, place incoming there (new side).
378  existing.first = incoming;
379  } else {
380  // Both slots are filled. Only upgrade missing owner info.
381  // (Do not overwrite an existing owner assignment.)
382  if (existing.first.index_unique < 0 && incoming.index_unique >= 0) {
383  existing.first.index_unique = incoming.index_unique;
384  existing.first.rank = incoming.rank;
385  } else if (existing.second.index_unique < 0 && incoming.index_unique >= 0) {
386  existing.second.index_unique = incoming.index_unique;
387  existing.second.rank = incoming.rank;
388  }
389  }
390  };
391 
392  // Merge both sides from the received pair.
393  ingest(rbuff_value[i].first);
394  ingest(rbuff_value[i].second);
395  } else {
396  // First time we see this key on the receiving rank.
397  rmap.insert({rbuff_key[i], rbuff_value[i]});
398  }
399  }
400 
401  for (size_t i = 0; i < nrecv; i++ ) {
402  auto it = rmap.find(rbuff_key[i]);
403  if(it != rmap.end()) rbuff_value[i] = it->second;
404  }
405 
406  // sending the assigned data back to original rank
407  grph.transpose();
408  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
409 
410  for (size_t i = 0; i < dsize; i++ ) {
411  auto it = map.find(sbuff_key[i]);
412  if (it != map.end()) {
413  // Overwrite with merged entry from exchange (rmap already combined both sides)
414  it->second = ibuff_value[i];
415  } else {
416  // Insert new entry received from another rank
417  // This ensures non-owner ranks get the complete mapping data
418  map.insert({sbuff_key[i], ibuff_value[i]});
419  }
420  }
421 }
422 
430 template<class K> inline
431 void assign_counter_vertices_tuple(hashmap::unordered_map<K, intersection_data>& map, const MPI_Comm comm)
432 {
433  size_t dsize = map.size();
434  vector<K> key_vec(dsize);
435  vector<intersection_data> value_vec(dsize);
436 
437  size_t idx = 0;
438  for (const auto& v : map) {
439  key_vec[idx] = v.first;
440  value_vec[idx] = v.second;
441  idx++;
442  }
443 
444  vector<int> perm, dest;
445  emi_select_merge_destinations(key_vec, dsize, sizeof(K) + sizeof(intersection_data), dest, comm,
446  "assign_counter_vertices_tuple",
447  [](const K& key) { return hashmap::hash_ops<K>::hash(key); });
448 
449  commgraph<size_t> grph;
450  grph.configure(dest, comm);
451  size_t nrecv = sum(grph.rcnt);
452 
453  interval(perm, 0, dsize);
454  binary_sort_copy(dest, perm);
455 
456  // fill send buffer and communicate
457  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
458  vector<intersection_data> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
459  for (size_t i = 0; i < dsize; i++) {
460  sbuff_key[i] = key_vec[perm[i]];
461  sbuff_value[i] = value_vec[perm[i]];
462  }
463  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
464  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
465 
467  for (size_t i = 0; i < nrecv; i++) {
468  auto it = rmap.find(rbuff_key[i]);
469  if (it != rmap.end()) {
470  intersection_data& map_val = it->second;
471  intersection_data& r_val = rbuff_value[i];
472 
473  int first_index = -1;
474  for (int j = 0; j < MAX_INTERSECTIONS; ++j) {
475  if (map_val.tags[j] == -1 && map_val.ranks[j] == -1) {
476  first_index = j;
477  break;
478  }
479  }
480 
481  if (first_index == -1) continue; // This means all possible intersections have already been assigned.
482 
483  for (int j = 0; j < MAX_INTERSECTIONS; ++j) {
484  if (r_val.tags[j] == -1) continue;
485 
486  int rtag = r_val.tags[j];
487  int rdata = r_val.data[j];
488  int rrank = r_val.ranks[j];
489 
490  bool exist_rtag = false;
491  for (int k = 0; k < first_index; ++k) {
492  if (map_val.tags[k] == rtag && map_val.ranks[k] == rrank) {
493  exist_rtag = true;
494  break;
495  }
496  }
497 
498  if (!exist_rtag && first_index < MAX_INTERSECTIONS) {
499  map_val.tags[first_index] = rtag;
500  map_val.data[first_index] = rdata;
501  map_val.ranks[first_index] = rrank;
502  first_index++;
503  }
504  }
505  }
506  else {
507  rmap.insert({ rbuff_key[i], rbuff_value[i] });
508  }
509  }
510 
511  for (size_t i = 0; i < nrecv; i++) {
512  auto it = rmap.find(rbuff_key[i]);
513  if (it != rmap.end()) rbuff_value[i] = it->second;
514  }
515 
516  // sending the assigned data back to original rank
517  grph.transpose();
518  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
519 
520  for (size_t i = 0; i < dsize; i++) {
521  auto it = map.find(sbuff_key[i]);
522  if (it != map.end()) it->second = ibuff_value[i];
523  }
524 }
525 
533 template<class K, class V> inline
534 void assign_dof_on_counter_face(hashmap::unordered_map<K, V> & map, const MPI_Comm comm)
535 {
536  size_t dsize = map.size();
537  vector<K> key_vec (dsize);
538  vector<V> value_vec (dsize);
539 
540  // make key and value vector of elements without counterparts
541  // which are from different ranks
542  size_t idx=0;
543  for(const auto & v : map) {
544  key_vec[idx] = v.first;
545  value_vec[idx] = v.second;
546  idx++;
547  }
548  dsize = idx;
549  key_vec.resize(dsize);
550  value_vec.resize(dsize);
551 
552  vector<int> perm, dest;
553  emi_select_merge_destinations(key_vec, dsize, sizeof(K) + sizeof(V), dest, comm,
554  "assign_dof_on_counter_face",
555  [](const K& key) { return key.first; });
556 
557  commgraph<size_t> grph;
558  grph.configure(dest, comm);
559  size_t nrecv = sum(grph.rcnt);
560 
561  interval(perm, 0, dsize);
562  binary_sort_copy(dest, perm);
563 
564  // fill send buffer and communicate
565  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
566  vector<V> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
567  for(size_t i=0; i<dsize; i++) {
568  sbuff_key[i] = key_vec[perm[i]];
569  sbuff_value[i] = value_vec[perm[i]];
570  }
571  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
572  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
573 
575 
576  for (size_t i = 0; i < nrecv; i++) {
577  auto it = rmap.find(rbuff_key[i]);
578  if(it != rmap.end()) {
579  if(rbuff_value[i] != -1 and it->second == -1){
580  it->second = rbuff_value[i];// add new indices
581  }
582  } else {
583  rmap.insert({rbuff_key[i], rbuff_value[i]});
584  }
585  }
586 
587  for (size_t i = 0; i < nrecv; i++ ) {
588  auto it = rmap.find(rbuff_key[i]);
589  if(it != rmap.end()) {
590  if(rbuff_value[i] == -1)
591  rbuff_value[i] = it->second;
592  }
593  }
594 
595  // sending the assigned data back to original rank
596  grph.transpose();
597  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
598 
599  for (size_t i = 0; i < dsize; i++ ) {
600  auto it = map.find(sbuff_key[i]);
601  if (it != map.end()) it->second = ibuff_value[i];
602  }
603 }
604 
613 template<class K, class V> inline
614 void assign_petsc_on_counter_face(hashmap::unordered_map<K, V> & map, const MPI_Comm comm)
615 {
616  size_t dsize = map.size();
617  vector<K> key_vec (dsize);
618  vector<V> value_vec (dsize);
619 
620  // make key and value vector of elements without counterparts
621  // which are from different ranks
622  size_t idx=0;
623  for(const auto & v : map) {
624  key_vec[idx] = v.first;
625  value_vec[idx] = v.second;
626  idx++;
627  }
628  dsize = idx;
629  key_vec.resize(dsize);
630  value_vec.resize(dsize);
631 
632  vector<int> perm, dest;
633  emi_select_merge_destinations(key_vec, dsize, sizeof(K) + sizeof(V), dest, comm,
634  "assign_petsc_on_counter_face",
635  [](const K& key) { return key.first; });
636 
637  commgraph<size_t> grph;
638  grph.configure(dest, comm);
639  size_t nrecv = sum(grph.rcnt);
640 
641  interval(perm, 0, dsize);
642  binary_sort_copy(dest, perm);
643 
644  // fill send buffer and communicate
645  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
646  vector<V> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
647  for(size_t i=0; i<dsize; i++) {
648  sbuff_key[i] = key_vec[perm[i]];
649  sbuff_value[i] = value_vec[perm[i]];
650  }
651  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
652  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
653 
655 
656  for (size_t i = 0; i < nrecv; i++) {
657  auto it = rmap.find(rbuff_key[i]);
658  if(it != rmap.end()) {
659  if(rbuff_value[i].second != -1 and it->second.second == -1){
660  it->second = rbuff_value[i];// add new indices
661  }
662  } else {
663  rmap.insert({rbuff_key[i], rbuff_value[i]});
664  }
665  }
666 
667  for (size_t i = 0; i < nrecv; i++ ) {
668  auto it = rmap.find(rbuff_key[i]);
669  if(it != rmap.end()) {
670  if(rbuff_value[i].second == -1)
671  rbuff_value[i].second = it->second.second;
672  }
673  }
674 
675  // sending the assigned data back to original rank
676  grph.transpose();
677  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
678 
679  for (size_t i = 0; i < dsize; i++ ) {
680  auto it = map.find(sbuff_key[i]);
681  if (it != map.end()) it->second = ibuff_value[i];
682  }
683 }
684 
692 template<class T>
693 void sort_surf_local_indices(tuple<T> & ref, tuple<T> & target)
694 {
695  vector<T> buff(2);
696 
697  buff[0] = ref.v1, buff[1] = ref.v2;
698  binary_sort(buff);
699 
700  target.v1 = buff[0], target.v2 = buff[1];
701 }
702 
710 template<class T>
711 void sort_surf_local_indices(triple<T> & ref, triple<T> & target)
712 {
713  vector<T> buff(3);
714 
715  buff[0] = ref.v1, buff[1] = ref.v2, buff[2] = ref.v3;
716  binary_sort(buff);
717 
718  target.v1 = buff[0], target.v2 = buff[1], target.v3 = buff[2];
719 }
720 
728 template<class T>
729 void sort_surf_local_indices(quadruple<T> & ref, quadruple<T> & target)
730 {
731  vector<T> buff(4);
732 
733  buff[0] = ref.v1, buff[1] = ref.v2, buff[2] = ref.v3, buff[3] = ref.v4;
734  binary_sort(buff);
735 
736  target.v1 = buff[0], target.v2 = buff[1], target.v3 = buff[2], target.v4 = buff[3];
737 }
738 
746 template<class K, class V> inline
747 void insert_surf_based_Tag(V & ref,
748  hashmap::unordered_map<K, std::pair<V, V>> & surfmap)
749 {
750  K surf;
751 
752  sort_surf_local_indices(ref.points, surf);
753 
754  auto it = surfmap.find(surf);
755  if (it != surfmap.end()) {
756  // Add as a counter face only if this face exists in surfMap and its tag number
757  // differs from the tag number of the first saved pair.
758  if ( it->second.first.tag != ref.tag )
759  it->second.second = ref;
760  } else {
761  std::pair<V, V> face;
762  face.first = ref;
763  surfmap.insert({surf,face});
764  }
765 }
766 
786 template<class T, class W, class V, class U> inline
787 void insert_surf_emi( int rank,
788  SF::vector<T> const & ref_con,
789  SF::vector<T> const & ptsData,
790  const T eidx, const T tag,
791  const std::vector<T> & line_con, const std::vector<T> & surf_con, const std::vector<T> & qsurf_con,
792  hashmap::unordered_map<tuple<T>, std::pair<W, W> > & line_face,
793  hashmap::unordered_map<triple<T>, std::pair<V, V> > & surfmap,
794  hashmap::unordered_map<quadruple<T>, std::pair<U, U> > & qsurfmap)
795 {
796  W line;
797  line.eidx = eidx;
798  line.tag = tag;
799  V face;
800  face.eidx = eidx;
801  face.tag = tag;
802  U qface;
803  qface.eidx = eidx;
804  qface.tag = tag;
805 
806  // Mark faces as potential candidates if all of their vertices lie on the interface.
807  for (size_t i = 0; i < line_con.size(); i += 2) {
808  if (ptsData[line_con[i ] - 1] > 0 &&
809  ptsData[line_con[i + 1] - 1] > 0 ) {
810  line.points.v1 = ref_con[line_con[i ] - 1];
811  line.points.v2 = ref_con[line_con[i + 1] - 1];
812  line.rank = rank;
813  insert_surf_based_Tag(line, line_face);
814  }
815  }
816 
817  for (size_t i = 0; i < surf_con.size(); i += 3) {
818  if (ptsData[surf_con[i ] - 1] > 0 &&
819  ptsData[surf_con[i + 1] - 1] > 0 &&
820  ptsData[surf_con[i + 2] - 1] > 0 ) {
821  face.points.v1 = ref_con[surf_con[i ] - 1];
822  face.points.v2 = ref_con[surf_con[i + 1] - 1];
823  face.points.v3 = ref_con[surf_con[i + 2] - 1];
824  face.rank = rank;
825  insert_surf_based_Tag(face, surfmap);
826  }
827  }
828 
829  for (size_t i = 0; i < qsurf_con.size(); i += 4) {
830  if (ptsData[qsurf_con[i ] - 1] > 0 &&
831  ptsData[qsurf_con[i + 1] - 1] > 0 &&
832  ptsData[qsurf_con[i + 2] - 1] > 0 &&
833  ptsData[qsurf_con[i + 3] - 1] > 0) {
834  qface.points.v1 = ref_con[qsurf_con[i ] - 1];
835  qface.points.v2 = ref_con[qsurf_con[i + 1] - 1];
836  qface.points.v3 = ref_con[qsurf_con[i + 2] - 1];
837  qface.points.v4 = ref_con[qsurf_con[i + 3] - 1];
838  qface.rank = rank;
839  insert_surf_based_Tag(qface, qsurfmap);
840  }
841  }
842 }
843 
850 template<class T, class V>
851 struct emi_face {
852  V points;
853  T eidx = -1;
854  T tag = 0;
855  T rank = -1;
856  T mem = 1; // default the face on memebarne, otherwise will be asssign to 2 for gap-junction face
857  bool mark_to_take = false; // if the face is selected on one of the ranks on the surface mesh with unique face, it will be true
858  // we assume that the face will be selected on the smaller rank. If both has the same rank, we choose the first pair.
859  T index_unique = -1; // index of face on emi_surface_unique_face_msh
860  T index_one = -1; // index of face on emi_surface_msh
861 };
862 
863 template<class T>
864 struct emi_unique_face {
865  T index_unique = -1;
866  T index_one = -1;
867  T rank = -1;
868 };
869 
896 template<class T>
897 struct emi_index_rank {
898  T index = -1;
899  T rank = -1;
900 };
913 template<class T, class S> inline
914 void compute_surface_with_tags(meshdata<T,S> & mesh,
915  const SF_nbr numbering,
916  hashmap::unordered_map<T,T> & vertex2ptsdata,
917  hashmap::unordered_set<int> & extra_tags,
918  hashmap::unordered_map<tuple<T>,
919  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
920  hashmap::unordered_map<triple<T>,
921  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_surf,
922  hashmap::unordered_map<quadruple<T>,
923  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_surf)
924 {
925  MPI_Comm comm = SF_COMM;
926  int size, rank;
927  MPI_Comm_size(comm, &size);
928  MPI_Comm_rank(comm, &rank);
929 
930  const T* con = mesh.con.data();
931  const T* nbr = mesh.get_numbering(numbering).data();
932  const SF::vector<mesh_int_t> & rnod = mesh.get_numbering(numbering);
933 
935 
936  for(size_t i=0; i<mesh.con.size(); i++){
937  g2ptsData[rnod[con[i]]] = vertex2ptsdata[rnod[con[i]]];
938  }
939 
940  const vector<T> & ref_eidx = mesh.get_numbering(NBR_ELEM_REF);
941 
942  // Collect potential faces whose three vertices lie on the EMI interfaces.
943  // The ptsData, computed earlier from the input mesh, is used to check
944  // the intersection of each vertex with different tag numbers.
945  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
946  T tag = mesh.tag[eidx];
947  T size_elem = mesh.dsp[eidx+1]-mesh.dsp[eidx];
948 
949  vector<T> dofvec(size_elem); // reference nodes of each element
950  vector<T> ptsDatavec(size_elem); // point data on the nodes of each element
951 
952  for (int n = mesh.dsp[eidx], i = 0; n < mesh.dsp[eidx+1];n++,i++)
953  {
954  dofvec[i] = rnod[con[n]];
955  ptsDatavec[i] = g2ptsData[rnod[con[n]]];
956  }
957  std::vector<T> surf_con ;
958  std::vector<T> qsurf_con;
959  std::vector<T> line_con;
960  switch(mesh.type[eidx]) {
961  // I need to find the right order of lines in 2D
962 
963  // I guess that in the 2D case we would define the interfaces as
964  // lines instead of surfaces.
965  case Tri: {
966  line_con = {1, 2,
967  2, 3,
968  3, 1};
969  surf_con = {};
970  qsurf_con = {};
971  break;
972  }
973 
974  case Quad: {
975  line_con = {1, 2,
976  2, 3,
977  3, 4,
978  4, 1};
979  surf_con = {};
980  qsurf_con = {};
981  break;
982  }
983 
984  case Tetra: {
985  // surfaces are (2,3,1) , (1,4,2) , (2,4,3) , (1,3,4)
986  line_con = {};
987  surf_con = {2,3,1,
988  1,4,2,
989  2,4,3,
990  1,3,4};
991  qsurf_con = {};
992  break;
993  }
994  case Pyramid: {
995  // surfaces are (1,5,2) , (2,5,3) , (3,5,4) , (4,5,1) , (1,2,3,4)
996  line_con = {};
997  surf_con = {1,5,2,
998  2,5,3,
999  3,5,4,
1000  4,5,1};
1001  qsurf_con = {1,2,3,4};
1002  break;
1003  }
1004  case Prism: {
1005  // surfaces are (1,2,3) , (4,5,6) , (1,2,6,4) , (2,3,5,6) , (3,1,4,5)
1006  line_con = {};
1007  surf_con = {1,2,3,
1008  4,5,6};
1009  qsurf_con = {1,2,6,4,
1010  2,3,5,6,
1011  3,1,4,5};
1012  break;
1013  }
1014  case Hexa: {
1015  // surfaces are (1,2,3,4) , (3,2,8,7) , (4,3,7,6) , (1,4,6,5) , (2,1,5,8), (6,7,8,5)
1016  line_con = {};
1017  surf_con = {};
1018  qsurf_con = {1,2,3,4,
1019  3,2,8,7,
1020  4,3,7,6,
1021  1,4,6,5,
1022  2,1,5,8,
1023  6,7,8,5 };
1024  break;
1025  }
1026  default:
1027  fprintf(stderr, "%s error: Unsupported element in surface computation!\n", __func__);
1028  exit(1);
1029  }
1030  insert_surf_emi(rank, dofvec, ptsDatavec, ref_eidx[eidx], mesh.tag[eidx], line_con, surf_con, qsurf_con, line_face, tri_surf, quad_surf);
1031  }
1032 
1033  // Find the counter faces corresponding to the selected potential faces.
1034  assign_counter_face(line_face, mesh.comm);
1035  assign_counter_face(tri_surf, mesh.comm);
1036  assign_counter_face(quad_surf, mesh.comm);
1037 
1038  // Removed potential faces without matching counter faces and
1039  // Removed potential faces where both associated tags belong to the extracellular domain.
1040  for(auto it = line_face.begin(); it != line_face.end(); ) {
1041  if( it->second.second.eidx == -1 ||(extra_tags.find(it->second.first.tag) != extra_tags.end() && extra_tags.find(it->second.second.tag) != extra_tags.end()))
1042  it = line_face.erase(it);
1043  else
1044  ++it;
1045  }
1046 
1047  for(auto it = tri_surf.begin(); it != tri_surf.end(); ) {
1048  if( it->second.second.eidx == -1 ||(extra_tags.find(it->second.first.tag) != extra_tags.end() && extra_tags.find(it->second.second.tag) != extra_tags.end()))
1049  it = tri_surf.erase(it);
1050  else
1051  ++it;
1052  }
1053 
1054  for(auto it = quad_surf.begin(); it != quad_surf.end(); ) {
1055  if( it->second.second.eidx == -1 ||(extra_tags.find(it->second.first.tag) != extra_tags.end() && extra_tags.find(it->second.second.tag) != extra_tags.end()))
1056  it = quad_surf.erase(it);
1057  else
1058  ++it;
1059  }
1060 }
1061 
1075 inline bool should_take_first(int tag1, int tag2,
1076  const hashmap::unordered_set<int>& extra_tags,
1077  const hashmap::unordered_set<int>& intra_tags)
1078 {
1079  bool tag1_is_intra = (intra_tags.find(tag1) != intra_tags.end());
1080  bool tag2_is_intra = (intra_tags.find(tag2) != intra_tags.end());
1081  bool tag1_is_extra = (extra_tags.find(tag1) != extra_tags.end());
1082  bool tag2_is_extra = (extra_tags.find(tag2) != extra_tags.end());
1083 
1084  // Case 1: Both are intra tags (gap junction) - take smaller tag
1085  if(tag1_is_intra && tag2_is_intra) {
1086  return tag1 < tag2;
1087  }
1088 
1089  // Case 2: One is extra, one is intra (membrane) - take intra
1090  if(tag1_is_intra && tag2_is_extra) {
1091  return true; // take tag1 (intra)
1092  }
1093  if(tag2_is_intra && tag1_is_extra) {
1094  return false; // take tag2 (intra)
1095  }
1096 
1097  // Default fallback: take smaller tag
1098  return tag1 < tag2;
1099 }
1100 
1122 inline
1123 void assign_map_between_elem_oneface_and_elem_uniqueFace(int rank,
1124  emi_unique_face <mesh_int_t> first,
1125  emi_unique_face <mesh_int_t> second,
1126  hashmap::unordered_map<mesh_int_t, std::pair<emi_index_rank<mesh_int_t>, emi_index_rank<mesh_int_t>>> & map_elem_uniqueFace_to_elem_oneface,
1127  hashmap::unordered_map<mesh_int_t, emi_index_rank<mesh_int_t>> & map_elem_oneface_to_elem_uniqueFace)
1128 {
1129  // Case 1: Both face and counter face belong to the same rank
1130  if(second.index_one==-1 && second.index_unique==-1){
1131  second.index_one = first.index_one;
1132  second.index_unique = first.index_unique;
1133  }
1134  // Case 2: Cross-rank faces - set up mapping for non-owner ranks
1135  else if(first.rank != second.rank) {
1136  // Find our local index_one (from whichever field matches our rank)
1137  int our_index_one = (first.rank == rank) ? first.index_one :
1138  (second.rank == rank) ? second.index_one : -1;
1139  // Find owner's index_unique (from whichever field has valid index_unique)
1140  int owner_index_unique = (first.index_unique >= 0) ? first.index_unique :
1141  (second.index_unique >= 0) ? second.index_unique : -1;
1142  int owner_rank = (first.index_unique >= 0) ? first.rank :
1143  (second.index_unique >= 0) ? second.rank : -1;
1144 
1145  // Only set mapping for non-owner ranks (owner already set it in STEP 3)
1146  if(our_index_one >= 0 && owner_index_unique >= 0 && owner_rank >= 0 && owner_rank != rank) {
1147  map_elem_oneface_to_elem_uniqueFace[our_index_one].index = owner_index_unique;
1148  map_elem_oneface_to_elem_uniqueFace[our_index_one].rank = owner_rank;
1149  }
1150  }
1151 
1152  // Populate map_elem_uniqueFace_to_elem_oneface (only on owner rank where index_unique >= 0)
1153  // Checked both first and second to find which one represents this rank as owner
1154  // After MPI exchange, the owner info could be in either field
1155  int local_index_unique = -1;
1156  emi_index_rank<mesh_int_t> value1, value2;
1157 
1158  if (first.index_unique >= 0 && first.rank == rank) {
1159  // This rank owns the unique face, info is in 'first'
1160  local_index_unique = first.index_unique;
1161  value1.index = first.index_one;
1162  value1.rank = first.rank;
1163  value2.index = second.index_one;
1164  value2.rank = second.rank;
1165  } else if (second.index_unique >= 0 && second.rank == rank) {
1166  // This rank owns the unique face, info is in 'second'
1167  local_index_unique = second.index_unique;
1168  value1.index = second.index_one;
1169  value1.rank = second.rank;
1170  value2.index = first.index_one;
1171  value2.rank = first.rank;
1172  }
1173 
1174  if (local_index_unique >= 0) {
1175  auto itr = map_elem_uniqueFace_to_elem_oneface.find(local_index_unique);
1176  if (itr == map_elem_uniqueFace_to_elem_oneface.end()) {
1177  std::pair <emi_index_rank<mesh_int_t>,emi_index_rank<mesh_int_t>> value;
1178  value = std::make_pair(value1, value2);
1179  map_elem_uniqueFace_to_elem_oneface.insert({local_index_unique, value});
1180  }
1181  }
1182 }
1183 
1206 template <class T, class Key, class PointsKey>
1207 inline void assign_ownership_rank_on_faces( const Key& key,
1208  std::pair<emi_face<T, PointsKey>, emi_face<T, PointsKey>>& v, // (face, counter-face)
1209  int rank,
1210  mesh_int_t idx_oneface,
1211  mesh_int_t& lsize_unique,
1212  const hashmap::unordered_set<int>& extra_tags,
1213  const hashmap::unordered_set<int>& intra_tags,
1214  hashmap::unordered_map<Key,std::pair<emi_unique_face<mesh_int_t>, emi_unique_face<mesh_int_t>>>& unique_face_to_elements,
1215  hashmap::unordered_map<mesh_int_t, emi_index_rank<mesh_int_t>>& map_elem_oneface_to_elem_uniqueFace)
1216 {
1217  const bool same_rank = (v.first.rank == v.second.rank);
1218 
1219  // Decide owner rank in an order-independent way:
1220  // - If one side is intra, always take the intra side
1221  // - If both intra, take smaller tag (tie -> smaller rank)
1222  // - Otherwise, take smaller tag (tie -> smaller rank)
1223  const bool first_intra = (intra_tags.find(v.first.tag) != intra_tags.end());
1224  const bool second_intra = (intra_tags.find(v.second.tag) != intra_tags.end());
1225 
1226  // Decide which side to take based purely on tags/intra (order-independent)
1227  bool take_first_by_tags = true;
1228  if (first_intra != second_intra) {
1229  take_first_by_tags = first_intra; // take intra side
1230  } else if (v.first.tag != v.second.tag) {
1231  take_first_by_tags = (v.first.tag < v.second.tag);
1232  } else {
1233  // Same tag (should be rare for interfaces); fall back to smaller rank to be deterministic
1234  take_first_by_tags = (v.first.rank <= v.second.rank);
1235  }
1236 
1237  int owner_rank = take_first_by_tags ? v.first.rank : v.second.rank;
1238 
1239  const bool i_am_owner = (rank == owner_rank);
1240 
1241  // If I'm the owner, mark the chosen side and set mapping now.
1242  if (i_am_owner) {
1243  auto& chosen = take_first_by_tags ? v.first : v.second;
1244 
1245  chosen.mark_to_take = true;
1246  chosen.index_unique = lsize_unique;
1247  chosen.index_one = idx_oneface;
1248 
1249  map_elem_oneface_to_elem_uniqueFace[idx_oneface].index = lsize_unique;
1250  map_elem_oneface_to_elem_uniqueFace[idx_oneface].rank = rank;
1251  }
1252 
1253  // Insert per-face record once:
1254  // first = "my side if I'm owner else -1" + rank info
1255  // second = placeholder for the other rank (or same rank if same_rank)
1256  auto it = unique_face_to_elements.find(key);
1257  if (it == unique_face_to_elements.end()) {
1258  emi_unique_face<mesh_int_t> first;
1259  emi_unique_face<mesh_int_t> second;
1260 
1261  if (i_am_owner) {
1262  // Owner knows its local unique index.
1263  first.index_unique = lsize_unique;
1264  first.index_one = idx_oneface;
1265  first.rank = rank;
1266 
1267  // "other side" rank:
1268  second.index_unique = -1;
1269  second.index_one = -1;
1270  second.rank = same_rank ? rank : (rank == v.first.rank ? v.second.rank : v.first.rank);
1271  } else {
1272  // Non-owner does not know unique index yet; MPI step will fill it.
1273  first.index_unique = -1;
1274  first.index_one = idx_oneface;
1275  first.rank = rank;
1276 
1277  second.index_unique = -1;
1278  second.index_one = -1;
1279  second.rank = owner_rank; // this is the rank that will eventually provide index_unique
1280  }
1281 
1282  unique_face_to_elements.insert({ key, std::make_pair(first, second) });
1283  } else {
1284  // Update existing entry to avoid losing owner info when non-owner inserted first.
1285  auto& existing = it->second;
1286 
1287  // Choose a slot to update for this rank.
1288  emi_unique_face<mesh_int_t>* slot = nullptr;
1289  if (existing.first.rank == rank || existing.first.rank < 0) {
1290  slot = &existing.first;
1291  } else if (existing.second.rank == rank || existing.second.rank < 0) {
1292  slot = &existing.second;
1293  } else if (existing.first.index_one == idx_oneface) {
1294  slot = &existing.first;
1295  } else if (existing.second.index_one == idx_oneface) {
1296  slot = &existing.second;
1297  } else if (existing.first.index_unique < 0) {
1298  slot = &existing.first;
1299  } else if (existing.second.index_unique < 0) {
1300  slot = &existing.second;
1301  }
1302 
1303  if (slot) {
1304  if (slot->rank < 0) slot->rank = rank;
1305  if (slot->index_one < 0) slot->index_one = idx_oneface;
1306  if (i_am_owner && slot->index_unique < 0) {
1307  slot->index_unique = lsize_unique;
1308  }
1309  }
1310  }
1311 
1312  // Only the owning rank increments lsize_unique
1313  if (i_am_owner) {
1314  lsize_unique += 1;
1315  }
1316 }
1317 
1318 
1377 template<class T, class S> inline
1378 void extract_face_based_tags(meshdata<T,S> & mesh,
1379  const SF_nbr numbering,
1380  hashmap::unordered_map<T,T> & vertex2ptsdata,
1381  hashmap::unordered_map<tuple<T> ,
1382  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
1383  hashmap::unordered_map<triple<T>,
1384  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
1385  hashmap::unordered_map<quadruple<T>,
1386  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face,
1387  hashmap::unordered_set<int> & extra_tags,
1388  hashmap::unordered_set<int> & intra_tags,
1389  meshdata<T,S> & surfmesh,
1390  meshdata<T,S> & surfmesh_w_counter,
1391  meshdata<T,S> & surfmesh_unique_face,
1392  hashmap::unordered_map<mesh_int_t, std::pair<emi_index_rank<mesh_int_t>, emi_index_rank<mesh_int_t>>> & map_elem_uniqueFace_to_elem_oneface,
1393  hashmap::unordered_map<mesh_int_t, emi_index_rank<mesh_int_t>> & map_elem_oneface_to_elem_uniqueFace)
1394 {
1395  // ============================================================================
1396  // STEP 1: Initialize surface meshes and extract interface faces
1397  // ============================================================================
1398  surfmesh.register_numbering(SF::NBR_REF);
1399  surfmesh_w_counter.register_numbering(SF::NBR_REF);
1400  surfmesh_unique_face.register_numbering(SF::NBR_REF);
1401 
1402  // Extract faces with counter faces on the interfaces to build surface mesh
1403 
1404  compute_surface_with_tags(mesh, numbering, vertex2ptsdata, extra_tags, line_face, tri_face, quad_face);
1405 
1406  int size, rank;
1407  MPI_Comm_size(mesh.comm, &size);
1408  MPI_Comm_rank(mesh.comm, &rank);
1409 
1410  long int g_num_line = line_face.size();
1411  long int g_num_tri = tri_face.size();
1412  long int g_num_quad = quad_face.size();
1413  MPI_Allreduce(MPI_IN_PLACE, &g_num_line, 1, MPI_LONG, MPI_SUM, mesh.comm);
1414  MPI_Allreduce(MPI_IN_PLACE, &g_num_tri, 1, MPI_LONG, MPI_SUM, mesh.comm);
1415  MPI_Allreduce(MPI_IN_PLACE, &g_num_quad, 1, MPI_LONG, MPI_SUM, mesh.comm);
1416 
1417  // Warn if any extra-extra faces remain (should be none after filtering)
1418  {
1419  long int local_extra_extra = 0;
1420  for (const auto& [key, v] : line_face) {
1421  if (extra_tags.find(v.first.tag) != extra_tags.end() &&
1422  extra_tags.find(v.second.tag) != extra_tags.end()) {
1423  local_extra_extra++;
1424  }
1425  }
1426  for (const auto& [key, v] : tri_face) {
1427  if (extra_tags.find(v.first.tag) != extra_tags.end() &&
1428  extra_tags.find(v.second.tag) != extra_tags.end()) {
1429  local_extra_extra++;
1430  }
1431  }
1432  for (const auto& [key, v] : quad_face) {
1433  if (extra_tags.find(v.first.tag) != extra_tags.end() &&
1434  extra_tags.find(v.second.tag) != extra_tags.end()) {
1435  local_extra_extra++;
1436  }
1437  }
1438  long int global_extra_extra = 0;
1439  MPI_Allreduce(&local_extra_extra, &global_extra_extra, 1, MPI_LONG, MPI_SUM, mesh.comm);
1440  if (global_extra_extra > 0 && rank == 0) {
1441  fprintf(stderr, "WARN: extra-extra faces detected before filtering: %ld\n", global_extra_extra);
1442  }
1443  }
1444 
1445  surfmesh.g_numelem = g_num_line + g_num_tri + g_num_quad;
1446  surfmesh.l_numelem = line_face.size() + tri_face.size() + quad_face.size();
1447 
1448  vector<T> cnt(surfmesh.l_numelem, 0);
1449  surfmesh.tag.resize(surfmesh.l_numelem);
1450  surfmesh.type.resize(surfmesh.l_numelem);
1451  surfmesh.con.resize(line_face.size() * 2 + tri_face.size() * 3 + quad_face.size() * 4);
1452 
1453  // ============================================================================
1454  // STEP 2: Sort face keys for deterministic ordering across MPI ranks
1455  // ============================================================================
1456  // To ensure deterministic order, we must sort the keys of the hashmaps before iterating.
1457  // This is crucial for MPI consistency - all ranks must process faces in the same order.
1458  vector<tuple<T>> line_keys;
1459  for(auto const& [key, val] : line_face) line_keys.push_back(key);
1460  std::sort(line_keys.begin(), line_keys.end());
1461 
1462  vector<triple<T>> tri_keys;
1463  for(auto const& [key, val] : tri_face) tri_keys.push_back(key);
1464  std::sort(tri_keys.begin(), tri_keys.end());
1465 
1466  vector<quadruple<T>> quad_keys;
1467  for(auto const& [key, val] : quad_face) quad_keys.push_back(key);
1468  std::sort(quad_keys.begin(), quad_keys.end());
1469 
1470  // Temporary hashmaps to track unique-face assignments before MPI communication
1471  hashmap::unordered_map<SF::tuple<mesh_int_t>, std::pair<emi_unique_face <mesh_int_t>,emi_unique_face<mesh_int_t>>> line_unique_face_to_elements;
1472  hashmap::unordered_map<SF::triple<mesh_int_t>, std::pair<emi_unique_face<mesh_int_t>,emi_unique_face<mesh_int_t>>> tri_unique_face_to_elements;
1473  hashmap::unordered_map<SF::quadruple<mesh_int_t>, std::pair<emi_unique_face<mesh_int_t>,emi_unique_face<mesh_int_t>>> quad_unique_face_to_elements;
1474 
1475  // Counters for local element indices:
1476  // - idx: index into surfmesh (one-side mesh)
1477  // - lsize_*: count of elements in surfmesh_w_counter
1478  // - lsize_unique_*: count of elements in surfmesh_unique_face
1479  mesh_int_t idx = 0, cidx = 0;
1480  mesh_int_t lsize_line = 0;
1481  mesh_int_t lsize_tri = 0;
1482  mesh_int_t lsize_quad = 0;
1483 
1484  mesh_int_t lsize_unique_line = 0;
1485  mesh_int_t lsize_unique_tri = 0;
1486  mesh_int_t lsize_unique_quad = 0;
1487 
1488  // ============================================================================
1489  // STEP 3: Process each face and determine unique-face ownership
1490  // ============================================================================
1491  // For each interface face, we:
1492  // 1. Add it to surfmesh
1493  for(auto const& key : line_keys) {
1494  auto& v = line_face.at(key);
1495  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
1496  if (!local_involved) {
1497  continue;
1498  }
1499  cnt[idx] = 2;
1500 
1501  // Mark gap-junction faces explicitly; membrane is the default.
1502  if ((intra_tags.find(v.first.tag) != intra_tags.end()) && (intra_tags.find(v.second.tag) != intra_tags.end())) {
1503  v.first.mem = 2; // to set the face located on gap-junction
1504  v.second.mem = 2; // to set the face located on gap-junction
1505  }
1506 
1507  // changing the strategy to select the face with the one with the same rank
1508  emi_face<T,tuple<T>> surf_neighbor =
1509  (v.first.rank == rank) ? v.first : v.second;
1510 
1511  surfmesh.type[idx] = Line;
1512  surfmesh.tag[idx] = surf_neighbor.tag;
1513  surfmesh.con[cidx + 0] = surf_neighbor.points.v1;
1514  surfmesh.con[cidx + 1] = surf_neighbor.points.v2;
1515 
1516  // If the ranks are equal, both faces will be added later.
1517  if(v.first.rank == v.second.rank)
1518  lsize_line+=2;
1519  else
1520  lsize_line+=1;
1521 
1522  assign_ownership_rank_on_faces<T, tuple<T>, tuple<T>>( key, v, rank,
1523  idx,
1524  lsize_unique_line,
1525  extra_tags, intra_tags,
1526  line_unique_face_to_elements,
1527  map_elem_oneface_to_elem_uniqueFace);
1528  idx += 1;
1529  cidx += 2;
1530  }
1531 
1532  for(auto const& key : tri_keys) {
1533  auto& v = tri_face.at(key);
1534  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
1535  if (!local_involved) {
1536  continue;
1537  }
1538  cnt[idx] = 3;
1539 
1540  if ((intra_tags.find(v.first.tag) != intra_tags.end()) && (intra_tags.find(v.second.tag) != intra_tags.end())) {
1541  v.first.mem = 2; // to set the face located on gap-junction
1542  v.second.mem = 2; // to set the face located on gap-junction
1543  }
1544 
1545  // changing the strategy to select the face with the one with the same rank
1546  emi_face<T,triple<T>> surf_neighbor =
1547  (v.first.rank == rank) ? v.first : v.second;
1548 
1549  surfmesh.type[idx] = Tri;
1550  surfmesh.tag[idx] = surf_neighbor.tag;
1551  surfmesh.con[cidx + 0] = surf_neighbor.points.v1;
1552  surfmesh.con[cidx + 1] = surf_neighbor.points.v2;
1553  surfmesh.con[cidx + 2] = surf_neighbor.points.v3;
1554 
1555  // If the ranks are equal, both faces will be added later.
1556  if(v.first.rank == v.second.rank)
1557  lsize_tri+=2;
1558  else
1559  lsize_tri+=1;
1560 
1561  assign_ownership_rank_on_faces<T, triple<T>, triple<T>>( key, v, rank,
1562  idx,
1563  lsize_unique_tri,
1564  extra_tags, intra_tags,
1565  tri_unique_face_to_elements,
1566  map_elem_oneface_to_elem_uniqueFace);
1567  idx += 1;
1568  cidx += 3;
1569  }
1570 
1571  for(auto const& key : quad_keys) {
1572  auto& v = quad_face.at(key);
1573  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
1574  if (!local_involved) {
1575  continue;
1576  }
1577  cnt[idx] = 4;
1578 
1579  if ((intra_tags.find(v.first.tag) != intra_tags.end()) && (intra_tags.find(v.second.tag) != intra_tags.end())) {
1580  v.first.mem = 2; // to set the face located on gap-junction
1581  v.second.mem = 2; // to set the face located on gap-junction
1582  }
1583 
1584  // changing the strategy to select the face with the one with the same rank
1585  emi_face<T,quadruple<T>> qsurf_neighbor =
1586  (v.first.rank == rank) ? v.first : v.second;
1587 
1588  surfmesh.type[idx] = Quad;
1589  surfmesh.tag[idx] = qsurf_neighbor.tag;
1590  surfmesh.con[cidx + 0] = qsurf_neighbor.points.v1;
1591  surfmesh.con[cidx + 1] = qsurf_neighbor.points.v2;
1592  surfmesh.con[cidx + 2] = qsurf_neighbor.points.v3;
1593  surfmesh.con[cidx + 3] = qsurf_neighbor.points.v4;
1594 
1595  // If the ranks are equal, both faces will be added later.
1596  if(v.first.rank == v.second.rank)
1597  lsize_quad+=2;
1598  else
1599  lsize_quad+=1;
1600 
1601  assign_ownership_rank_on_faces<T, quadruple<T>, quadruple<T>>( key, v, rank,
1602  idx,
1603  lsize_unique_quad,
1604  extra_tags, intra_tags,
1605  quad_unique_face_to_elements,
1606  map_elem_oneface_to_elem_uniqueFace);
1607  idx += 1;
1608  cidx += 4;
1609  }
1610 
1611  surfmesh.l_numelem = idx;
1612  surfmesh.tag.resize(idx);
1613  surfmesh.type.resize(idx);
1614  cnt.resize(idx);
1615  dsp_from_cnt(cnt, surfmesh.dsp);
1616  surfmesh.con.resize(cidx);
1617 
1618  long int g_num_surf = surfmesh.l_numelem;
1619  MPI_Allreduce(MPI_IN_PLACE, &g_num_surf, 1, MPI_LONG, MPI_SUM, mesh.comm);
1620  surfmesh.g_numelem = g_num_surf;
1621 
1622  if(rank ==0){
1623  std::cout << "surfmesh.g_numelem: " << surfmesh.g_numelem <<std::endl;
1624  std::cout << "surfmesh.l_numelem: " << surfmesh.l_numelem <<std::endl;
1625  }
1626 
1627  {
1628  long int g_num_w_counter_line = lsize_line;
1629  long int g_num_w_counter_tri = lsize_tri;
1630  long int g_num_w_counter_quad = lsize_quad;
1631 
1632  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_counter_line, 1, MPI_LONG, MPI_SUM, SF_COMM);
1633  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_counter_tri, 1, MPI_LONG, MPI_SUM, SF_COMM);
1634  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_counter_quad, 1, MPI_LONG, MPI_SUM, SF_COMM);
1635 
1636  // Assign global and local elements for the surface mesh used in ionic computation.
1637  surfmesh_w_counter.g_numelem = g_num_w_counter_line+g_num_w_counter_tri+g_num_w_counter_quad;
1638  surfmesh_w_counter.l_numelem = lsize_line + lsize_tri + lsize_quad;
1639 
1640  surfmesh_w_counter.tag.resize(surfmesh_w_counter.l_numelem);
1641  surfmesh_w_counter.type.resize(surfmesh_w_counter.l_numelem);
1642  surfmesh_w_counter.con.resize(lsize_line * 2 + lsize_tri * 3 + lsize_quad * 4);
1643 
1644  if(rank ==0){
1645  std::cout << "surfmesh_w_counter.g_numelem: " << surfmesh_w_counter.g_numelem <<std::endl;
1646  std::cout << "surfmesh_w_counter.l_numelem: " << surfmesh_w_counter.l_numelem <<std::endl;
1647  }
1648  }
1649 
1650  {
1651  long int g_num_w_unique_line = lsize_unique_line;
1652  long int g_num_w_unique_tri = lsize_unique_tri;
1653  long int g_num_w_unique_quad = lsize_unique_quad;
1654 
1655  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_unique_line, 1, MPI_LONG, MPI_SUM, SF_COMM);
1656  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_unique_tri, 1, MPI_LONG, MPI_SUM, SF_COMM);
1657  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_unique_quad, 1, MPI_LONG, MPI_SUM, SF_COMM);
1658 
1659  // Assign global and local elements for the surface mesh used in ionic computation.
1660  surfmesh_unique_face.g_numelem = g_num_w_unique_line+g_num_w_unique_tri+g_num_w_unique_quad;
1661  surfmesh_unique_face.l_numelem = lsize_unique_line + lsize_unique_tri + lsize_unique_quad;
1662 
1663  surfmesh_unique_face.tag.resize(surfmesh_unique_face.l_numelem);
1664  surfmesh_unique_face.type.resize(surfmesh_unique_face.l_numelem);
1665  surfmesh_unique_face.con.resize(lsize_unique_line * 2 + lsize_unique_tri * 3 + lsize_unique_quad * 4);
1666 
1667  if(rank ==0){
1668  std::cout << "surfmesh_unique_face.g_numelem: " << surfmesh_unique_face.g_numelem <<std::endl;
1669  std::cout << "surfmesh_unique_face.l_numelem: " << surfmesh_unique_face.l_numelem <<std::endl;
1670  }
1671  }
1672 
1673  // ============================================================================
1674  // STEP 4: MPI communication to share unique-face ownership between ranks
1675  // ============================================================================
1676  // When faces span multiple ranks, only one rank "owns" the unique-face.
1677  // assign_unique_first_face communicates the owner's index_unique to non-owner ranks
1678  // so they can set up their map_elem_oneface_to_elem_uniqueFace mappings correctly.
1679 
1680  // DEBUG: Before exchange
1681  #ifdef EMI_DEBUG_MESH
1682  {
1683  fprintf(stderr, "RANK %d BEFORE exchange: line_hash=%zu, tri_hash=%zu, quad_hash=%zu (expected unique total=%ld)\n",
1684  rank, line_unique_face_to_elements.size(), tri_unique_face_to_elements.size(),
1685  quad_unique_face_to_elements.size(), lsize_unique_line + lsize_unique_tri + lsize_unique_quad);
1686  fflush(stderr);
1687  }
1688  #endif
1689 
1690  if(size>1)
1691  {
1692  // Find the counter faces corresponding to the selected potential faces.
1693  assign_unique_first_face(line_unique_face_to_elements, mesh.comm);
1694  assign_unique_first_face(tri_unique_face_to_elements, mesh.comm);
1695  assign_unique_first_face(quad_unique_face_to_elements, mesh.comm);
1696  }
1697 
1698  // DEBUG: After exchange
1699  #ifdef EMI_DEBUG_MESH
1700  {
1701  int valid_line = 0, valid_tri = 0, valid_quad = 0;
1702  for (const auto& [key, val] : line_unique_face_to_elements) {
1703  if (val.first.index_unique >= 0 || val.second.index_unique >= 0) valid_line++;
1704  }
1705  for (const auto& [key, val] : tri_unique_face_to_elements) {
1706  if (val.first.index_unique >= 0 || val.second.index_unique >= 0) valid_tri++;
1707  }
1708  for (const auto& [key, val] : quad_unique_face_to_elements) {
1709  if (val.first.index_unique >= 0 || val.second.index_unique >= 0) valid_quad++;
1710  }
1711  fprintf(stderr, "RANK %d AFTER exchange: line_hash=%zu, tri_hash=%zu, quad_hash=%zu\n",
1712  rank, line_unique_face_to_elements.size(), tri_unique_face_to_elements.size(),
1713  quad_unique_face_to_elements.size());
1714  fprintf(stderr, "RANK %d VALID entries (index_unique>=0): line=%d, tri=%d, quad=%d\n",
1715  rank, valid_line, valid_tri, valid_quad);
1716  fflush(stderr);
1717  }
1718  #endif
1719 
1720  // ============================================================================
1721  // STEP 5: Finalize mappings after MPI communication
1722  // ============================================================================
1723  for(auto it = line_unique_face_to_elements.begin(); it != line_unique_face_to_elements.end(); ++it) {
1724 
1725  auto& first = it->second.first;
1726  auto& second = it->second.second;
1727 
1728  assign_map_between_elem_oneface_and_elem_uniqueFace(rank, first, second, map_elem_uniqueFace_to_elem_oneface, map_elem_oneface_to_elem_uniqueFace);
1729  }
1730 
1731  for(auto it = tri_unique_face_to_elements.begin(); it != tri_unique_face_to_elements.end(); ++it) {
1732  auto& first = it->second.first;
1733  auto& second = it->second.second;
1734 
1735  assign_map_between_elem_oneface_and_elem_uniqueFace(rank, first, second, map_elem_uniqueFace_to_elem_oneface, map_elem_oneface_to_elem_uniqueFace);
1736  }
1737 
1738  for(auto it = quad_unique_face_to_elements.begin(); it != quad_unique_face_to_elements.end(); ++it) {
1739  auto& first = it->second.first;
1740  auto& second = it->second.second;
1741 
1742  assign_map_between_elem_oneface_and_elem_uniqueFace(rank, first, second, map_elem_uniqueFace_to_elem_oneface, map_elem_oneface_to_elem_uniqueFace);
1743  }
1744 }
1745 
1755 template<class T> inline
1756 void complete_map_vertex_to_dof_with_counter_face(hashmap::unordered_map<tuple<T>,
1757  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
1758  hashmap::unordered_map<triple<T>,
1759  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
1760  hashmap::unordered_map<quadruple<T>,
1761  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face,
1762  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
1763  mesh_int_t> & map_vertex_tag_to_dof)
1764 {
1765  MPI_Comm comm = SF_COMM;
1766  int size, rank;
1767  MPI_Comm_size(comm, &size);
1768  MPI_Comm_rank(comm, &rank);
1769 
1770  // look for counter part line_inteface to add to the map key<vertex,tag> -> value<dof>
1771  for(const auto & v : line_face) {
1772 
1773  // if the ranks are equal, then it has already added to the map
1774  if(v.second.first.rank == v.second.second.rank){
1775  continue;
1776  }
1777  // otherwise we look for the counter part with different rank
1778  emi_face<T,tuple<T>> surf_neighbor =
1779  (v.second.first.rank != rank) ? v.second.first : v.second.second;
1780 
1781  {
1782  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1783  Index_tag_old = std::make_pair(surf_neighbor.points.v1,surf_neighbor.tag);
1784  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1785  }
1786 
1787  {
1788  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1789  Index_tag_old = std::make_pair(surf_neighbor.points.v2,surf_neighbor.tag);
1790  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1791  }
1792  }
1793 
1794  // look for counter part line_inteface to add to the map key<vertex,tag> -> value<dof>
1795  for(const auto & v : tri_face) {
1796 
1797  // if the ranks are equal, then it has already added to the map
1798  if(v.second.first.rank == v.second.second.rank){
1799  continue;
1800  }
1801  // otherwise we look for the counter part with different rank
1802  emi_face<T,triple<T>> surf_neighbor =
1803  (v.second.first.rank != rank) ? v.second.first : v.second.second;
1804 
1805  {
1806  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1807  Index_tag_old = std::make_pair(surf_neighbor.points.v1,surf_neighbor.tag);
1808  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1809  }
1810 
1811  {
1812  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1813  Index_tag_old = std::make_pair(surf_neighbor.points.v2,surf_neighbor.tag);
1814  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1815  }
1816 
1817  {
1818  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1819  Index_tag_old = std::make_pair(surf_neighbor.points.v3,surf_neighbor.tag);
1820  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1821  }
1822  }
1823 
1824  // look for counter part line_inteface to add to the map key<vertex,tag> -> value<dof>
1825  for(const auto & v : quad_face) {
1826 
1827  // if the ranks are equal, then it has already added to the map
1828  if(v.second.first.rank == v.second.second.rank){
1829  continue;
1830  }
1831  // otherwise we look for the counter part with different rank
1832  emi_face<T,quadruple<T>> qsurf_neighbor =
1833  (v.second.first.rank != rank) ? v.second.first : v.second.second;
1834 
1835  {
1836  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1837  Index_tag_old = std::make_pair(qsurf_neighbor.points.v1,qsurf_neighbor.tag);
1838  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1839  }
1840 
1841  {
1842  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1843  Index_tag_old = std::make_pair(qsurf_neighbor.points.v2,qsurf_neighbor.tag);
1844  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1845  }
1846 
1847  {
1848  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1849  Index_tag_old = std::make_pair(qsurf_neighbor.points.v3,qsurf_neighbor.tag);
1850  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1851  }
1852 
1853  {
1854  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1855  Index_tag_old = std::make_pair(qsurf_neighbor.points.v4,qsurf_neighbor.tag);
1856  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1857  }
1858  }
1859 
1860  // assign the counter faces which we haven't assigned it yet the new dof
1861  assign_dof_on_counter_face(map_vertex_tag_to_dof,comm);
1862 }
1863 
1893 template<class T, class S> inline
1894 void compute_map_vertex_to_dof(meshdata<T,S> & mesh,
1895  const SF_nbr numbering,
1896  hashmap::unordered_map<T,T> & vertex2ptsdata,
1897  hashmap::unordered_set<int> & extra_tags,
1898  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
1899  mesh_int_t> & map_vertex_tag_to_dof)
1900 {
1901  MPI_Comm comm = mesh.comm;
1902  int size, rank;
1903  MPI_Comm_size(comm, &size);
1904  MPI_Comm_rank(comm, &rank);
1905 
1906  hashmap::unordered_map<mesh_int_t, intersection_data> map_vertex_to_tags_data_ranks;
1907 
1908  const T* con = mesh.con.data();
1909  const T* nbr = mesh.get_numbering(numbering).data();
1910  const SF::vector<T> & rnod = mesh.get_numbering(numbering);
1911 
1912  hashmap::unordered_map<T,T> g2ptsData;
1913 
1914  for(size_t i=0; i<mesh.con.size(); i++)
1915  {
1916  g2ptsData[rnod[con[i]]] = vertex2ptsdata[rnod[con[i]]];
1917  }
1918 
1919  //-----------------------------------------------------------------
1920  // Step 1: Local Information Gathering
1921  //-----------------------------------------------------------------
1922  // - Identify Simple Cases: For any vertex that is not on an interface (i.e., it's completely inside the extracellular space or completely
1923  // inside a unique myocyte), the mapping is simple. The new DOF index is the same as the original vertex index. The function populates the
1924  // output map, map_vertex_tag_to_dof, with these direct mappings, e.g., (vertex {1}, tag_3) -> 1.
1925  // - Collect Intersection Data: This is the most important part of this step. For every vertex, it builds a list of all the regions (tags) and
1926  // MPI ranks that share it, based on the elements the current process knows about. This information is stored in a temporary map called
1927  // map_vertex_to_tags_data_ranks.
1928  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++)
1929  {
1930  T tag = mesh.tag[eidx];
1931  auto pos_extra = extra_tags.find(tag);
1932  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
1933  {
1934  T gIndex = rnod[con[n]];
1935  T data_on_gIndex = g2ptsData[rnod[con[n]]];
1936 
1937  // This handles inner DoFs, either within the extracellular domain or vertices belonging to the intracellular domain.
1938  // and marked as a visited dofs and saved key<gIndex,tag)> -> val<gIndex>
1939  if((pos_extra != extra_tags.end()) || data_on_gIndex==0){
1940  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = std::make_pair(gIndex,tag);
1941  if (map_vertex_tag_to_dof.find(Index_tag_old) == map_vertex_tag_to_dof.end() )
1942  {
1943  map_vertex_tag_to_dof.insert({Index_tag_old,gIndex});
1944  }
1945  }
1946 
1947  // When gIndex is located on an interface — membrane, gap junction, or their intersection.
1948  auto it = map_vertex_to_tags_data_ranks.find(gIndex);
1949  if (it != map_vertex_to_tags_data_ranks.end() )
1950  {
1951  intersection_data& value = it->second;
1952  bool check_tag_rank = true;
1953  for(int i=0; i<MAX_INTERSECTIONS; ++i) {
1954  if(value.tags[i] == tag && value.ranks[i] == rank) {
1955  check_tag_rank = false;
1956  break;
1957  }
1958  }
1959 
1960  // Insert the tag and corresponding rank if they haven’t been inserted yet.
1961  if(check_tag_rank){
1962  for(int i=0; i<MAX_INTERSECTIONS; ++i) {
1963  if(value.tags[i] == -1 && value.ranks[i] == -1) {
1964  value.tags[i] = tag;
1965  value.data[i] = data_on_gIndex;
1966  value.ranks[i] = rank;
1967  break;
1968  }
1969  }
1970  }
1971  }
1972  else
1973  {
1974  // insert gIndex to the map_vertex_to_tags_data_ranks
1975  intersection_data value;
1976  value.tags[0] = tag;
1977  value.data[0] = data_on_gIndex;
1978  value.ranks[0] = rank;
1979  map_vertex_to_tags_data_ranks.insert({gIndex,value});
1980  }
1981  }
1982  }
1983  // until here: map_vertex_to_tags_data_ranks knows about the interfaces for the elements it owns,
1984  // but it doesn't know if a vertex it owns is also part of an interface on another process.
1985 
1986  //-----------------------------------------------------------------
1987  // Step 2: Globalize Information with MPI
1988  //-----------------------------------------------------------------
1989  // - assign_counter_vertices_tuple: This helper function is called to manage the communication. It takes the partial
1990  // map_vertex_to_tags_data_ranks from each process and uses MPI (specifically MPI_Exchange, which is likely a wrapper around routines like
1991  // MPI_Alltoallv) to send and receive data.
1992  // - Build Global Map: After the communication, the partial maps are merged. The result is that the map_vertex_to_tags_data_ranks on every
1993  // process now contains the complete sharing information for every vertex in the entire mesh. Each process now knows the full list of tags
1994  // and ranks associated with any given vertex.
1995  assign_counter_vertices_tuple(map_vertex_to_tags_data_ranks, comm);
1996 
1997  //-----------------------------------------------------------------
1998  // Step 3: Identify and Count DOFs to Be Created
1999  //-----------------------------------------------------------------
2000  // Now that every process has the same global information, they can independently and deterministically decide which new DOFs to create.
2001  // - Analyze Intersections: The code iterates through the now-global map_vertex_to_tags_data_ranks. For each vertex, it inspects the list of
2002  // tags that share it to classify the type of interface.
2003  // - Apply Rules: Based on the interface type, it decides which side needs a new, unique DOF index.
2004  // * Membrane (Myocyte-Extracellular): The extracellular DOF keeps the original vertex index. The myocyte DOF is marked as needing a new, unique index.
2005  // * Gap Junction (Myocyte-Myocyte): To separate the two myocytes, one of them needs a new index. A deterministic rule (e.g., the myocyte
2006  // with the higher tag number gets the new index) is applied to ensure all processes make the same decision.
2007  // * Complex Junctions: For even more complex intersections (e.g., where multiple myocytes and the extracellular space meet), similar
2008  // rules are applied to create the necessary new DOFs.
2009  // - Count Local Contribution: Each process counts how many new DOFs (shift) it is responsible for creating based on these rules.
2010  int shift = 0;
2012  for(const auto & key : map_vertex_to_tags_data_ranks)
2013  {
2014  T gIndex = key.first;
2015  const intersection_data& value = key.second;
2016  T data_on_gIndex = g2ptsData[gIndex];
2017 
2018  // find the first -1 in value
2019  int first_index = 0;
2020  for (int i = 0; i < MAX_INTERSECTIONS; ++i) {
2021  if(value.tags[i] == -1 && value.ranks[i] == -1) {
2022  first_index = i;
2023  break;
2024  }else{
2025  // check if all the intersection_data has the same g2ptsData
2026  if(value.data[i]!=data_on_gIndex)
2027  {
2028  // Error for unhandled cases.
2029  fprintf(stderr,
2030  "WARN: Due to an inconsistency in ptsData computation at gIndex %d with tag %d, "
2031  "the program is stopped while preparing the data structure to introduce DOFs.\n",
2032  gIndex, value.tags[i]);
2033  fprintf(stderr, "\n");
2034  EXIT(1);
2035  }
2036  }
2037  if (i == MAX_INTERSECTIONS -1) first_index = MAX_INTERSECTIONS;
2038  }
2039 
2040  if(data_on_gIndex==1){ // membrane
2041  // check the mumber of tag numbers,
2042  int count_intra_tags = 0;
2043  int count_extra_tags = 0;
2044  for (int i = 0; i < first_index; ++i) {
2045  if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])==extra_tags.end()){
2046  count_intra_tags += 1;
2047  }else if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])!=extra_tags.end()){
2048  count_extra_tags += 1;
2049  }
2050  }
2051 
2052  // Check if the collected data at gIndex is indeed on the membrane.
2053  // A valid membrane vertex must intersect exactly one intracellular tag and one or more extracellular tags.
2054  if(count_intra_tags>1 || count_extra_tags == 0)
2055  {
2056  // Error for unhandled cases.
2057  fprintf(stderr,
2058  "WARN: Due to an issue in defining ptsData for gIndex %d, , on the membrane"
2059  "the program is stopped while preparing the data structure to introduce DOFs.\n",
2060  gIndex);
2061  fprintf(stderr, "\n");
2062  EXIT(1);
2063  }
2064 
2065  int index_intra = -1;
2066  // Find the intracellular tag number. Since data_on_gIndex == 1, there will be only one tag corresponding to the intracellular domain.
2067  for (int i = 0; i < first_index; ++i) {
2068  if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])==extra_tags.end()){
2069  index_intra = i;
2070  break;
2071  }
2072  }
2073  if(index_intra == -1) {
2074  fprintf(stderr, "WARN: unhandled case in membrane while decoupling: gIndex %d", gIndex);
2075  fprintf(stderr, "\n");
2076  EXIT(1);
2077  }
2078  if(index_intra != -1) {
2079  T tag_myocyte = value.tags[index_intra];
2080  T rank_myocyte = value.ranks[index_intra];
2081 
2082  if(rank_myocyte == rank) { // mark intracellular tag with gIndex
2083  std::pair <mesh_int_t,mesh_int_t> Index_tag_next = std::make_pair(gIndex,tag_myocyte);
2084  if (map_mark_new_dofs.find(Index_tag_next) == map_mark_new_dofs.end()) {
2085  map_mark_new_dofs.insert({Index_tag_next,true}); //mark a new index
2086  shift++;
2087  }
2088  }
2089  }
2090  }
2091  else if(data_on_gIndex==2){ // gap junction
2092  // count the number of extra tags and intra tags
2093  int count_intra_tags = 0;
2094  int count_extra_tags = 0;
2095  for (int i = 0; i < first_index; ++i) {
2096  if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])==extra_tags.end()){
2097  count_intra_tags += 1;
2098  }else if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])!=extra_tags.end()){
2099  count_extra_tags += 1;
2100  }
2101  }
2102  if(count_intra_tags!=2 || count_extra_tags != 0){
2103  // Error for unhandled cases.
2104  fprintf(stderr,
2105  "WARN: Due to an issue in defining ptsData for gIndex %d, on the gapjunction"
2106  "the program is stopped while preparing the data structure to introduce DOFs.\n",
2107  gIndex);
2108  fprintf(stderr, "\n");
2109  EXIT(1);
2110  }
2111 
2112  T tag1 = -1;
2113  T tag2 = -1;
2114  T rank1 = -1;
2115  T rank2 = -1;
2116  bool tag1_checked = false;
2117  bool tag2_checked = false;
2118  // find the tag1 and tag2 on the gap junction
2119  for (int i = 0; i < first_index; ++i)
2120  {
2121  auto pos_extra = extra_tags.find(value.tags[i]);
2122  if(pos_extra == extra_tags.end()) {
2123  if(tag1==-1){
2124  tag1 = value.tags[i];
2125  rank1 = value.ranks[i];
2126  tag1_checked = true;
2127  }else if(tag1!=-1){
2128  tag2 = value.tags[i];
2129  rank2 = value.ranks[i];
2130  tag2_checked = true;
2131  break;
2132  }
2133  }
2134  }
2135 
2136  if(!tag1_checked && !tag2_checked || (rank1==rank2 && rank1!=rank) ) {
2137  fprintf(stderr, "WARN: unhandled case in gap junction while decoupling: gIndex %d with tags %d:%d ", gIndex, tag1, tag2);
2138  fprintf(stderr, "\n");
2139  EXIT(1);
2140  }
2141 
2142  if(rank1==rank2 && (rank1 == rank)){ // If both tag numbers belong to the same rank, keep the smaller tag’s vertex index and reserve a new index for the larger tag.
2143  T smallest_tag = tag1 < tag2? tag1:tag2;
2144  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = std::make_pair(gIndex,smallest_tag);
2145  if (map_vertex_tag_to_dof.find(Index_tag_old) == map_vertex_tag_to_dof.end() ) {
2146  map_vertex_tag_to_dof.insert({Index_tag_old,gIndex});
2147  }
2148 
2149  T bigger_tag = tag1 < tag2? tag2:tag1;
2150  std::pair <mesh_int_t,mesh_int_t> Index_tag_new = std::make_pair(gIndex,bigger_tag);
2151  if (map_mark_new_dofs.find(Index_tag_new) == map_mark_new_dofs.end() ) {
2152  map_mark_new_dofs.insert({Index_tag_new,true});
2153  shift++;
2154  }
2155  }
2156  else if(rank1!=rank2 && (rank1 == rank)) { // If one belong to the current rank, keep the smaller tag’s vertex index else reserve a new index for the larger tag.
2157  if(tag1<tag2) {
2158  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = std::make_pair(gIndex,tag1);
2159  if (map_vertex_tag_to_dof.find(Index_tag_old) == map_vertex_tag_to_dof.end() ) {
2160  map_vertex_tag_to_dof.insert({Index_tag_old,gIndex});
2161  }
2162  } else {
2163  std::pair <mesh_int_t,mesh_int_t> Index_tag_new = std::make_pair(gIndex,tag1);
2164  if (map_mark_new_dofs.find(Index_tag_new) == map_mark_new_dofs.end() ) {
2165  map_mark_new_dofs.insert({Index_tag_new,true});
2166  shift++;
2167  }
2168  }
2169  }
2170  else if(rank1!=rank2 && (rank2 == rank)) { // If one belong to the current rank, keep the smaller tag’s vertex index else reserve a new index for the larger tag.
2171  if(tag2<tag1) {
2172  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = std::make_pair(gIndex,tag2);
2173  if (map_vertex_tag_to_dof.find(Index_tag_old) == map_vertex_tag_to_dof.end() ) {
2174  map_vertex_tag_to_dof.insert({Index_tag_old,gIndex});
2175  }
2176  } else {
2177  std::pair <mesh_int_t,mesh_int_t> Index_tag_new = std::make_pair(gIndex,tag2);
2178  if (map_mark_new_dofs.find(Index_tag_new) == map_mark_new_dofs.end() ) {
2179  map_mark_new_dofs.insert({Index_tag_new,true});
2180  shift++;
2181  }
2182  }
2183  }
2184  }
2185  else if(data_on_gIndex==3)
2186  {
2187  // count the number of extra tags and intra tags
2188  int count_intra_tags = 0;
2189  int count_extra_tags = 0;
2190  for (int i = 0; i < first_index; ++i) {
2191  if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])==extra_tags.end()){
2192  count_intra_tags += 1;
2193  }else if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])!=extra_tags.end()){
2194  count_extra_tags += 1;
2195  }
2196  }
2197  if(count_intra_tags!=2 || count_extra_tags == 0){
2198  // Error for unhandled cases.
2199  fprintf(stderr,
2200  "WARN: Due to an issue in defining ptsData for gIndex %d, on intersection between membrane and gapjunction"
2201  "the program is stopped while preparing the data structure to introduce DOFs.\n",
2202  gIndex);
2203  fprintf(stderr, "\n");
2204  EXIT(1);
2205  }
2206 
2207  // All extracellular tags have already been added in the previous loop that iterated over the mesh.
2208  // on the interface between gap juntion and membarne
2209  for (int i = 0; i < first_index; ++i) {
2210  // find the first two tag number which belong to itracellular tags, add for indexing
2211  auto pos_extra = extra_tags.find(value.tags[i]);
2212  if(value.tags[i]!=-1 && value.ranks[i] == rank && pos_extra == extra_tags.end()) { // myocytes, for sure we should have two myocytes
2213  std::pair <mesh_int_t,mesh_int_t> Index_tag_new = std::make_pair(gIndex,value.tags[i]);
2214  if (map_mark_new_dofs.find(Index_tag_new) == map_mark_new_dofs.end() )
2215  {
2216  map_mark_new_dofs.insert({Index_tag_new,true}); //mark a new index
2217  shift++;
2218  }
2219  }
2220  }
2221  }
2222  }
2223 
2224  vector<size_t> dsp_dof(size);
2225  MPI_Allgather(&shift, sizeof(size_t), MPI_BYTE, dsp_dof.data(), sizeof(size_t), MPI_BYTE, comm);
2226 
2227  //-----------------------------------------------------------------
2228  // Step 4: Assign Unique Global Indices in Parallel
2229  //-----------------------------------------------------------------
2230  // The final step is to assign the new indices without any conflicts between processes.
2231  // - Calculate Global Offset (`start`): An MPI_Allgather is used to share the shift counts among all processes. Each process then calculates
2232  // its own unique start index for the block of new DOFs it will create. This is done by taking the total number of vertices in the original
2233  // mesh and adding the shift counts from all processes with a smaller rank. This is a parallel prefix sum, and it guarantees that the ranges
2234  // of new indices created by different processes will not overlap.
2235  // - Assign New Indices: Finally, the function loops through the vertices it marked for new DOF creation and assigns them a unique global
2236  // index (newIndex = start + count). This final mapping, (original_vertex_index, tag) -> new_unique_dof_index, is inserted into the output
2237  // map map_vertex_tag_to_dof.
2238  T start;
2239  if(rank==0){
2240  start = mesh.g_numpts;
2241  }
2242  else
2243  {
2244  start = mesh.g_numpts;
2245  for (int r = 0; r < rank; ++r)
2246  {
2247  start+= dsp_dof[r];
2248  }
2249  }
2250 
2251  int count = 0;
2252  auto lexicographic_comp_pair = [](const std::pair<mesh_int_t,mesh_int_t> & a, const std::pair<mesh_int_t,mesh_int_t> & b)
2253  {
2254  if (a.second == b.second) return a.first < b.first;
2255  return a.second < b.second;
2256  };
2257 
2258  map_mark_new_dofs.sort(lexicographic_comp_pair);
2259 
2260  for(const auto & entry : map_mark_new_dofs)
2261  {
2262  T newIndex = start+count;
2263  map_vertex_tag_to_dof.insert({entry.first,newIndex});
2264  count++;
2265  }
2266 }
2267 
2274 template<class K, class intersection_indices> inline
2275 void assign_counter_dofs(hashmap::unordered_map<K, intersection_indices>& map, const MPI_Comm comm)
2276 {
2277  size_t dsize = map.size();
2278  vector<K> key_vec(dsize);
2279  vector<intersection_indices> value_vec(dsize);
2280  intersection_indices indices;
2281  size_t idx = 0;
2282  for (const auto& v : map) {
2283  key_vec[idx] = v.first;
2284  value_vec[idx] = v.second;
2285  idx++;
2286  }
2287 
2288  vector<int> perm, dest;
2289  emi_select_merge_destinations(key_vec, dsize, sizeof(K) + sizeof(intersection_indices), dest, comm,
2290  "assign_counter_dofs",
2291  [](const K& key) { return hashmap::hash_ops<K>::hash(key); });
2292 
2293  commgraph<size_t> grph;
2294  grph.configure(dest, comm);
2295  size_t nrecv = sum(grph.rcnt);
2296 
2297  interval(perm, 0, dsize);
2298  binary_sort_copy(dest, perm);
2299 
2300  // fill send buffer and communicate
2301  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
2302  vector<intersection_indices> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
2303  for (size_t i = 0; i < dsize; i++) {
2304  sbuff_key[i] = key_vec[perm[i]];
2305  sbuff_value[i] = value_vec[perm[i]];
2306  }
2307  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
2308  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
2309 
2311  for (size_t i = 0; i < nrecv; i++)
2312  {
2313  auto it = rmap.find(rbuff_key[i]);
2314  if (it != rmap.end())
2315  {
2316  intersection_indices& map_indices = it->second;
2317  intersection_indices& r_indices = rbuff_value[i];
2318 
2319  for (int r_indices : r_indices.indices)
2320  {
2321  if (r_indices == -1) continue;
2322  bool found = false;
2323  for (int m_index : map_indices.indices)
2324  {
2325  if (m_index == r_indices) {
2326  found = true;
2327  break;
2328  }
2329  }
2330  if (!found)
2331  {
2332  for (size_t k = 0; k < MAX_INTERSECTIONS; ++k) {
2333  if (map_indices.indices[k] == -1) {
2334  map_indices.indices[k] = r_indices;
2335  break;
2336  }
2337  }
2338  }
2339  }
2340  }
2341  else
2342  {
2343  rmap.insert({ rbuff_key[i], rbuff_value[i] });
2344  }
2345  }
2346 
2347  for (size_t i = 0; i < nrecv; i++) {
2348  auto it = rmap.find(rbuff_key[i]);
2349  if (it != rmap.end()) rbuff_value[i] = it->second;
2350  }
2351 
2352  // sending the assigned data back to original rank
2353  grph.transpose();
2354  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
2355 
2356  for (size_t i = 0; i < dsize; i++) {
2357  auto it = map.find(sbuff_key[i]);
2358  if (it != map.end()) it->second = ibuff_value[i];
2359  }
2360 }
2371 template<class T, class S> inline
2372 void update_emi_mesh_with_dofs(meshdata<T,S> & mesh,
2373  const SF_nbr numbering,
2374  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
2375  mesh_int_t> & map_vertex_tag_to_dof,
2377 {
2378  // map between old index to set of new indices
2380  MPI_Comm comm = mesh.comm;
2381  mesh.globalize(numbering);
2382  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
2383  T tag = mesh.tag[eidx];
2384  T size_elem = mesh.dsp[eidx+1]-mesh.dsp[eidx];
2385 
2386  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
2387  {
2388  T gIndex = mesh.con[n];
2389 
2390  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2391  Index_tag_old = std::make_pair(gIndex,tag);
2392 
2393  auto it_new = map_vertex_tag_to_dof.find(Index_tag_old);
2394  if (it_new != map_vertex_tag_to_dof.end() )
2395  {
2396  mesh_int_t dof_new = (*it_new).second;
2397  dof2vertex.insert({dof_new, gIndex});
2398  mesh.con[n] = dof_new;
2399 
2400  auto it = vertex2dof.find(gIndex);
2401  if (it != vertex2dof.end())
2402  {
2403  intersection_indices& indices = it->second;
2404  bool found = false;
2405  for(T t : indices.indices) {
2406  if(t == dof_new) {
2407  found = true;
2408  break;
2409  }
2410  }
2411  if(!found) {
2412  // the default size is 3, since the maximum intersection shoudl be on the intersection between the gap junction and memeranre, and since the 2 indiecs for myocyesand one indices for extracellular.
2413  for(size_t i=0; i<MAX_INTERSECTIONS; ++i) {
2414  if(indices.indices[i] == -1) {
2415  indices.indices[i] = dof_new;
2416  break;
2417  }
2418  }
2419  }
2420  }else{
2421  intersection_indices indices;
2422  indices.indices[0] = dof_new;
2423  vertex2dof.insert({gIndex,indices});
2424  }
2425  }
2426  }
2427  }
2428 
2429  assign_counter_dofs(vertex2dof, comm);
2430 
2431  // Now update dof2vertex using vertex2dof
2432  for (const auto& [old_idx, indices] : vertex2dof) {
2433  for (int new_idx : indices.indices) {
2434  if (new_idx != -1) {
2435  // Avoid overwriting unless necessary
2436  if (dof2vertex.find(new_idx) == dof2vertex.end()) {
2437  dof2vertex.insert({new_idx, old_idx});
2438  }
2439  }
2440  }
2441  }
2442 
2443  mesh.localize(numbering);
2444 }
2445 
2457 template<class T, class S> inline
2458 void update_map_indices_to_petsc(meshdata<T,S> & mesh,
2459  const SF_nbr numbering_ref,const SF_nbr numbering_petsc,
2460  hashmap::unordered_set<int> extra_tags,
2461  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
2462  std::pair<mesh_int_t,mesh_int_t>> & map_vertex_tag_to_dof_petsc,
2464  SF::vector<mesh_int_t> & elemTag_emi_mesh)
2465 {
2466  const SF::vector<SF_int> & ref_nbr = mesh.get_numbering(numbering_ref);
2467  const SF::vector<SF_int> & petsc_nbr = mesh.get_numbering(numbering_petsc);
2468 
2469  elemTag_emi_mesh.resize(mesh.l_numelem);
2470  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
2471  T tag = mesh.tag[eidx];
2472  elemTag_emi_mesh[eidx] = 1; // assigned 1 for for extracellular located on extracellular side
2473  if(extra_tags.find(tag) == extra_tags.end())
2474  elemTag_emi_mesh[eidx] = 2; // assigned 2 for for intracellular located on myocyte
2475  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
2476  {
2477  T l_Indx = mesh.con[n];
2478  T petsc_Idx = petsc_nbr[l_Indx];
2479  T n_Indx = ref_nbr[l_Indx];
2480  T o_Indx = dof2vertex[n_Indx];
2481 
2482  std::pair <mesh_int_t,mesh_int_t> oIndx_tag;
2483  oIndx_tag = std::make_pair(o_Indx,tag);
2484 
2485  auto it_new = map_vertex_tag_to_dof_petsc.find(oIndx_tag);
2486  if (it_new != map_vertex_tag_to_dof_petsc.end() )
2487  {
2488  std::pair <mesh_int_t,mesh_int_t> nIndx_petsc;
2489  nIndx_petsc = std::make_pair(n_Indx,petsc_Idx);
2490 
2491  (*it_new).second = nIndx_petsc;
2492  }
2493  }
2494  }
2495 }
2496 
2509 template<class T, class S>
2510 inline void update_faces_on_surface_mesh_after_decoupling_with_dofs(meshdata<T, S> & mesh,
2511  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
2512  mesh_int_t> & map_vertex_tag_to_dof,
2513  hashmap::unordered_map<tuple<T>,
2514  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
2515  hashmap::unordered_map<triple<T>,
2516  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
2517  hashmap::unordered_map<quadruple<T>,
2518  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face)
2519 {
2520  MPI_Comm comm = SF_COMM;
2521  int size, rank;
2522  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
2523 
2524  SF::vector<mesh_int_t> idxbuff(mesh.con);
2525  binary_sort(idxbuff); unique_resize(idxbuff);
2526 
2529  for(size_t i=0; i<idxbuff.size(); i++){
2530  g2l[idxbuff[i]] = i;
2531  l2g[i] = idxbuff[i];
2532  }
2533 
2534  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
2535  T tag = mesh.tag[eidx];
2536  T size_elem = mesh.dsp[eidx+1]-mesh.dsp[eidx];
2537  std::vector<int> elem_nodes;
2538  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
2539  {
2540  T gIndex = mesh.con[n];
2541  elem_nodes.push_back(gIndex);
2542  }
2543  std::sort(elem_nodes.begin(),elem_nodes.end());
2544  if(elem_nodes.size()==2){
2545  tuple<T> key;
2546  key.v1 = elem_nodes[0];
2547  key.v2 = elem_nodes[1];
2548  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>> value = line_face[key];
2549 
2550  line_face.erase(key);
2551  std::vector<int> new_nodes(2);
2552 
2553  auto tag_key = (value.first.rank == rank) ? value.first.tag : value.second.tag;
2554  {
2555  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2556  Index_tag_old = std::make_pair(elem_nodes[0],tag_key);
2557  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2558  new_nodes[0] = Index_new;
2559  }
2560  {
2561  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2562  Index_tag_old = std::make_pair(elem_nodes[1],tag_key);
2563  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2564  new_nodes[1] = Index_new;
2565  }
2566  std::sort(new_nodes.begin(),new_nodes.end());
2567  // update key
2568  tuple<T> new_key;
2569  new_key.v1 = new_nodes[0];
2570  new_key.v2 = new_nodes[1];
2571 
2572  // first
2573  {
2574  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v1, value.first.tag)];
2575  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v2, value.first.tag)];
2576  value.first.points.v1 = Index_new1;
2577  value.first.points.v2 = Index_new2;
2578  }
2579  // second
2580  {
2581  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v1, value.second.tag)];
2582  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v2, value.second.tag)];
2583  value.second.points.v1 = Index_new1;
2584  value.second.points.v2 = Index_new2;
2585  }
2586  line_face.insert({new_key,value});
2587  }
2588  if(elem_nodes.size()==3){
2589  triple<T> key;
2590  key.v1 = elem_nodes[0];
2591  key.v2 = elem_nodes[1];
2592  key.v3 = elem_nodes[2];
2593  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>> value = tri_face[key];
2594  tri_face.erase(key);
2595  std::vector<int> new_nodes(3);
2596 
2597  auto tag_key = (value.first.rank == rank) ? value.first.tag : value.second.tag;
2598  {
2599  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2600  Index_tag_old = std::make_pair(elem_nodes[0],tag_key);
2601  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2602  new_nodes[0] = Index_new;
2603  }
2604  {
2605  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2606  Index_tag_old = std::make_pair(elem_nodes[1],tag_key);
2607  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2608  new_nodes[1] = Index_new;
2609  }
2610  {
2611  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2612  Index_tag_old = std::make_pair(elem_nodes[2],tag_key);
2613  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2614  new_nodes[2] = Index_new;
2615  }
2616  std::sort(new_nodes.begin(),new_nodes.end());
2617  // update key
2618  triple<T> new_key;
2619  new_key.v1 = new_nodes[0];
2620  new_key.v2 = new_nodes[1];
2621  new_key.v3 = new_nodes[2];
2622 
2623  // first
2624  {
2625  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v1, value.first.tag)];
2626  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v2, value.first.tag)];
2627  mesh_int_t Index_new3 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v3, value.first.tag)];
2628  value.first.points.v1 = Index_new1;
2629  value.first.points.v2 = Index_new2;
2630  value.first.points.v3 = Index_new3;
2631  }
2632  // second
2633  {
2634  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v1, value.second.tag)];
2635  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v2, value.second.tag)];
2636  mesh_int_t Index_new3 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v3, value.second.tag)];
2637  value.second.points.v1 = Index_new1;
2638  value.second.points.v2 = Index_new2;
2639  value.second.points.v3 = Index_new3;
2640  }
2641  tri_face.insert({new_key,value});
2642  }
2643  if(elem_nodes.size()==4){
2644  quadruple<T> key;
2645  key.v1 = elem_nodes[0];
2646  key.v2 = elem_nodes[1];
2647  key.v3 = elem_nodes[2];
2648  key.v4 = elem_nodes[3];
2649  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>> value = quad_face[key];
2650  quad_face.erase(key);
2651  std::vector<int> new_nodes(4);
2652 
2653  auto tag_key = (value.first.rank == rank) ? value.first.tag : value.second.tag;
2654  {
2655  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2656  Index_tag_old = std::make_pair(elem_nodes[0],tag_key);
2657  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2658  new_nodes[0] = Index_new;
2659  }
2660  {
2661  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2662  Index_tag_old = std::make_pair(elem_nodes[1],tag_key);
2663  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2664  new_nodes[1] = Index_new;
2665  }
2666  {
2667  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2668  Index_tag_old = std::make_pair(elem_nodes[2],tag_key);
2669  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2670  new_nodes[2] = Index_new;
2671  }
2672  {
2673  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2674  Index_tag_old = std::make_pair(elem_nodes[3],tag_key);
2675  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2676  new_nodes[3] = Index_new;
2677  }
2678  std::sort(new_nodes.begin(),new_nodes.end());
2679  // update key
2680  quadruple<T> new_key;
2681  new_key.v1 = new_nodes[0];
2682  new_key.v2 = new_nodes[1];
2683  new_key.v3 = new_nodes[2];
2684  new_key.v4 = new_nodes[3];
2685 
2686  // first
2687  {
2688  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v1, value.first.tag)];
2689  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v2, value.first.tag)];
2690  mesh_int_t Index_new3 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v3, value.first.tag)];
2691  mesh_int_t Index_new4 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v4, value.first.tag)];
2692  value.first.points.v1 = Index_new1;
2693  value.first.points.v2 = Index_new2;
2694  value.first.points.v3 = Index_new3;
2695  value.first.points.v4 = Index_new4;
2696  }
2697  // second
2698  {
2699  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v1, value.second.tag)];
2700  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v2, value.second.tag)];
2701  mesh_int_t Index_new3 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v3, value.second.tag)];
2702  mesh_int_t Index_new4 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v4, value.second.tag)];
2703  value.second.points.v1 = Index_new1;
2704  value.second.points.v2 = Index_new2;
2705  value.second.points.v3 = Index_new3;
2706  value.second.points.v4 = Index_new4;
2707  }
2708  quad_face.insert({new_key,value});
2709  }
2710  }
2711 }
2712 
2724 template<class T, class S> inline
2725 void compute_surface_mesh_with_counter_face(meshdata<T, S> & mesh, const SF_nbr numbering,
2726  hashmap::unordered_map<tuple<T>,
2727  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
2728  hashmap::unordered_map<triple<T>,
2729  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
2730  hashmap::unordered_map<quadruple<T>,
2731  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face)
2732 {
2733  mesh.register_numbering(numbering);
2734 
2735  MPI_Comm comm = SF_COMM;
2736  int size, rank;
2737  MPI_Comm_size(comm, &size);
2738  MPI_Comm_rank(comm, &rank);
2739 
2740  // To ensure deterministic order, we must sort the keys of the hashmaps before iterating.
2741  vector<tuple<T>> line_keys;
2742  for(auto const& [key, val] : line_face) line_keys.push_back(key);
2743  std::sort(line_keys.begin(), line_keys.end());
2744 
2745  vector<triple<T>> tri_keys;
2746  for(auto const& [key, val] : tri_face) tri_keys.push_back(key);
2747  std::sort(tri_keys.begin(), tri_keys.end());
2748 
2749  vector<quadruple<T>> quad_keys;
2750  for(auto const& [key, val] : quad_face) quad_keys.push_back(key);
2751  std::sort(quad_keys.begin(), quad_keys.end());
2752 
2753  vector<T> cnt(mesh.l_numelem);
2754  size_t idx = 0, cidx = 0;
2755 
2756  // Add faces to the surface mesh used in ionic computation.
2757  // For each {line,tri,quad}_face, add faces only if one of the pair matches the current rank,
2758  // or if both faces have the same rank, add both.
2759  for(const auto & key : line_keys) {
2760  const auto & v = line_face.at(key);
2761  bool both_faces = (v.first.rank == v.second.rank) ? true: false;
2762 
2763  cnt[idx] = 2;
2764 
2765  emi_face<T,tuple<T>> surf_first;
2766  emi_face<T,tuple<T>> surf_second;
2767  if(v.first.rank == rank)
2768  {
2769  surf_first = v.first;
2770  surf_second = v.second;
2771  }
2772  else
2773  {
2774  surf_first = v.second;
2775  surf_second = v.first;
2776  }
2777 
2778  mesh.type[idx] = Line;
2779  mesh.tag[idx] = surf_first.tag;
2780  mesh.con[cidx + 0] = surf_first.points.v1;
2781  mesh.con[cidx + 1] = surf_first.points.v2;
2782  idx += 1;
2783  cidx += 2;
2784 
2785  if(both_faces){
2786 
2787  cnt[idx] = 2;
2788  mesh.type[idx] = Line;
2789  mesh.tag[idx] = surf_second.tag;
2790  mesh.con[cidx + 0] = surf_second.points.v1;
2791  mesh.con[cidx + 1] = surf_second.points.v2;
2792 
2793  idx += 1;
2794  cidx += 2;
2795  }
2796  }
2797 
2798  for(const auto & key : tri_keys) {
2799  const auto & v = tri_face.at(key);
2800  bool both_faces = (v.first.rank == v.second.rank) ? true: false;
2801 
2802  cnt[idx] = 3;
2803 
2804  emi_face<T,triple<T>> surf_first;
2805  emi_face<T,triple<T>> surf_second;
2806  if(v.first.rank == rank)
2807  {
2808  surf_first = v.first;
2809  surf_second = v.second;
2810  }
2811  else
2812  {
2813  surf_first = v.second;
2814  surf_second = v.first;
2815  }
2816 
2817  mesh.type[idx] = Tri;
2818  mesh.tag[idx] = surf_first.tag;
2819  mesh.con[cidx + 0] = surf_first.points.v1;
2820  mesh.con[cidx + 1] = surf_first.points.v2;
2821  mesh.con[cidx + 2] = surf_first.points.v3;
2822 
2823  idx += 1;
2824  cidx += 3;
2825 
2826  if(both_faces){
2827 
2828  cnt[idx] = 3;
2829 
2830  mesh.type[idx] = Tri;
2831  mesh.tag[idx] = surf_second.tag;
2832  mesh.con[cidx + 0] = surf_second.points.v1;
2833  mesh.con[cidx + 1] = surf_second.points.v2;
2834  mesh.con[cidx + 2] = surf_second.points.v3;
2835 
2836  idx += 1;
2837  cidx += 3;
2838  }
2839  }
2840 
2841  for(const auto & key : quad_keys) {
2842  const auto & v = quad_face.at(key);
2843  bool both_faces = (v.first.rank == v.second.rank) ? true: false;
2844  cnt[idx] = 4;
2845 
2846  emi_face<T,quadruple<T>> surf_first;
2847  emi_face<T,quadruple<T>> surf_second;
2848  if(v.first.rank == rank)
2849  {
2850  surf_first = v.first;
2851  surf_second = v.second;
2852  }
2853  else
2854  {
2855  surf_first = v.second;
2856  surf_second = v.first;
2857  }
2858 
2859  mesh.type[idx] = Quad;
2860  mesh.tag[idx] = surf_first.tag;
2861  mesh.con[cidx + 0] = surf_first.points.v1;
2862  mesh.con[cidx + 1] = surf_first.points.v2;
2863  mesh.con[cidx + 2] = surf_first.points.v3;
2864  mesh.con[cidx + 3] = surf_first.points.v4;
2865 
2866  idx += 1;
2867  cidx += 4;
2868 
2869  if(both_faces){
2870 
2871  cnt[idx] = 4;
2872 
2873  mesh.type[idx] = Quad;
2874  mesh.tag[idx] = surf_second.tag;
2875  mesh.con[cidx + 0] = surf_second.points.v1;
2876  mesh.con[cidx + 1] = surf_second.points.v2;
2877  mesh.con[cidx + 2] = surf_second.points.v3;
2878  mesh.con[cidx + 3] = surf_second.points.v4;
2879 
2880  idx += 1;
2881  cidx += 4;
2882  }
2883  }
2884  dsp_from_cnt(cnt, mesh.dsp);
2885 }
2886 
2917 template<class T, class S> inline
2918 void compute_surface_mesh_with_unique_face(meshdata<T, S> & mesh, const SF_nbr numbering,
2919  hashmap::unordered_map<tuple<T>,
2920  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
2921  hashmap::unordered_map<triple<T>,
2922  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
2923  hashmap::unordered_map<quadruple<T>,
2924  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face,
2925  hashmap::unordered_map<mesh_int_t, std::pair<mesh_int_t, mesh_int_t>> & map_elem_uniqueFace_to_tags)
2926 {
2927  mesh.register_numbering(numbering);
2928 
2929  MPI_Comm comm = SF_COMM;
2930  int size, rank;
2931  MPI_Comm_size(comm, &size);
2932  MPI_Comm_rank(comm, &rank);
2933 
2934  // To ensure deterministic order, we must sort the keys of the hashmaps before iterating.
2935  vector<tuple<T>> line_keys;
2936  for(auto const& [key, val] : line_face) line_keys.push_back(key);
2937  std::sort(line_keys.begin(), line_keys.end());
2938 
2939  vector<triple<T>> tri_keys;
2940  for(auto const& [key, val] : tri_face) tri_keys.push_back(key);
2941  std::sort(tri_keys.begin(), tri_keys.end());
2942 
2943  vector<quadruple<T>> quad_keys;
2944  for(auto const& [key, val] : quad_face) quad_keys.push_back(key);
2945  std::sort(quad_keys.begin(), quad_keys.end());
2946 
2947  vector<T> cnt(mesh.l_numelem);
2948  size_t idx = 0, cidx = 0;
2949 
2950  // Add faces to the surface mesh used in ionic computation.
2951  // For each {line,tri,quad}_face, add faces only if one of the pair matches the current rank,
2952  // or if both faces have the same rank, add both.
2953  for(const auto & key : line_keys) {
2954  const auto & v = line_face.at(key);
2955 
2956  bool to_take_face = false;
2957 
2958  emi_face<T,tuple<T>> surf_take;
2959  std::pair <mesh_int_t,mesh_int_t> tag_pairs_take;
2960  {
2961  if(v.first.rank == rank && v.first.mark_to_take == true)
2962  {
2963  surf_take = v.first;
2964  to_take_face = true;
2965  tag_pairs_take = std::make_pair(v.first.tag,v.second.tag);
2966  }
2967  else if(v.second.rank == rank && v.second.mark_to_take == true)
2968  {
2969  surf_take = v.second;
2970  to_take_face = true;
2971  tag_pairs_take = std::make_pair(v.second.tag,v.first.tag);
2972  }
2973 
2974  if(to_take_face)
2975  {
2976  map_elem_uniqueFace_to_tags.insert({idx,tag_pairs_take});
2977  cnt[idx] = 2;
2978  mesh.type[idx] = Line;
2979  mesh.tag[idx] = surf_take.tag;
2980  mesh.con[cidx + 0] = surf_take.points.v1;
2981  mesh.con[cidx + 1] = surf_take.points.v2;
2982  idx += 1;
2983  cidx += 2;
2984  }
2985  }
2986  }
2987 
2988  for(const auto & key : tri_keys) {
2989  const auto & v = tri_face.at(key);
2990 
2991  bool to_take_face = false;
2992 
2993  emi_face<T,triple<T>> surf_take;
2994  std::pair <mesh_int_t,mesh_int_t> tag_pairs_take;
2995  {
2996  if(v.first.rank == rank && v.first.mark_to_take == true)
2997  {
2998  surf_take = v.first;
2999  to_take_face = true;
3000  tag_pairs_take = std::make_pair(v.first.tag,v.second.tag);
3001  }
3002  else if(v.second.rank == rank && v.second.mark_to_take == true)
3003  {
3004  surf_take = v.second;
3005  to_take_face = true;
3006  tag_pairs_take = std::make_pair(v.second.tag,v.first.tag);
3007  }
3008 
3009  if(to_take_face)
3010  {
3011  map_elem_uniqueFace_to_tags.insert({idx,tag_pairs_take});
3012  cnt[idx] = 3;
3013  mesh.type[idx] = Tri;
3014  mesh.tag[idx] = surf_take.tag;
3015  mesh.con[cidx + 0] = surf_take.points.v1;
3016  mesh.con[cidx + 1] = surf_take.points.v2;
3017  mesh.con[cidx + 2] = surf_take.points.v3;
3018 
3019  idx += 1;
3020  cidx += 3;
3021  }
3022  }
3023  }
3024 
3025  for(const auto & key : quad_keys) {
3026  const auto & v = quad_face.at(key);
3027  cnt[idx] = 4;
3028 
3029  bool to_take_face = false;
3030 
3031  emi_face<T,quadruple<T>> surf_take;
3032  std::pair <mesh_int_t,mesh_int_t> tag_pairs_take;
3033  {
3034  if(v.first.rank == rank && v.first.mark_to_take == true)
3035  {
3036  surf_take = v.first;
3037  to_take_face = true;
3038  tag_pairs_take = std::make_pair(v.first.tag,v.second.tag);
3039  }
3040  else if(v.second.rank == rank && v.second.mark_to_take == true)
3041  {
3042  surf_take = v.second;
3043  to_take_face = true;
3044  tag_pairs_take = std::make_pair(v.second.tag,v.first.tag);
3045  }
3046 
3047  if(to_take_face)
3048  {
3049  map_elem_uniqueFace_to_tags.insert({idx,tag_pairs_take});
3050  cnt[idx] = 4;
3051  mesh.type[idx] = Quad;
3052  mesh.tag[idx] = surf_take.tag;
3053  mesh.con[cidx + 0] = surf_take.points.v1;
3054  mesh.con[cidx + 1] = surf_take.points.v2;
3055  mesh.con[cidx + 2] = surf_take.points.v3;
3056  mesh.con[cidx + 3] = surf_take.points.v4;
3057 
3058  idx += 1;
3059  cidx += 4;
3060  }
3061  }
3062  }
3063  dsp_from_cnt(cnt, mesh.dsp);
3064 }
3065 
3085 template<class T> inline
3086 void create_reverse_elem_mapping_between_surface_meshes(
3087  hashmap::unordered_map<tuple<T>, std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>>& line_face,
3088  hashmap::unordered_map<triple<T>, std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>>& tri_face,
3089  hashmap::unordered_map<quadruple<T>, std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>>& quad_face,
3090  SF::vector<T>& vec_one_to_both_face,
3091  MPI_Comm comm)
3092 {
3093  int rank;
3094  MPI_Comm_rank(comm, &rank);
3095 
3096  // To ensure deterministic order, we must sort the keys of the hashmaps before iterating.
3097  vector<tuple<T>> line_keys;
3098  for(auto const& [key, val] : line_face) line_keys.push_back(key);
3099  std::sort(line_keys.begin(), line_keys.end());
3100 
3101  vector<triple<T>> tri_keys;
3102  for(auto const& [key, val] : tri_face) tri_keys.push_back(key);
3103  std::sort(tri_keys.begin(), tri_keys.end());
3104 
3105  vector<quadruple<T>> quad_keys;
3106  for(auto const& [key, val] : quad_face) quad_keys.push_back(key);
3107  std::sort(quad_keys.begin(), quad_keys.end());
3108 
3109  size_t surf_elem_idx = 0;
3110  size_t w_counter_elem_idx = 0;
3111 
3112  size_t lsize_line = 0;
3113  size_t lsize_tri = 0;
3114  size_t lsize_quad = 0;
3115 
3116  // Count how many both-face elements are local on this rank.
3117  // If both faces are on the same rank, the both-face mesh has two entries.
3118  for(const auto& key : line_keys) {
3119  const auto& v = line_face.at(key);
3120  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3121  if (!local_involved) {
3122  continue;
3123  }
3124  if(v.first.rank == v.second.rank) lsize_line+=2;
3125  else lsize_line+=1;
3126  }
3127  for(const auto& key : tri_keys) {
3128  const auto& v = tri_face.at(key);
3129  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3130  if (!local_involved) {
3131  continue;
3132  }
3133  if(v.first.rank == v.second.rank) lsize_tri+=2;
3134  else lsize_tri+=1;
3135  }
3136  for(const auto& key : quad_keys) {
3137  const auto& v = quad_face.at(key);
3138  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3139  if (!local_involved) {
3140  continue;
3141  }
3142  if(v.first.rank == v.second.rank) lsize_quad+=2;
3143  else lsize_quad+=1;
3144  }
3145 
3146  vec_one_to_both_face.resize(lsize_line + lsize_tri + lsize_quad);
3147 
3148  // Fill mapping in deterministic order:
3149  // each both-face element gets the index of its corresponding one-face element.
3150  for (const auto& key : line_keys) {
3151  const auto& v = line_face.at(key);
3152  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3153  if (!local_involved) {
3154  continue;
3155  }
3156  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3157  if (v.first.rank == v.second.rank) {
3158  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3159  }
3160  surf_elem_idx++;
3161  }
3162 
3163  for (const auto& key : tri_keys) {
3164  const auto& v = tri_face.at(key);
3165  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3166  if (!local_involved) {
3167  continue;
3168  }
3169  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3170  if (v.first.rank == v.second.rank) {
3171  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3172  }
3173  surf_elem_idx++;
3174  }
3175 
3176  for (const auto& key : quad_keys) {
3177  const auto& v = quad_face.at(key);
3178  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3179  if (!local_involved) {
3180  continue;
3181  }
3182  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3183  if (v.first.rank == v.second.rank) {
3184  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3185  }
3186  surf_elem_idx++;
3187  }
3188 }
3189 
3198 template<class T> inline
3199 void added_counter_faces_to_map(hashmap::unordered_map<tuple<T>,
3200  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
3201  hashmap::unordered_map<triple<T>,
3202  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
3203  hashmap::unordered_map<quadruple<T>,
3204  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face)
3205 {
3206  // Buffer inserts not to mutate an unordered_map while iterating over it and later add to unordered_map
3207  std::vector<std::pair<tuple<T>, std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>>> line_inserts;
3208  for(const auto & v : line_face) {
3209  SF::emi_face<mesh_int_t,SF::tuple<mesh_int_t>> surf_first = v.second.first;
3210  SF::emi_face<mesh_int_t,SF::tuple<mesh_int_t>> surf_second = v.second.second;
3211 
3212  std::vector<mesh_int_t> elem_nodes_first; elem_nodes_first.resize(2);
3213  std::vector<mesh_int_t> elem_nodes_second; elem_nodes_second.resize(2);
3214 
3215  SF::tuple<mesh_int_t> key_first;
3216  elem_nodes_first[0] = surf_first.points.v1;
3217  elem_nodes_first[1] = surf_first.points.v2;
3218  std::sort(elem_nodes_first.begin(),elem_nodes_first.end());
3219  key_first.v1 = elem_nodes_first[0];
3220  key_first.v2 = elem_nodes_first[1];
3221 
3222  SF::tuple<mesh_int_t> key_second;
3223  elem_nodes_second[0] = surf_second.points.v1;
3224  elem_nodes_second[1] = surf_second.points.v2;
3225  std::sort(elem_nodes_second.begin(),elem_nodes_second.end());
3226  key_second.v1 = elem_nodes_second[0];
3227  key_second.v2 = elem_nodes_second[1];
3228 
3229  const bool same_key_first = (v.first.v1 == key_first.v1 && v.first.v2 == key_first.v2);
3230  if (!same_key_first && line_face.find(key_first) == line_face.end()) {
3231  line_inserts.push_back({key_first, v.second});
3232  }
3233  const bool same_key_second = (v.first.v1 == key_second.v1 && v.first.v2 == key_second.v2);
3234  if (!same_key_second && line_face.find(key_second) == line_face.end()) {
3235  line_inserts.push_back({key_second, v.second});
3236  }
3237  }
3238  for (const auto & entry : line_inserts) line_face.insert(entry);
3239 
3240  std::vector<std::pair<triple<T>, std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>>> tri_inserts;
3241  for(const auto & v : tri_face) {
3242  SF::emi_face<mesh_int_t,SF::triple<mesh_int_t>> surf_first = v.second.first;
3243  SF::emi_face<mesh_int_t,SF::triple<mesh_int_t>> surf_second = v.second.second;
3244 
3245  std::vector<mesh_int_t> elem_nodes_first; elem_nodes_first.resize(3);
3246  std::vector<mesh_int_t> elem_nodes_second; elem_nodes_second.resize(3);
3247 
3248  SF::triple<mesh_int_t> key_first;
3249  elem_nodes_first[0] = surf_first.points.v1;
3250  elem_nodes_first[1] = surf_first.points.v2;
3251  elem_nodes_first[2] = surf_first.points.v3;
3252  std::sort(elem_nodes_first.begin(),elem_nodes_first.end());
3253  key_first.v1 = elem_nodes_first[0];
3254  key_first.v2 = elem_nodes_first[1];
3255  key_first.v3 = elem_nodes_first[2];
3256 
3257  SF::triple<mesh_int_t> key_second;
3258  elem_nodes_second[0] = surf_second.points.v1;
3259  elem_nodes_second[1] = surf_second.points.v2;
3260  elem_nodes_second[2] = surf_second.points.v3;
3261  std::sort(elem_nodes_second.begin(),elem_nodes_second.end());
3262  key_second.v1 = elem_nodes_second[0];
3263  key_second.v2 = elem_nodes_second[1];
3264  key_second.v3 = elem_nodes_second[2];
3265 
3266  const bool same_key_first = (v.first.v1 == key_first.v1 &&
3267  v.first.v2 == key_first.v2 &&
3268  v.first.v3 == key_first.v3);
3269  if (!same_key_first && tri_face.find(key_first) == tri_face.end()) {
3270  tri_inserts.push_back({key_first, v.second});
3271  }
3272  const bool same_key_second = (v.first.v1 == key_second.v1 &&
3273  v.first.v2 == key_second.v2 &&
3274  v.first.v3 == key_second.v3);
3275  if (!same_key_second && tri_face.find(key_second) == tri_face.end()) {
3276  tri_inserts.push_back({key_second, v.second});
3277  }
3278  }
3279  for (const auto & entry : tri_inserts) tri_face.insert(entry);
3280 
3281  std::vector<std::pair<quadruple<T>, std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>>> quad_inserts;
3282  for(const auto & v : quad_face) {
3283  SF::emi_face<mesh_int_t,SF::quadruple<mesh_int_t>> surf_first = v.second.first;
3284  SF::emi_face<mesh_int_t,SF::quadruple<mesh_int_t>> surf_second = v.second.second;
3285 
3286  std::vector<mesh_int_t> elem_nodes_first; elem_nodes_first.resize(4);
3287  std::vector<mesh_int_t> elem_nodes_second; elem_nodes_second.resize(4);
3288 
3289  SF::quadruple<mesh_int_t> key_first;
3290  elem_nodes_first[0] = surf_first.points.v1;
3291  elem_nodes_first[1] = surf_first.points.v2;
3292  elem_nodes_first[2] = surf_first.points.v3;
3293  elem_nodes_first[3] = surf_first.points.v4;
3294  std::sort(elem_nodes_first.begin(),elem_nodes_first.end());
3295  key_first.v1 = elem_nodes_first[0];
3296  key_first.v2 = elem_nodes_first[1];
3297  key_first.v3 = elem_nodes_first[2];
3298  key_first.v4 = elem_nodes_first[3];
3299 
3300  SF::quadruple<mesh_int_t> key_second;
3301  elem_nodes_second[0] = surf_second.points.v1;
3302  elem_nodes_second[1] = surf_second.points.v2;
3303  elem_nodes_second[2] = surf_second.points.v3;
3304  elem_nodes_second[3] = surf_second.points.v4;
3305  std::sort(elem_nodes_second.begin(),elem_nodes_second.end());
3306  key_second.v1 = elem_nodes_second[0];
3307  key_second.v2 = elem_nodes_second[1];
3308  key_second.v3 = elem_nodes_second[2];
3309  key_second.v4 = elem_nodes_second[3];
3310 
3311  const bool same_key_first = (v.first.v1 == key_first.v1 &&
3312  v.first.v2 == key_first.v2 &&
3313  v.first.v3 == key_first.v3 &&
3314  v.first.v4 == key_first.v4);
3315  if (!same_key_first && quad_face.find(key_first) == quad_face.end()) {
3316  quad_inserts.push_back({key_first, v.second});
3317  }
3318  const bool same_key_second = (v.first.v1 == key_second.v1 &&
3319  v.first.v2 == key_second.v2 &&
3320  v.first.v3 == key_second.v3 &&
3321  v.first.v4 == key_second.v4);
3322  if (!same_key_second && quad_face.find(key_second) == quad_face.end()) {
3323  quad_inserts.push_back({key_second, v.second});
3324  }
3325  }
3326  for (const auto & entry : quad_inserts) quad_face.insert(entry);
3327 }
3328 
3335 template<class K> inline
3336 void assign_counter_tags(hashmap::unordered_map<K, intersection_tags>& map, const MPI_Comm comm)
3337 {
3338  size_t dsize = map.size();
3339  vector<K> key_vec(dsize);
3340  vector<intersection_tags> value_vec(dsize);
3341 
3342  size_t idx = 0;
3343  for (const auto& v : map) {
3344  key_vec[idx] = v.first;
3345  value_vec[idx] = v.second;
3346  idx++;
3347  }
3348 
3349  vector<int> perm, dest;
3350  emi_select_merge_destinations(key_vec, dsize, sizeof(K) + sizeof(intersection_tags), dest, comm,
3351  "assign_counter_tags",
3352  [](const K& key) { return hashmap::hash_ops<K>::hash(key); });
3353 
3354  commgraph<size_t> grph;
3355  grph.configure(dest, comm);
3356  size_t nrecv = sum(grph.rcnt);
3357 
3358  interval(perm, 0, dsize);
3359  binary_sort_copy(dest, perm);
3360 
3361  // fill send buffer and communicate
3362  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
3363  vector<intersection_tags> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
3364  for (size_t i = 0; i < dsize; i++) {
3365  sbuff_key[i] = key_vec[perm[i]];
3366  sbuff_value[i] = value_vec[perm[i]];
3367  }
3368  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
3369  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
3370 
3372  for (size_t i = 0; i < nrecv; i++)
3373  {
3374  auto it = rmap.find(rbuff_key[i]);
3375  if (it != rmap.end())
3376  {
3377  intersection_tags& map_tags = it->second;
3378  intersection_tags& r_tags = rbuff_value[i];
3379 
3380  for (int r_tag : r_tags.tags)
3381  {
3382  if (r_tag == -1) continue;
3383  bool found = false;
3384  for (int m_tag : map_tags.tags)
3385  {
3386  if (m_tag == r_tag) {
3387  found = true;
3388  break;
3389  }
3390  }
3391  if (!found)
3392  {
3393  for (size_t k = 0; k < MAX_INTERSECTIONS; ++k) {
3394  if (map_tags.tags[k] == -1) {
3395  map_tags.tags[k] = r_tag;
3396  break;
3397  }
3398  }
3399  }
3400  }
3401  }
3402  else
3403  {
3404  rmap.insert({ rbuff_key[i], rbuff_value[i] });
3405  }
3406  }
3407 
3408  for (size_t i = 0; i < nrecv; i++) {
3409  auto it = rmap.find(rbuff_key[i]);
3410  if (it != rmap.end()) rbuff_value[i] = it->second;
3411  }
3412 
3413  // sending the assigned data back to original rank
3414  grph.transpose();
3415  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
3416 
3417  for (size_t i = 0; i < dsize; i++) {
3418  auto it = map.find(sbuff_key[i]);
3419  if (it != map.end()) it->second = ibuff_value[i];
3420  }
3421 }
3422 
3444 template<class T, class S> inline
3445 void compute_ptsdata_from_original_mesh(meshdata<T,S> & mesh,
3446  const SF_nbr numbering,
3447  hashmap::unordered_map<T,T> & vertex2ptsdata,
3448  hashmap::unordered_set<int> extra_tags,
3449  hashmap::unordered_set<int> intra_tags)
3450 {
3451  MPI_Comm comm = mesh.comm;
3452  int size, rank;
3453  MPI_Comm_size(comm, &size);
3454  MPI_Comm_rank(comm, &rank);
3455 
3456  const T* con = mesh.con.data();
3457  const SF::vector<T> & rnod = mesh.get_numbering(numbering);
3458 
3459  // Step 1: Gather all unique region tags associated with each vertex.
3460  // This map is local to each process and will store a list of tags for each vertex
3461  // based on the elements owned by the current process.
3463  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++)
3464  {
3465  T tag = mesh.tag[eidx];
3466  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
3467  {
3468  T gIndex = rnod[con[n]];
3469 
3470  auto it_new = map_index_to_tags.find(gIndex);
3471  if (it_new != map_index_to_tags.end())
3472  {
3473  // If vertex is already in the map, add the new tag if it's not already present.
3474  intersection_tags& tags = it_new->second;
3475  bool found = false;
3476  for(T t : tags.tags) {
3477  if(t == tag) {
3478  found = true;
3479  break;
3480  }
3481  }
3482  if(!found) {
3483  // This is an intracellular tag. Find an empty slot and add it.
3484  for(int i=0; i <MAX_INTERSECTIONS; ++i) {
3485  if (tags.tags[i] == -1) {
3486  tags.tags[i] = tag;
3487  break;
3488  }
3489  }
3490  }
3491  }
3492  else{
3493  // If vertex is not in the map, add it with the current tag.
3494  intersection_tags tags;
3495  tags.tags[0] = tag;
3496  map_index_to_tags.insert({gIndex,tags});
3497  }
3498  }
3499  }
3500 
3501  // Step 2: Globalize the tag information.
3502  // After this call, `map_index_to_tags` on every process will contain the
3503  // complete list of unique tags for every vertex in the entire mesh.
3504  assign_counter_tags(map_index_to_tags, comm);
3505 
3506  // Step 3: Classify each vertex based on its global list of tags.
3507  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++)
3508  {
3509  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
3510  {
3511  T gIndex = rnod[con[n]];
3512  intersection_tags& vec_tags = map_index_to_tags[gIndex];
3513  int count_intra_tags = 0;
3514  int count_extra_tags = 0;
3515 
3516  // Count the number of unique intracellular and extracellular tags for the vertex.
3517  int count_tags = 0;
3518  for(T t : vec_tags.tags) {
3519  if(intra_tags.count(t)) count_intra_tags++;
3520  if(extra_tags.count(t)) count_extra_tags++;
3521  }
3522 
3523  // Apply classification rules based on the tag counts.
3524  if(count_extra_tags>=1 && count_intra_tags==0){
3525  vertex2ptsdata[gIndex] = 0; // Interior to an extracellular region
3526  }else if(count_extra_tags==0 && count_intra_tags==1){
3527  vertex2ptsdata[gIndex] = 0; // Interior to a myocyte
3528  }else if(count_extra_tags>=1 && count_intra_tags==1){
3529  vertex2ptsdata[gIndex] = 1; // On a membrane (1 myocyte, 1 extracellular)
3530  }else if(count_extra_tags==0 && count_intra_tags==2){
3531  vertex2ptsdata[gIndex] = 2; // On a gap junction (2 myocytes)
3532  }else if(count_extra_tags>=1 && count_intra_tags==2){
3533  vertex2ptsdata[gIndex] = 3; // On a complex junction (2 myocytes, 1 extracellular)
3534  }
3535  else if(count_intra_tags>2){
3536  // Error for unhandled cases.
3537  fprintf(stderr, "More than two intracellular are connected %d. Tags:", gIndex);
3538  for(T tag : vec_tags.tags) if(tag != -1) fprintf(stderr, " %d", tag);
3539  fprintf(stderr, "\n");
3540  EXIT(1);
3541  }
3542  else{
3543  // Error for unhandled cases.
3544  fprintf(stderr, "WARN: unhandled case in ptsData computation for gIndex %d. Tags:", gIndex);
3545  for(T tag : vec_tags.tags) if(tag != -1) fprintf(stderr, " %d", tag);
3546  fprintf(stderr, "\n");
3547  EXIT(1);
3548  }
3549  }
3550  }
3551 }
3552 
3553 
3569 template <class T, class S, class V, class emi_index_rank>
3570 inline void assemble_map_both_to_unique(SF::abstract_matrix<T, S>& op,
3572  std::pair<emi_index_rank, emi_index_rank>>& map,
3573  const SF::meshdata<T, V>& unique_mesh,
3574  const SF::meshdata<T, V>& both_mesh)
3575 {
3576  // Recover global row/column ownership so locally assembled entries can be
3577  // inserted directly with global matrix indices.
3578  int rank;
3579  MPI_Comm_rank(both_mesh.comm, &rank);
3580 
3581  SF::vector<long int> layout_both;
3582  SF::layout_from_count<long int>(both_mesh.l_numelem, layout_both, both_mesh.comm);
3583 
3584  SF::vector<long int> layout_unique;
3585  SF::layout_from_count<long int>(unique_mesh.l_numelem, layout_unique, unique_mesh.comm);
3586 
3587  SF::vector<SF_int> row_idx(1), col_idx(1);
3588  SF::dmat<SF_real> ebuff(1, 1);
3589 
3590  // Assemble one restriction row per locally owned unique face.
3591  for (const auto& [local_unique_idx, both_pair] : map) {
3592  const auto& first_both = both_pair.first;
3593  const auto& second_both = both_pair.second;
3594 
3595  // The unique-face row is owned by this rank, so its global row index comes
3596  // from the unique-face layout plus the local unique-face element id.
3597  T global_unique_row = layout_unique[rank] + local_unique_idx;
3598  row_idx[0] = global_unique_row;
3599 
3600  // Resolve the candidate both-face columns and discard invalid sides.
3601  bool first_valid = (first_both.index >= 0 && first_both.rank >= 0);
3602  bool second_valid = (second_both.index >= 0 && second_both.rank >= 0);
3603 
3604  T global_both_col_first = -1;
3605  T global_both_col_second = -1;
3606  if (first_valid) {
3607  global_both_col_first = layout_both[first_both.rank] + first_both.index;
3608  }
3609  if (second_valid) {
3610  global_both_col_second = layout_both[second_both.rank] + second_both.index;
3611  }
3612 
3613  // If both sides collapse to the same global both-face column, do not insert
3614  // the same contribution twice; treat it as a single-sided restriction row.
3615  if (first_valid && second_valid && global_both_col_first == global_both_col_second) {
3616  second_valid = false;
3617  }
3618 
3619  // Count how many distinct both-face columns contribute to this unique face.
3620  int count = 0;
3621  if (first_valid) count++;
3622  if (second_valid) count++;
3623 
3624  if (count == 0) continue; // No data to average
3625 
3626  // Use arithmetic averaging when both sides exist; otherwise preserve the
3627  // single available value without halving it.
3628  double weight = 0.5;
3629  if (count == 1) weight = 1.0;
3630 
3631  ebuff.assign(1, 1, weight);
3632 
3633  // Insert the first contributing both-face column.
3634  if (first_valid) {
3635  col_idx[0] = global_both_col_first;
3636  op.set_values(row_idx, col_idx, ebuff.data(), false);
3637  }
3638 
3639  // Insert the second contributing both-face column when present.
3640  if (second_valid) {
3641  col_idx[0] = global_both_col_second;
3642  op.set_values(row_idx, col_idx, ebuff.data(), false);
3643  }
3644  }
3645 
3646  op.finish_assembly();
3647 }
3648 
3665 template<class T>
3666 inline void restrict_to_membrane(vector<T> & v, hashmap::unordered_map<T,T> & dof2ptsData, const meshdata<int, float> & mesh)
3667 {
3668  const SF::vector<mesh_int_t> & rnod = mesh.get_numbering(SF::NBR_REF);
3671  for(size_t i=0; i<rnod.size(); i++){
3672  g2l[rnod[i]] = i;
3673  l2g[i] = rnod[i];
3674  }
3675 
3676  size_t widx = 0;
3677 
3678  for(size_t i=0; i<v.size(); i++){
3679  T g = l2g[v[i]];
3680 
3681  if (dof2ptsData[g] > 0) { // dof2ptsData requires global index
3682  v[widx++] = v[i];
3683  }
3684  }
3685 
3686  v.resize(widx);
3687 }
3688 
3689 }
3690 
3691 #endif
3692 #endif
int mesh_int_t
Definition: SF_container.h:46
#define SF_COMM
the default SlimFem MPI communicator
Definition: SF_globals.h:27
Functions related to EMI mesh IO.
Functions handling a distributed mesh.
Definition: mesher.cc:237
virtual void finish_assembly()=0
virtual void set_values(const vector< T > &row_idx, const vector< T > &col_idx, const vector< S > &vals, bool add)=0
Dense matrix class.
Definition: dense_mat.hpp:43
The mesh storage class. It contains both element and vertex data.
Definition: SF_container.h:396
size_t l_numelem
local number of elements
Definition: SF_container.h:399
MPI_Comm comm
the parallel mesh is defined on a MPI world
Definition: SF_container.h:404
A vector storing arbitrary data.
Definition: SF_vector.h:43
size_t size() const
The current size of the vector.
Definition: SF_vector.h:104
void resize(size_t n)
Resize a vector.
Definition: SF_vector.h:209
iterator find(const K &key)
Search for key. Return iterator.
Definition: hashmap.hpp:637
void sort(Compare comp=Compare())
Sort data entries.
Definition: hashmap.hpp:703
void insert(InputIterator first, InputIterator last)
Insert Iterator range.
Definition: hashmap.hpp:583
T & at(const K &key)
Data access by key.
Definition: hashmap.hpp:657
size_t size() const
Definition: hashmap.hpp:731
iterator find(const K &key)
Definition: hashmap.hpp:1092
hm_int count(const K &key) const
Definition: hashmap.hpp:1078
@ Line
Definition: filament.h:16
@ Quad
Definition: filament.h:16
@ Tri
Definition: filament.h:16
Definition: dense_mat.hpp:34
void dsp_from_cnt(const vector< T > &cnt, vector< T > &dsp)
Compute displacements from counts.
Definition: SF_vector.h:310
void interval(vector< T > &vec, size_t start, size_t end)
Create an integer interval between start and end.
Definition: SF_vector.h:350
void binary_sort_copy(vector< T > &_V, vector< S > &_W)
Definition: SF_sort.h:296
T sum(const vector< T > &vec)
Compute sum of a vector's entries.
Definition: SF_vector.h:340
void unique_resize(vector< T > &_P)
Definition: SF_sort.h:348
void count(const vector< T > &data, vector< S > &cnt)
Count number of occurrences of indices.
Definition: SF_vector.h:332
void MPI_Exchange(commgraph< T > &grph, vector< S > &send, vector< S > &recv, MPI_Comm comm)
Exchange data in parallel over MPI.
Definition: SF_network.h:47
void binary_sort(vector< T > &_V)
Definition: SF_sort.h:284
@ Tri
Definition: SF_container.h:60
@ Prism
Definition: SF_container.h:58
@ Pyramid
Definition: SF_container.h:57
@ Tetra
Definition: SF_container.h:54
@ Quad
Definition: SF_container.h:59
@ Hexa
Definition: SF_container.h:55
SF_nbr
Enumeration encoding the different supported numberings.
Definition: SF_container.h:200
@ NBR_ELEM_REF
The element numbering of the reference mesh (the one stored on HD).
Definition: SF_container.h:204
@ NBR_REF
The nodal numbering of the reference mesh (the one stored on HD).
Definition: SF_container.h:201
void dsp_from_cnt(const std::vector< T > &cnt, std::vector< T > &dsp)
Compute displacements from counts.
Definition: kdpart.hpp:140
constexpr T min(T a, T b)
Definition: ion_type.h:33
constexpr T max(T a, T b)
Definition: ion_type.h:31
static hm_uint hash(const T &a)
Definition: hashmap.hpp:92