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 "SF_mesh_utils.h"
32 #include "SF_mesh_io_emi.h"
33 
34 namespace SF {
35 
36 
37 const int MAX_INTERSECTIONS = 20;
46 struct intersection_tags {
54  int tags[MAX_INTERSECTIONS];
55  intersection_tags() {
56  std::fill_n(tags, MAX_INTERSECTIONS, -1);
57  }
58 
59 };
60 
61 struct intersection_data {
62  int tags[MAX_INTERSECTIONS];
63  int data[MAX_INTERSECTIONS];
64  int ranks[MAX_INTERSECTIONS];
65 
66  intersection_data() {
67  std::fill_n(tags, MAX_INTERSECTIONS, -1);
68  std::fill_n(data, MAX_INTERSECTIONS, -1);
69  std::fill_n(ranks, MAX_INTERSECTIONS, -1);
70  }
71 };
72 
73 struct intersection_indices {
74  int indices[MAX_INTERSECTIONS];
75  intersection_indices() {
76  std::fill_n(indices, MAX_INTERSECTIONS, -1);
77  }
78 };
87 template<class K, class V> inline
88 void assign_counter_face(hashmap::unordered_map<K, V> & map, const MPI_Comm comm)
89 {
90  size_t dsize = map.size();
91  vector<K> key_vec (dsize);
92  vector<V> value_vec (dsize);
93 
94  // make key and value vector of elements without counterparts
95  // which are from different ranks
96  size_t idx=0;
97  for(const auto & v : map) {
98  if (v.second.second.eidx == -1){
99  key_vec[idx] = v.first;
100  value_vec[idx] = v.second;
101  idx++;
102  }
103  }
104  dsize = idx;
105  key_vec.resize(dsize);
106  value_vec.resize(dsize);
107 
108  // inspired in get_hashmap_duplicates
109  int size, rank;
110  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
111 
112  vector<int> perm, dest(dsize);
113 
114  for(size_t i=0; i<dsize; i++)
115  dest[i] = hashmap::hash_ops<K>::hash(key_vec[i]) % size;
116 
117  commgraph<size_t> grph;
118  grph.configure(dest, comm);
119  size_t nrecv = sum(grph.rcnt);
120 
121  interval(perm, 0, dsize);
122  binary_sort_copy(dest, perm);
123 
124  // fill send buffer and communicate
125  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
126  vector<V> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
127  for(size_t i=0; i<dsize; i++) {
128  sbuff_key[i] = key_vec[perm[i]];
129  sbuff_value[i] = value_vec[perm[i]];
130  }
131  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
132  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
133 
134  // Merge values for identical keys on the receiving side.
135  // Each merged entry should contain both sides and preserve owner index_unique.
137 
138  // Merge all received entries for the same key.
139  for (size_t i = 0; i < nrecv; i++) {
140  auto it = rmap.find(rbuff_key[i]);
141  if(it != rmap.end()) {
142  // Add counter faces only if the tag numbers of the two faces are different.
143  // That is, the counter face must have a different tag number from the current face.
144  // However, in some special cases, both faces may belong to the extracellular domain,
145  // which are removed later.
146  if (it->second.first.tag != rbuff_value[i].first.tag) {
147  it->second.second = rbuff_value[i].first;
148  }
149  } else {
150  rmap.insert({rbuff_key[i], rbuff_value[i]});
151  }
152  }
153 
154  for (size_t i = 0; i < nrecv; i++ ) {
155  auto it = rmap.find(rbuff_key[i]);
156  if(it != rmap.end()) rbuff_value[i] = it->second;
157  }
158 
159  // Send merged values back to the originating ranks.
160  grph.transpose();
161  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
162 
163  for (size_t i = 0; i < dsize; i++ ) {
164  auto it = map.find(sbuff_key[i]);
165  if (it != map.end()) it->second = ibuff_value[i];
166  }
167 }
168 
177 template<class K, class V> inline
178 void assign_unique_first_face(hashmap::unordered_map<K, V> & map, const MPI_Comm comm)
179 {
180  // Exchange and merge interface-face pairs across ranks.
181  // Each key represents one physical interface; after this call, both sides
182  // (first/second) carry rank/index data, and the owner side has index_unique set.
183  //
184  // Steps:
185  // 1) Collect only partial entries (missing counter-face) that need MPI exchange.
186  // 2) Hash-partition these keys to a merge owner rank.
187  // 3) Exchange keys/values to owners and merge both sides for each key.
188  // 4) Send merged entries back so every rank has complete data for its keys.
189  size_t dsize = map.size();
190  vector<K> key_vec (dsize);
191  vector<V> value_vec (dsize);
192 
193  // Collect only entries without counterparts (cross-rank interfaces).
194  size_t idx=0;
195  for(const auto & v : map) {
196  if (v.second.second.index_unique == -1 && v.second.second.index_one == -1){
197  key_vec[idx] = v.first;
198  value_vec[idx] = v.second;
199  idx++;
200  }
201  }
202  dsize = idx;
203  key_vec.resize(dsize);
204  value_vec.resize(dsize);
205 
206  // Hash-partition keys to decide which rank will merge each key.
207  int size, rank;
208  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
209 
210  vector<int> perm, dest(dsize);
211 
212  for(size_t i=0; i<dsize; i++)
213  dest[i] = hashmap::hash_ops<K>::hash(key_vec[i]) % size;
214 
215  commgraph<size_t> grph;
216  grph.configure(dest, comm);
217  size_t nrecv = sum(grph.rcnt);
218 
219  interval(perm, 0, dsize);
220  binary_sort_copy(dest, perm);
221 
222  // Send keys and partial values to the owner rank for each key.
223  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
224  vector<V> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
225  for(size_t i=0; i<dsize; i++) {
226  sbuff_key[i] = key_vec[perm[i]];
227  sbuff_value[i] = value_vec[perm[i]];
228  }
229  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
230  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
231 
232  // rmap accumulates merged entries per key on the receiving rank.
233  // Multiple ranks may send the same key; we combine them into one pair.
235 
236  // Merge all received entries for the same key.
237  for (size_t i = 0; i < nrecv; i++) {
238  auto it = rmap.find(rbuff_key[i]);
239  if(it != rmap.end()) {
240  // Merge entries without losing existing owner info.
241  auto& existing = it->second;
242  using Face = std::decay_t<decltype(existing.first)>;
243 
244  auto update_slot = [&](Face& slot, const Face& src) {
245  // Merge strategy for one side:
246  // - index_unique is only set by the owner and must never be overwritten by a non-owner.
247  // - rank/index_one are filled lazily when missing.
248  if (src.index_unique >= 0) {
249  slot.index_unique = src.index_unique;
250  slot.rank = src.rank;
251  } else {
252  if (slot.rank < 0) slot.rank = src.rank;
253  }
254  if (slot.index_one < 0 && src.index_one >= 0) slot.index_one = src.index_one;
255  };
256 
257  auto same_face = [&](const Face& a, const Face& b) {
258  // Use (rank,index_one) to identify a side uniquely (index_one is local).
259  return (a.rank >= 0 && b.rank >= 0 &&
260  a.index_one >= 0 && b.index_one >= 0 &&
261  a.rank == b.rank && a.index_one == b.index_one);
262  };
263 
264  auto ingest = [&](const Face& incoming) {
265  // Ignore empty placeholders (no useful info).
266  if (incoming.index_one < 0 && incoming.index_unique < 0 && incoming.rank < 0) return;
267  // Try to match incoming with an existing side (same rank or same index_one).
268  if (same_face(existing.first, incoming) || existing.first.rank == incoming.rank) {
269  update_slot(existing.first, incoming);
270  } else if (same_face(existing.second, incoming) || existing.second.rank == incoming.rank) {
271  update_slot(existing.second, incoming);
272  } else if (existing.second.index_one == -1 && existing.second.index_unique == -1) {
273  // If the second slot is empty, place incoming there (new side).
274  existing.second = incoming;
275  } else if (existing.first.index_one == -1 && existing.first.index_unique == -1) {
276  // If the first slot is empty, place incoming there (new side).
277  existing.first = incoming;
278  } else {
279  // Both slots are filled. Only upgrade missing owner info.
280  // (Do not overwrite an existing owner assignment.)
281  if (existing.first.index_unique < 0 && incoming.index_unique >= 0) {
282  existing.first.index_unique = incoming.index_unique;
283  existing.first.rank = incoming.rank;
284  } else if (existing.second.index_unique < 0 && incoming.index_unique >= 0) {
285  existing.second.index_unique = incoming.index_unique;
286  existing.second.rank = incoming.rank;
287  }
288  }
289  };
290 
291  // Merge both sides from the received pair.
292  ingest(rbuff_value[i].first);
293  ingest(rbuff_value[i].second);
294  } else {
295  // First time we see this key on the receiving rank.
296  rmap.insert({rbuff_key[i], rbuff_value[i]});
297  }
298  }
299 
300  for (size_t i = 0; i < nrecv; i++ ) {
301  auto it = rmap.find(rbuff_key[i]);
302  if(it != rmap.end()) rbuff_value[i] = it->second;
303  }
304 
305  // sending the assigned data back to original rank
306  grph.transpose();
307  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
308 
309  for (size_t i = 0; i < dsize; i++ ) {
310  auto it = map.find(sbuff_key[i]);
311  if (it != map.end()) {
312  // Overwrite with merged entry from exchange (rmap already combined both sides)
313  it->second = ibuff_value[i];
314  } else {
315  // Insert new entry received from another rank
316  // This ensures non-owner ranks get the complete mapping data
317  map.insert({sbuff_key[i], ibuff_value[i]});
318  }
319  }
320 }
321 
329 template<class K> inline
330 void assign_counter_vertices_tuple(hashmap::unordered_map<K, intersection_data>& map, const MPI_Comm comm)
331 {
332  size_t dsize = map.size();
333  vector<K> key_vec(dsize);
334  vector<intersection_data> value_vec(dsize);
335 
336  size_t idx = 0;
337  for (const auto& v : map) {
338  key_vec[idx] = v.first;
339  value_vec[idx] = v.second;
340  idx++;
341  }
342 
343  int size, rank;
344  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
345 
346  vector<int> perm, dest(dsize);
347  for (size_t i = 0; i < dsize; i++)
348  dest[i] = hashmap::hash_ops<K>::hash(key_vec[i]) % size;
349 
350  commgraph<size_t> grph;
351  grph.configure(dest, comm);
352  size_t nrecv = sum(grph.rcnt);
353 
354  interval(perm, 0, dsize);
355  binary_sort_copy(dest, perm);
356 
357  // fill send buffer and communicate
358  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
359  vector<intersection_data> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
360  for (size_t i = 0; i < dsize; i++) {
361  sbuff_key[i] = key_vec[perm[i]];
362  sbuff_value[i] = value_vec[perm[i]];
363  }
364  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
365  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
366 
368  for (size_t i = 0; i < nrecv; i++) {
369  auto it = rmap.find(rbuff_key[i]);
370  if (it != rmap.end()) {
371  intersection_data& map_val = it->second;
372  intersection_data& r_val = rbuff_value[i];
373 
374  int first_index = -1;
375  for (int j = 0; j < MAX_INTERSECTIONS; ++j) {
376  if (map_val.tags[j] == -1 && map_val.ranks[j] == -1) {
377  first_index = j;
378  break;
379  }
380  }
381 
382  if (first_index == -1) continue; // This means all possible intersections have already been assigned.
383 
384  for (int j = 0; j < MAX_INTERSECTIONS; ++j) {
385  if (r_val.tags[j] == -1) continue;
386 
387  int rtag = r_val.tags[j];
388  int rdata = r_val.data[j];
389  int rrank = r_val.ranks[j];
390 
391  bool exist_rtag = false;
392  for (int k = 0; k < first_index; ++k) {
393  if (map_val.tags[k] == rtag && map_val.ranks[k] == rrank) {
394  exist_rtag = true;
395  break;
396  }
397  }
398 
399  if (!exist_rtag && first_index < MAX_INTERSECTIONS) {
400  map_val.tags[first_index] = rtag;
401  map_val.data[first_index] = rdata;
402  map_val.ranks[first_index] = rrank;
403  first_index++;
404  }
405  }
406  }
407  else {
408  rmap.insert({ rbuff_key[i], rbuff_value[i] });
409  }
410  }
411 
412  for (size_t i = 0; i < nrecv; i++) {
413  auto it = rmap.find(rbuff_key[i]);
414  if (it != rmap.end()) rbuff_value[i] = it->second;
415  }
416 
417  // sending the assigned data back to original rank
418  grph.transpose();
419  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
420 
421  for (size_t i = 0; i < dsize; i++) {
422  auto it = map.find(sbuff_key[i]);
423  if (it != map.end()) it->second = ibuff_value[i];
424  }
425 }
426 
434 template<class K, class V> inline
435 void assign_dof_on_counter_face(hashmap::unordered_map<K, V> & map, const MPI_Comm comm)
436 {
437 
438  int size, rank;
439  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
440 
441  size_t dsize = map.size();
442  vector<K> key_vec (dsize);
443  vector<V> value_vec (dsize);
444 
445  // make key and value vector of elements without counterparts
446  // which are from different ranks
447  size_t idx=0;
448  for(const auto & v : map) {
449  key_vec[idx] = v.first;
450  value_vec[idx] = v.second;
451  idx++;
452  }
453  dsize = idx;
454  key_vec.resize(dsize);
455  value_vec.resize(dsize);
456 
457  vector<int> perm, dest(dsize);
458 
459  for(size_t i=0; i<dsize; i++){
460  dest[i] = hashmap::hash_ops<typename K::first_type>::hash(key_vec[i].first) % size;
461  }
462 
463  commgraph<size_t> grph;
464  grph.configure(dest, comm);
465  size_t nrecv = sum(grph.rcnt);
466 
467  interval(perm, 0, dsize);
468  binary_sort_copy(dest, perm);
469 
470  // fill send buffer and communicate
471  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
472  vector<V> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
473  for(size_t i=0; i<dsize; i++) {
474  sbuff_key[i] = key_vec[perm[i]];
475  sbuff_value[i] = value_vec[perm[i]];
476  }
477  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
478  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
479 
481 
482  for (size_t i = 0; i < nrecv; i++) {
483  auto it = rmap.find(rbuff_key[i]);
484  if(it != rmap.end()) {
485  if(rbuff_value[i] != -1 and it->second == -1){
486  it->second = rbuff_value[i];// add new indices
487  }
488  } else {
489  rmap.insert({rbuff_key[i], rbuff_value[i]});
490  }
491  }
492 
493  for (size_t i = 0; i < nrecv; i++ ) {
494  auto it = rmap.find(rbuff_key[i]);
495  if(it != rmap.end()) {
496  if(rbuff_value[i] == -1)
497  rbuff_value[i] = it->second;
498  }
499  }
500 
501  // sending the assigned data back to original rank
502  grph.transpose();
503  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
504 
505  for (size_t i = 0; i < dsize; i++ ) {
506  auto it = map.find(sbuff_key[i]);
507  if (it != map.end()) it->second = ibuff_value[i];
508  }
509 }
510 
519 template<class K, class V> inline
520 void assign_petsc_on_counter_face(hashmap::unordered_map<K, V> & map, const MPI_Comm comm)
521 {
522 
523  int size, rank;
524  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
525 
526  size_t dsize = map.size();
527  vector<K> key_vec (dsize);
528  vector<V> value_vec (dsize);
529 
530  // make key and value vector of elements without counterparts
531  // which are from different ranks
532  size_t idx=0;
533  for(const auto & v : map) {
534  key_vec[idx] = v.first;
535  value_vec[idx] = v.second;
536  idx++;
537  }
538  dsize = idx;
539  key_vec.resize(dsize);
540  value_vec.resize(dsize);
541 
542  vector<int> perm, dest(dsize);
543  for(size_t i=0; i<dsize; i++){
544  dest[i] = hashmap::hash_ops<typename K::first_type>::hash(key_vec[i].first) % size;
545  }
546 
547  commgraph<size_t> grph;
548  grph.configure(dest, comm);
549  size_t nrecv = sum(grph.rcnt);
550 
551  interval(perm, 0, dsize);
552  binary_sort_copy(dest, perm);
553 
554  // fill send buffer and communicate
555  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
556  vector<V> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
557  for(size_t i=0; i<dsize; i++) {
558  sbuff_key[i] = key_vec[perm[i]];
559  sbuff_value[i] = value_vec[perm[i]];
560  }
561  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
562  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
563 
565 
566  for (size_t i = 0; i < nrecv; i++) {
567  auto it = rmap.find(rbuff_key[i]);
568  if(it != rmap.end()) {
569  if(rbuff_value[i].second != -1 and it->second.second == -1){
570  it->second = rbuff_value[i];// add new indices
571  }
572  } else {
573  rmap.insert({rbuff_key[i], rbuff_value[i]});
574  }
575  }
576 
577  for (size_t i = 0; i < nrecv; i++ ) {
578  auto it = rmap.find(rbuff_key[i]);
579  if(it != rmap.end()) {
580  if(rbuff_value[i].second == -1)
581  rbuff_value[i].second = it->second.second;
582  }
583  }
584 
585  // sending the assigned data back to original rank
586  grph.transpose();
587  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
588 
589  for (size_t i = 0; i < dsize; i++ ) {
590  auto it = map.find(sbuff_key[i]);
591  if (it != map.end()) it->second = ibuff_value[i];
592  }
593 }
594 
602 template<class T>
603 void sort_surf_local_indices(tuple<T> & ref, tuple<T> & target)
604 {
605  vector<T> buff(2);
606 
607  buff[0] = ref.v1, buff[1] = ref.v2;
608  binary_sort(buff);
609 
610  target.v1 = buff[0], target.v2 = buff[1];
611 }
612 
620 template<class T>
621 void sort_surf_local_indices(triple<T> & ref, triple<T> & target)
622 {
623  vector<T> buff(3);
624 
625  buff[0] = ref.v1, buff[1] = ref.v2, buff[2] = ref.v3;
626  binary_sort(buff);
627 
628  target.v1 = buff[0], target.v2 = buff[1], target.v3 = buff[2];
629 }
630 
638 template<class T>
639 void sort_surf_local_indices(quadruple<T> & ref, quadruple<T> & target)
640 {
641  vector<T> buff(4);
642 
643  buff[0] = ref.v1, buff[1] = ref.v2, buff[2] = ref.v3, buff[3] = ref.v4;
644  binary_sort(buff);
645 
646  target.v1 = buff[0], target.v2 = buff[1], target.v3 = buff[2], target.v4 = buff[3];
647 }
648 
656 template<class K, class V> inline
657 void insert_surf_based_Tag(V & ref,
658  hashmap::unordered_map<K, std::pair<V, V>> & surfmap)
659 {
660  K surf;
661 
662  sort_surf_local_indices(ref.points, surf);
663 
664  auto it = surfmap.find(surf);
665  if (it != surfmap.end()) {
666  // Add as a counter face only if this face exists in surfMap and its tag number
667  // differs from the tag number of the first saved pair.
668  if ( it->second.first.tag != ref.tag )
669  it->second.second = ref;
670  } else {
671  std::pair<V, V> face;
672  face.first = ref;
673  surfmap.insert({surf,face});
674  }
675 }
676 
696 template<class T, class W, class V, class U> inline
697 void insert_surf_emi( int rank,
698  SF::vector<T> const & ref_con,
699  SF::vector<T> const & ptsData,
700  const T eidx, const T tag,
701  const std::vector<T> & line_con, const std::vector<T> & surf_con, const std::vector<T> & qsurf_con,
702  hashmap::unordered_map<tuple<T>, std::pair<W, W> > & line_face,
703  hashmap::unordered_map<triple<T>, std::pair<V, V> > & surfmap,
704  hashmap::unordered_map<quadruple<T>, std::pair<U, U> > & qsurfmap)
705 {
706  W line;
707  line.eidx = eidx;
708  line.tag = tag;
709  V face;
710  face.eidx = eidx;
711  face.tag = tag;
712  U qface;
713  qface.eidx = eidx;
714  qface.tag = tag;
715 
716  // Mark faces as potential candidates if all of their vertices lie on the interface.
717  for (size_t i = 0; i < line_con.size(); i += 2) {
718  if (ptsData[line_con[i ] - 1] > 0 &&
719  ptsData[line_con[i + 1] - 1] > 0 ) {
720  line.points.v1 = ref_con[line_con[i ] - 1];
721  line.points.v2 = ref_con[line_con[i + 1] - 1];
722  line.rank = rank;
723  insert_surf_based_Tag(line, line_face);
724  }
725  }
726 
727  for (size_t i = 0; i < surf_con.size(); i += 3) {
728  if (ptsData[surf_con[i ] - 1] > 0 &&
729  ptsData[surf_con[i + 1] - 1] > 0 &&
730  ptsData[surf_con[i + 2] - 1] > 0 ) {
731  face.points.v1 = ref_con[surf_con[i ] - 1];
732  face.points.v2 = ref_con[surf_con[i + 1] - 1];
733  face.points.v3 = ref_con[surf_con[i + 2] - 1];
734  face.rank = rank;
735  insert_surf_based_Tag(face, surfmap);
736  }
737  }
738 
739  for (size_t i = 0; i < qsurf_con.size(); i += 4) {
740  if (ptsData[qsurf_con[i ] - 1] > 0 &&
741  ptsData[qsurf_con[i + 1] - 1] > 0 &&
742  ptsData[qsurf_con[i + 2] - 1] > 0 &&
743  ptsData[qsurf_con[i + 3] - 1] > 0) {
744  qface.points.v1 = ref_con[qsurf_con[i ] - 1];
745  qface.points.v2 = ref_con[qsurf_con[i + 1] - 1];
746  qface.points.v3 = ref_con[qsurf_con[i + 2] - 1];
747  qface.points.v4 = ref_con[qsurf_con[i + 3] - 1];
748  qface.rank = rank;
749  insert_surf_based_Tag(qface, qsurfmap);
750  }
751  }
752 }
753 
760 template<class T, class V>
761 struct emi_face {
762  V points;
763  T eidx = -1;
764  T tag = 0;
765  T rank = -1;
766  T mem = 1; // default the face on memebarne, otherwise will be asssign to 2 for gap-junction face
767  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
768  // we assume that the face will be selected on the smaller rank. If both has the same rank, we choose the first pair.
769  T index_unique = -1; // index of face on emi_surface_unique_face_msh
770  T index_one = -1; // index of face on emi_surface_msh
771 };
772 
773 template<class T>
774 struct emi_unique_face {
775  T index_unique = -1;
776  T index_one = -1;
777  T rank = -1;
778 };
779 
806 template<class T>
807 struct emi_index_rank {
808  T index = -1;
809  T rank = -1;
810 };
823 template<class T, class S> inline
824 void compute_surface_with_tags(meshdata<T,S> & mesh,
825  const SF_nbr numbering,
826  hashmap::unordered_map<T,T> & vertex2ptsdata,
827  hashmap::unordered_set<int> & extra_tags,
828  hashmap::unordered_map<tuple<T>,
829  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
830  hashmap::unordered_map<triple<T>,
831  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_surf,
832  hashmap::unordered_map<quadruple<T>,
833  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_surf)
834 {
835  MPI_Comm comm = SF_COMM;
836  int size, rank;
837  MPI_Comm_size(comm, &size);
838  MPI_Comm_rank(comm, &rank);
839 
840  const T* con = mesh.con.data();
841  const T* nbr = mesh.get_numbering(numbering).data();
842  const SF::vector<mesh_int_t> & rnod = mesh.get_numbering(numbering);
843 
845 
846  for(size_t i=0; i<mesh.con.size(); i++){
847  g2ptsData[rnod[con[i]]] = vertex2ptsdata[rnod[con[i]]];
848  }
849 
850  const vector<T> & ref_eidx = mesh.get_numbering(NBR_ELEM_REF);
851 
852  // Collect potential faces whose three vertices lie on the EMI interfaces.
853  // The ptsData, computed earlier from the input mesh, is used to check
854  // the intersection of each vertex with different tag numbers.
855  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
856  T tag = mesh.tag[eidx];
857  T size_elem = mesh.dsp[eidx+1]-mesh.dsp[eidx];
858 
859  vector<T> dofvec(size_elem); // reference nodes of each element
860  vector<T> ptsDatavec(size_elem); // point data on the nodes of each element
861 
862  for (int n = mesh.dsp[eidx], i = 0; n < mesh.dsp[eidx+1];n++,i++)
863  {
864  dofvec[i] = rnod[con[n]];
865  ptsDatavec[i] = g2ptsData[rnod[con[n]]];
866  }
867  std::vector<T> surf_con ;
868  std::vector<T> qsurf_con;
869  std::vector<T> line_con;
870  switch(mesh.type[eidx]) {
871  // I need to find the right order of lines in 2D
872 
873  // I guess that in the 2D case we would define the interfaces as
874  // lines instead of surfaces.
875  case Tri: {
876  line_con = {1, 2,
877  2, 3,
878  3, 1};
879  surf_con = {};
880  qsurf_con = {};
881  break;
882  }
883 
884  case Quad: {
885  line_con = {1, 2,
886  2, 3,
887  3, 4,
888  4, 1};
889  surf_con = {};
890  qsurf_con = {};
891  break;
892  }
893 
894  case Tetra: {
895  // surfaces are (2,3,1) , (1,4,2) , (2,4,3) , (1,3,4)
896  line_con = {};
897  surf_con = {2,3,1,
898  1,4,2,
899  2,4,3,
900  1,3,4};
901  qsurf_con = {};
902  break;
903  }
904  case Pyramid: {
905  // surfaces are (1,5,2) , (2,5,3) , (3,5,4) , (4,5,1) , (1,2,3,4)
906  line_con = {};
907  surf_con = {1,5,2,
908  2,5,3,
909  3,5,4,
910  4,5,1};
911  qsurf_con = {1,2,3,4};
912  break;
913  }
914  case Prism: {
915  // surfaces are (1,2,3) , (4,5,6) , (1,2,6,4) , (2,3,5,6) , (3,1,4,5)
916  line_con = {};
917  surf_con = {1,2,3,
918  4,5,6};
919  qsurf_con = {1,2,6,4,
920  2,3,5,6,
921  3,1,4,5};
922  break;
923  }
924  case Hexa: {
925  // 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)
926  line_con = {};
927  surf_con = {};
928  qsurf_con = {1,2,3,4,
929  3,2,8,7,
930  4,3,7,6,
931  1,4,6,5,
932  2,1,5,8,
933  6,7,8,5 };
934  break;
935  }
936  default:
937  fprintf(stderr, "%s error: Unsupported element in surface computation!\n", __func__);
938  exit(1);
939  }
940  insert_surf_emi(rank, dofvec, ptsDatavec, ref_eidx[eidx], mesh.tag[eidx], line_con, surf_con, qsurf_con, line_face, tri_surf, quad_surf);
941  }
942 
943  // Find the counter faces corresponding to the selected potential faces.
944  assign_counter_face(line_face, mesh.comm);
945  assign_counter_face(tri_surf, mesh.comm);
946  assign_counter_face(quad_surf, mesh.comm);
947 
948  // Removed potential faces without matching counter faces and
949  // Removed potential faces where both associated tags belong to the extracellular domain.
950  for(auto it = line_face.begin(); it != line_face.end(); ) {
951  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()))
952  it = line_face.erase(it);
953  else
954  ++it;
955  }
956 
957  for(auto it = tri_surf.begin(); it != tri_surf.end(); ) {
958  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()))
959  it = tri_surf.erase(it);
960  else
961  ++it;
962  }
963 
964  for(auto it = quad_surf.begin(); it != quad_surf.end(); ) {
965  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()))
966  it = quad_surf.erase(it);
967  else
968  ++it;
969  }
970 }
971 
985 inline bool should_take_first(int tag1, int tag2,
986  const hashmap::unordered_set<int>& extra_tags,
987  const hashmap::unordered_set<int>& intra_tags)
988 {
989  bool tag1_is_intra = (intra_tags.find(tag1) != intra_tags.end());
990  bool tag2_is_intra = (intra_tags.find(tag2) != intra_tags.end());
991  bool tag1_is_extra = (extra_tags.find(tag1) != extra_tags.end());
992  bool tag2_is_extra = (extra_tags.find(tag2) != extra_tags.end());
993 
994  // Case 1: Both are intra tags (gap junction) - take smaller tag
995  if(tag1_is_intra && tag2_is_intra) {
996  return tag1 < tag2;
997  }
998 
999  // Case 2: One is extra, one is intra (membrane) - take intra
1000  if(tag1_is_intra && tag2_is_extra) {
1001  return true; // take tag1 (intra)
1002  }
1003  if(tag2_is_intra && tag1_is_extra) {
1004  return false; // take tag2 (intra)
1005  }
1006 
1007  // Default fallback: take smaller tag
1008  return tag1 < tag2;
1009 }
1010 
1032 inline
1033 void assign_map_between_elem_oneface_and_elem_uniqueFace(int rank,
1034  emi_unique_face <mesh_int_t> first,
1035  emi_unique_face <mesh_int_t> second,
1036  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,
1037  hashmap::unordered_map<mesh_int_t, emi_index_rank<mesh_int_t>> & map_elem_oneface_to_elem_uniqueFace)
1038 {
1039  // Case 1: Both face and counter face belong to the same rank
1040  if(second.index_one==-1 && second.index_unique==-1){
1041  second.index_one = first.index_one;
1042  second.index_unique = first.index_unique;
1043  }
1044  // Case 2: Cross-rank faces - set up mapping for non-owner ranks
1045  else if(first.rank != second.rank) {
1046  // Find our local index_one (from whichever field matches our rank)
1047  int our_index_one = (first.rank == rank) ? first.index_one :
1048  (second.rank == rank) ? second.index_one : -1;
1049  // Find owner's index_unique (from whichever field has valid index_unique)
1050  int owner_index_unique = (first.index_unique >= 0) ? first.index_unique :
1051  (second.index_unique >= 0) ? second.index_unique : -1;
1052  int owner_rank = (first.index_unique >= 0) ? first.rank :
1053  (second.index_unique >= 0) ? second.rank : -1;
1054 
1055  // Only set mapping for non-owner ranks (owner already set it in STEP 3)
1056  if(our_index_one >= 0 && owner_index_unique >= 0 && owner_rank >= 0 && owner_rank != rank) {
1057  map_elem_oneface_to_elem_uniqueFace[our_index_one].index = owner_index_unique;
1058  map_elem_oneface_to_elem_uniqueFace[our_index_one].rank = owner_rank;
1059  }
1060  }
1061 
1062  // Populate map_elem_uniqueFace_to_elem_oneface (only on owner rank where index_unique >= 0)
1063  // Checked both first and second to find which one represents this rank as owner
1064  // After MPI exchange, the owner info could be in either field
1065  int local_index_unique = -1;
1066  emi_index_rank<mesh_int_t> value1, value2;
1067 
1068  if (first.index_unique >= 0 && first.rank == rank) {
1069  // This rank owns the unique face, info is in 'first'
1070  local_index_unique = first.index_unique;
1071  value1.index = first.index_one;
1072  value1.rank = first.rank;
1073  value2.index = second.index_one;
1074  value2.rank = second.rank;
1075  } else if (second.index_unique >= 0 && second.rank == rank) {
1076  // This rank owns the unique face, info is in 'second'
1077  local_index_unique = second.index_unique;
1078  value1.index = second.index_one;
1079  value1.rank = second.rank;
1080  value2.index = first.index_one;
1081  value2.rank = first.rank;
1082  }
1083 
1084  if (local_index_unique >= 0) {
1085  auto itr = map_elem_uniqueFace_to_elem_oneface.find(local_index_unique);
1086  if (itr == map_elem_uniqueFace_to_elem_oneface.end()) {
1087  std::pair <emi_index_rank<mesh_int_t>,emi_index_rank<mesh_int_t>> value;
1088  value = std::make_pair(value1, value2);
1089  map_elem_uniqueFace_to_elem_oneface.insert({local_index_unique, value});
1090  }
1091  }
1092 }
1093 
1116 template <class T, class Key, class PointsKey>
1117 inline void assign_ownership_rank_on_faces( const Key& key,
1118  std::pair<emi_face<T, PointsKey>, emi_face<T, PointsKey>>& v, // (face, counter-face)
1119  int rank,
1120  mesh_int_t idx_oneface,
1121  mesh_int_t& lsize_unique,
1122  const hashmap::unordered_set<int>& extra_tags,
1123  const hashmap::unordered_set<int>& intra_tags,
1124  hashmap::unordered_map<Key,std::pair<emi_unique_face<mesh_int_t>, emi_unique_face<mesh_int_t>>>& unique_face_to_elements,
1125  hashmap::unordered_map<mesh_int_t, emi_index_rank<mesh_int_t>>& map_elem_oneface_to_elem_uniqueFace)
1126 {
1127  const bool same_rank = (v.first.rank == v.second.rank);
1128 
1129  // Decide owner rank in an order-independent way:
1130  // - If one side is intra, always take the intra side
1131  // - If both intra, take smaller tag (tie -> smaller rank)
1132  // - Otherwise, take smaller tag (tie -> smaller rank)
1133  const bool first_intra = (intra_tags.find(v.first.tag) != intra_tags.end());
1134  const bool second_intra = (intra_tags.find(v.second.tag) != intra_tags.end());
1135 
1136  // Decide which side to take based purely on tags/intra (order-independent)
1137  bool take_first_by_tags = true;
1138  if (first_intra != second_intra) {
1139  take_first_by_tags = first_intra; // take intra side
1140  } else if (v.first.tag != v.second.tag) {
1141  take_first_by_tags = (v.first.tag < v.second.tag);
1142  } else {
1143  // Same tag (should be rare for interfaces); fall back to smaller rank to be deterministic
1144  take_first_by_tags = (v.first.rank <= v.second.rank);
1145  }
1146 
1147  int owner_rank = take_first_by_tags ? v.first.rank : v.second.rank;
1148 
1149  const bool i_am_owner = (rank == owner_rank);
1150 
1151  // If I'm the owner, mark the chosen side and set mapping now.
1152  if (i_am_owner) {
1153  auto& chosen = take_first_by_tags ? v.first : v.second;
1154 
1155  chosen.mark_to_take = true;
1156  chosen.index_unique = lsize_unique;
1157  chosen.index_one = idx_oneface;
1158 
1159  map_elem_oneface_to_elem_uniqueFace[idx_oneface].index = lsize_unique;
1160  map_elem_oneface_to_elem_uniqueFace[idx_oneface].rank = rank;
1161  }
1162 
1163  // Insert per-face record once:
1164  // first = "my side if I'm owner else -1" + rank info
1165  // second = placeholder for the other rank (or same rank if same_rank)
1166  auto it = unique_face_to_elements.find(key);
1167  if (it == unique_face_to_elements.end()) {
1168  emi_unique_face<mesh_int_t> first;
1169  emi_unique_face<mesh_int_t> second;
1170 
1171  if (i_am_owner) {
1172  // Owner knows its local unique index.
1173  first.index_unique = lsize_unique;
1174  first.index_one = idx_oneface;
1175  first.rank = rank;
1176 
1177  // "other side" rank:
1178  second.index_unique = -1;
1179  second.index_one = -1;
1180  second.rank = same_rank ? rank : (rank == v.first.rank ? v.second.rank : v.first.rank);
1181  } else {
1182  // Non-owner does not know unique index yet; MPI step will fill it.
1183  first.index_unique = -1;
1184  first.index_one = idx_oneface;
1185  first.rank = rank;
1186 
1187  second.index_unique = -1;
1188  second.index_one = -1;
1189  second.rank = owner_rank; // this is the rank that will eventually provide index_unique
1190  }
1191 
1192  unique_face_to_elements.insert({ key, std::make_pair(first, second) });
1193  } else {
1194  // Update existing entry to avoid losing owner info when non-owner inserted first.
1195  auto& existing = it->second;
1196 
1197  // Choose a slot to update for this rank.
1198  emi_unique_face<mesh_int_t>* slot = nullptr;
1199  if (existing.first.rank == rank || existing.first.rank < 0) {
1200  slot = &existing.first;
1201  } else if (existing.second.rank == rank || existing.second.rank < 0) {
1202  slot = &existing.second;
1203  } else if (existing.first.index_one == idx_oneface) {
1204  slot = &existing.first;
1205  } else if (existing.second.index_one == idx_oneface) {
1206  slot = &existing.second;
1207  } else if (existing.first.index_unique < 0) {
1208  slot = &existing.first;
1209  } else if (existing.second.index_unique < 0) {
1210  slot = &existing.second;
1211  }
1212 
1213  if (slot) {
1214  if (slot->rank < 0) slot->rank = rank;
1215  if (slot->index_one < 0) slot->index_one = idx_oneface;
1216  if (i_am_owner && slot->index_unique < 0) {
1217  slot->index_unique = lsize_unique;
1218  }
1219  }
1220  }
1221 
1222  // Only the owning rank increments lsize_unique
1223  if (i_am_owner) {
1224  lsize_unique += 1;
1225  }
1226 }
1227 
1228 
1289 template<class T, class S> inline
1290 void extract_face_based_tags(meshdata<T,S> & mesh,
1291  const SF_nbr numbering,
1292  hashmap::unordered_map<T,T> & vertex2ptsdata,
1293  hashmap::unordered_map<tuple<T> ,
1294  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
1295  hashmap::unordered_map<triple<T>,
1296  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
1297  hashmap::unordered_map<quadruple<T>,
1298  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face,
1299  hashmap::unordered_set<int> & extra_tags,
1300  hashmap::unordered_set<int> & intra_tags,
1301  std::set<std::pair<mesh_int_t,mesh_int_t>> &membraneFace_default,
1302  std::set<std::pair<mesh_int_t,mesh_int_t>> &gapjunctionFace_default,
1303  meshdata<T,S> & surfmesh,
1304  meshdata<T,S> & surfmesh_w_counter,
1305  meshdata<T,S> & surfmesh_unique_face,
1306  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,
1307  hashmap::unordered_map<mesh_int_t, emi_index_rank<mesh_int_t>> & map_elem_oneface_to_elem_uniqueFace)
1308 {
1309  // ============================================================================
1310  // STEP 1: Initialize surface meshes and extract interface faces
1311  // ============================================================================
1312  surfmesh.register_numbering(SF::NBR_REF);
1313  surfmesh_w_counter.register_numbering(SF::NBR_REF);
1314  surfmesh_unique_face.register_numbering(SF::NBR_REF);
1315 
1316  // Extract faces with counter faces on the interfaces to build surface mesh
1317 
1318  compute_surface_with_tags(mesh, numbering, vertex2ptsdata, extra_tags, line_face, tri_face, quad_face);
1319 
1320  int size, rank;
1321  MPI_Comm_size(mesh.comm, &size);
1322  MPI_Comm_rank(mesh.comm, &rank);
1323 
1324  long int g_num_line = line_face.size();
1325  long int g_num_tri = tri_face.size();
1326  long int g_num_quad = quad_face.size();
1327  MPI_Allreduce(MPI_IN_PLACE, &g_num_line, 1, MPI_LONG, MPI_SUM, mesh.comm);
1328  MPI_Allreduce(MPI_IN_PLACE, &g_num_tri, 1, MPI_LONG, MPI_SUM, mesh.comm);
1329  MPI_Allreduce(MPI_IN_PLACE, &g_num_quad, 1, MPI_LONG, MPI_SUM, mesh.comm);
1330 
1331  // Warn if any extra-extra faces remain (should be none after filtering)
1332  {
1333  long int local_extra_extra = 0;
1334  for (const auto& [key, v] : line_face) {
1335  if (extra_tags.find(v.first.tag) != extra_tags.end() &&
1336  extra_tags.find(v.second.tag) != extra_tags.end()) {
1337  local_extra_extra++;
1338  }
1339  }
1340  for (const auto& [key, v] : tri_face) {
1341  if (extra_tags.find(v.first.tag) != extra_tags.end() &&
1342  extra_tags.find(v.second.tag) != extra_tags.end()) {
1343  local_extra_extra++;
1344  }
1345  }
1346  for (const auto& [key, v] : quad_face) {
1347  if (extra_tags.find(v.first.tag) != extra_tags.end() &&
1348  extra_tags.find(v.second.tag) != extra_tags.end()) {
1349  local_extra_extra++;
1350  }
1351  }
1352  long int global_extra_extra = 0;
1353  MPI_Allreduce(&local_extra_extra, &global_extra_extra, 1, MPI_LONG, MPI_SUM, mesh.comm);
1354  if (global_extra_extra > 0 && rank == 0) {
1355  fprintf(stderr, "WARN: extra-extra faces detected before filtering: %ld\n", global_extra_extra);
1356  }
1357  }
1358 
1359  surfmesh.g_numelem = g_num_line + g_num_tri + g_num_quad;
1360  surfmesh.l_numelem = line_face.size() + tri_face.size() + quad_face.size();
1361 
1362  vector<T> cnt(surfmesh.l_numelem, 0);
1363  surfmesh.tag.resize(surfmesh.l_numelem);
1364  surfmesh.type.resize(surfmesh.l_numelem);
1365  surfmesh.con.resize(line_face.size() * 2 + tri_face.size() * 3 + quad_face.size() * 4);
1366 
1367  // ============================================================================
1368  // STEP 2: Sort face keys for deterministic ordering across MPI ranks
1369  // ============================================================================
1370  // To ensure deterministic order, we must sort the keys of the hashmaps before iterating.
1371  // This is crucial for MPI consistency - all ranks must process faces in the same order.
1372  vector<tuple<T>> line_keys;
1373  for(auto const& [key, val] : line_face) line_keys.push_back(key);
1374  std::sort(line_keys.begin(), line_keys.end());
1375 
1376  vector<triple<T>> tri_keys;
1377  for(auto const& [key, val] : tri_face) tri_keys.push_back(key);
1378  std::sort(tri_keys.begin(), tri_keys.end());
1379 
1380  vector<quadruple<T>> quad_keys;
1381  for(auto const& [key, val] : quad_face) quad_keys.push_back(key);
1382  std::sort(quad_keys.begin(), quad_keys.end());
1383 
1384  // Temporary hashmaps to track unique-face assignments before MPI communication
1385  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;
1386  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;
1387  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;
1388 
1389  // Counters for local element indices:
1390  // - idx: index into surfmesh (one-side mesh)
1391  // - lsize_*: count of elements in surfmesh_w_counter
1392  // - lsize_unique_*: count of elements in surfmesh_unique_face
1393  mesh_int_t idx = 0, cidx = 0;
1394  mesh_int_t lsize_line = 0;
1395  mesh_int_t lsize_tri = 0;
1396  mesh_int_t lsize_quad = 0;
1397 
1398  mesh_int_t lsize_unique_line = 0;
1399  mesh_int_t lsize_unique_tri = 0;
1400  mesh_int_t lsize_unique_quad = 0;
1401 
1402  // ============================================================================
1403  // STEP 3: Process each face and determine unique-face ownership
1404  // ============================================================================
1405  // For each interface face, we:
1406  // 1. Add it to surfmesh
1407  for(auto const& key : line_keys) {
1408  auto& v = line_face.at(key);
1409  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
1410  if (!local_involved) {
1411  continue;
1412  }
1413  cnt[idx] = 2;
1414 
1415  // check if the face is located on gapjuntion
1416  if ((intra_tags.find(v.first.tag) != intra_tags.end()) && (intra_tags.find(v.second.tag) != intra_tags.end()))
1417  {
1418  gapjunctionFace_default.insert(std::make_pair(v.first.tag,v.second.tag));
1419  gapjunctionFace_default.insert(std::make_pair(v.second.tag,v.first.tag));
1420  v.first.mem = 2; // to set the face located on gap-junction
1421  v.second.mem = 2; // to set the face located on gap-junction
1422  }else{
1423  membraneFace_default.insert(std::make_pair(v.first.tag,v.second.tag));
1424  membraneFace_default.insert(std::make_pair(v.second.tag,v.first.tag));
1425  // The default value of .mem = 1 is already set for membrane.
1426  }
1427 
1428  // changing the strategy to select the face with the one with the same rank
1429  emi_face<T,tuple<T>> surf_neighbor =
1430  (v.first.rank == rank) ? v.first : v.second;
1431 
1432  surfmesh.type[idx] = Line;
1433  surfmesh.tag[idx] = surf_neighbor.tag;
1434  surfmesh.con[cidx + 0] = surf_neighbor.points.v1;
1435  surfmesh.con[cidx + 1] = surf_neighbor.points.v2;
1436 
1437  // If the ranks are equal, both faces will be added later.
1438  if(v.first.rank == v.second.rank)
1439  lsize_line+=2;
1440  else
1441  lsize_line+=1;
1442 
1443  assign_ownership_rank_on_faces<T, tuple<T>, tuple<T>>( key, v, rank,
1444  idx,
1445  lsize_unique_line,
1446  extra_tags, intra_tags,
1447  line_unique_face_to_elements,
1448  map_elem_oneface_to_elem_uniqueFace);
1449  idx += 1;
1450  cidx += 2;
1451  }
1452 
1453  for(auto const& key : tri_keys) {
1454  auto& v = tri_face.at(key);
1455  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
1456  if (!local_involved) {
1457  continue;
1458  }
1459  cnt[idx] = 3;
1460 
1461  // check if the face is located on gapjuntion
1462  if ((intra_tags.find(v.first.tag) != intra_tags.end()) && (intra_tags.find(v.second.tag) != intra_tags.end()))
1463  {
1464  gapjunctionFace_default.insert(std::make_pair(v.first.tag,v.second.tag));
1465  gapjunctionFace_default.insert(std::make_pair(v.second.tag,v.first.tag));
1466  v.first.mem = 2; // to set the face located on gap-junction
1467  v.second.mem = 2; // to set the face located on gap-junction
1468  }else{
1469  membraneFace_default.insert(std::make_pair(v.first.tag,v.second.tag));
1470  membraneFace_default.insert(std::make_pair(v.second.tag,v.first.tag));
1471  // The default value of .mem = 1 is already set for membrane.
1472  }
1473 
1474  // changing the strategy to select the face with the one with the same rank
1475  emi_face<T,triple<T>> surf_neighbor =
1476  (v.first.rank == rank) ? v.first : v.second;
1477 
1478  surfmesh.type[idx] = Tri;
1479  surfmesh.tag[idx] = surf_neighbor.tag;
1480  surfmesh.con[cidx + 0] = surf_neighbor.points.v1;
1481  surfmesh.con[cidx + 1] = surf_neighbor.points.v2;
1482  surfmesh.con[cidx + 2] = surf_neighbor.points.v3;
1483 
1484  // If the ranks are equal, both faces will be added later.
1485  if(v.first.rank == v.second.rank)
1486  lsize_tri+=2;
1487  else
1488  lsize_tri+=1;
1489 
1490  assign_ownership_rank_on_faces<T, triple<T>, triple<T>>( key, v, rank,
1491  idx,
1492  lsize_unique_tri,
1493  extra_tags, intra_tags,
1494  tri_unique_face_to_elements,
1495  map_elem_oneface_to_elem_uniqueFace);
1496  idx += 1;
1497  cidx += 3;
1498  }
1499 
1500  for(auto const& key : quad_keys) {
1501  auto& v = quad_face.at(key);
1502  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
1503  if (!local_involved) {
1504  continue;
1505  }
1506  cnt[idx] = 4;
1507 
1508  // check if the face is located on gapjuntion
1509  if ((intra_tags.find(v.first.tag) != intra_tags.end()) && (intra_tags.find(v.second.tag) != intra_tags.end()))
1510  {
1511  gapjunctionFace_default.insert(std::make_pair(v.first.tag,v.second.tag));
1512  gapjunctionFace_default.insert(std::make_pair(v.second.tag,v.first.tag));
1513  v.first.mem = 2; // to set the face located on gap-junction
1514  v.second.mem = 2; // to set the face located on gap-junction
1515  }else{
1516  membraneFace_default.insert(std::make_pair(v.first.tag,v.second.tag));
1517  membraneFace_default.insert(std::make_pair(v.second.tag,v.first.tag));
1518  // The default value of .mem = 1 is already set for membrane.
1519  }
1520 
1521  // changing the strategy to select the face with the one with the same rank
1522  emi_face<T,quadruple<T>> qsurf_neighbor =
1523  (v.first.rank == rank) ? v.first : v.second;
1524 
1525  surfmesh.type[idx] = Quad;
1526  surfmesh.tag[idx] = qsurf_neighbor.tag;
1527  surfmesh.con[cidx + 0] = qsurf_neighbor.points.v1;
1528  surfmesh.con[cidx + 1] = qsurf_neighbor.points.v2;
1529  surfmesh.con[cidx + 2] = qsurf_neighbor.points.v3;
1530  surfmesh.con[cidx + 3] = qsurf_neighbor.points.v4;
1531 
1532  // If the ranks are equal, both faces will be added later.
1533  if(v.first.rank == v.second.rank)
1534  lsize_quad+=2;
1535  else
1536  lsize_quad+=1;
1537 
1538  assign_ownership_rank_on_faces<T, quadruple<T>, quadruple<T>>( key, v, rank,
1539  idx,
1540  lsize_unique_quad,
1541  extra_tags, intra_tags,
1542  quad_unique_face_to_elements,
1543  map_elem_oneface_to_elem_uniqueFace);
1544  idx += 1;
1545  cidx += 4;
1546  }
1547 
1548  surfmesh.l_numelem = idx;
1549  surfmesh.tag.resize(idx);
1550  surfmesh.type.resize(idx);
1551  cnt.resize(idx);
1552  dsp_from_cnt(cnt, surfmesh.dsp);
1553  surfmesh.con.resize(cidx);
1554 
1555  long int g_num_surf = surfmesh.l_numelem;
1556  MPI_Allreduce(MPI_IN_PLACE, &g_num_surf, 1, MPI_LONG, MPI_SUM, mesh.comm);
1557  surfmesh.g_numelem = g_num_surf;
1558 
1559  if(rank ==0){
1560  std::cout << "surfmesh.g_numelem: " << surfmesh.g_numelem <<std::endl;
1561  std::cout << "surfmesh.l_numelem: " << surfmesh.l_numelem <<std::endl;
1562  }
1563 
1564  {
1565  long int g_num_w_counter_line = lsize_line;
1566  long int g_num_w_counter_tri = lsize_tri;
1567  long int g_num_w_counter_quad = lsize_quad;
1568 
1569  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_counter_line, 1, MPI_LONG, MPI_SUM, SF_COMM);
1570  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_counter_tri, 1, MPI_LONG, MPI_SUM, SF_COMM);
1571  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_counter_quad, 1, MPI_LONG, MPI_SUM, SF_COMM);
1572 
1573  // Assign global and local elements for the surface mesh used in ionic computation.
1574  surfmesh_w_counter.g_numelem = g_num_w_counter_line+g_num_w_counter_tri+g_num_w_counter_quad;
1575  surfmesh_w_counter.l_numelem = lsize_line + lsize_tri + lsize_quad;
1576 
1577  surfmesh_w_counter.tag.resize(surfmesh_w_counter.l_numelem);
1578  surfmesh_w_counter.type.resize(surfmesh_w_counter.l_numelem);
1579  surfmesh_w_counter.con.resize(lsize_line * 2 + lsize_tri * 3 + lsize_quad * 4);
1580 
1581  if(rank ==0){
1582  std::cout << "surfmesh_w_counter.g_numelem: " << surfmesh_w_counter.g_numelem <<std::endl;
1583  std::cout << "surfmesh_w_counter.l_numelem: " << surfmesh_w_counter.l_numelem <<std::endl;
1584  }
1585  }
1586 
1587  {
1588  long int g_num_w_unique_line = lsize_unique_line;
1589  long int g_num_w_unique_tri = lsize_unique_tri;
1590  long int g_num_w_unique_quad = lsize_unique_quad;
1591 
1592  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_unique_line, 1, MPI_LONG, MPI_SUM, SF_COMM);
1593  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_unique_tri, 1, MPI_LONG, MPI_SUM, SF_COMM);
1594  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_unique_quad, 1, MPI_LONG, MPI_SUM, SF_COMM);
1595 
1596  // Assign global and local elements for the surface mesh used in ionic computation.
1597  surfmesh_unique_face.g_numelem = g_num_w_unique_line+g_num_w_unique_tri+g_num_w_unique_quad;
1598  surfmesh_unique_face.l_numelem = lsize_unique_line + lsize_unique_tri + lsize_unique_quad;
1599 
1600  surfmesh_unique_face.tag.resize(surfmesh_unique_face.l_numelem);
1601  surfmesh_unique_face.type.resize(surfmesh_unique_face.l_numelem);
1602  surfmesh_unique_face.con.resize(lsize_unique_line * 2 + lsize_unique_tri * 3 + lsize_unique_quad * 4);
1603 
1604  if(rank ==0){
1605  std::cout << "surfmesh_unique_face.g_numelem: " << surfmesh_unique_face.g_numelem <<std::endl;
1606  std::cout << "surfmesh_unique_face.l_numelem: " << surfmesh_unique_face.l_numelem <<std::endl;
1607  }
1608  }
1609 
1610  // ============================================================================
1611  // STEP 4: MPI communication to share unique-face ownership between ranks
1612  // ============================================================================
1613  // When faces span multiple ranks, only one rank "owns" the unique-face.
1614  // assign_unique_first_face communicates the owner's index_unique to non-owner ranks
1615  // so they can set up their map_elem_oneface_to_elem_uniqueFace mappings correctly.
1616 
1617  // DEBUG: Before exchange
1618  #ifdef EMI_DEBUG_MESH
1619  {
1620  fprintf(stderr, "RANK %d BEFORE exchange: line_hash=%zu, tri_hash=%zu, quad_hash=%zu (expected unique total=%ld)\n",
1621  rank, line_unique_face_to_elements.size(), tri_unique_face_to_elements.size(),
1622  quad_unique_face_to_elements.size(), lsize_unique_line + lsize_unique_tri + lsize_unique_quad);
1623  fflush(stderr);
1624  }
1625  #endif
1626 
1627  if(size>1)
1628  {
1629  // Find the counter faces corresponding to the selected potential faces.
1630  assign_unique_first_face(line_unique_face_to_elements, mesh.comm);
1631  assign_unique_first_face(tri_unique_face_to_elements, mesh.comm);
1632  assign_unique_first_face(quad_unique_face_to_elements, mesh.comm);
1633  }
1634 
1635  // DEBUG: After exchange
1636  #ifdef EMI_DEBUG_MESH
1637  {
1638  int valid_line = 0, valid_tri = 0, valid_quad = 0;
1639  for (const auto& [key, val] : line_unique_face_to_elements) {
1640  if (val.first.index_unique >= 0 || val.second.index_unique >= 0) valid_line++;
1641  }
1642  for (const auto& [key, val] : tri_unique_face_to_elements) {
1643  if (val.first.index_unique >= 0 || val.second.index_unique >= 0) valid_tri++;
1644  }
1645  for (const auto& [key, val] : quad_unique_face_to_elements) {
1646  if (val.first.index_unique >= 0 || val.second.index_unique >= 0) valid_quad++;
1647  }
1648  fprintf(stderr, "RANK %d AFTER exchange: line_hash=%zu, tri_hash=%zu, quad_hash=%zu\n",
1649  rank, line_unique_face_to_elements.size(), tri_unique_face_to_elements.size(),
1650  quad_unique_face_to_elements.size());
1651  fprintf(stderr, "RANK %d VALID entries (index_unique>=0): line=%d, tri=%d, quad=%d\n",
1652  rank, valid_line, valid_tri, valid_quad);
1653  fflush(stderr);
1654  }
1655  #endif
1656 
1657  // ============================================================================
1658  // STEP 5: Finalize mappings after MPI communication
1659  // ============================================================================
1660  for(auto it = line_unique_face_to_elements.begin(); it != line_unique_face_to_elements.end(); ++it) {
1661 
1662  auto& first = it->second.first;
1663  auto& second = it->second.second;
1664 
1665  assign_map_between_elem_oneface_and_elem_uniqueFace(rank, first, second, map_elem_uniqueFace_to_elem_oneface, map_elem_oneface_to_elem_uniqueFace);
1666  }
1667 
1668  for(auto it = tri_unique_face_to_elements.begin(); it != tri_unique_face_to_elements.end(); ++it) {
1669  auto& first = it->second.first;
1670  auto& second = it->second.second;
1671 
1672  assign_map_between_elem_oneface_and_elem_uniqueFace(rank, first, second, map_elem_uniqueFace_to_elem_oneface, map_elem_oneface_to_elem_uniqueFace);
1673  }
1674 
1675  for(auto it = quad_unique_face_to_elements.begin(); it != quad_unique_face_to_elements.end(); ++it) {
1676  auto& first = it->second.first;
1677  auto& second = it->second.second;
1678 
1679  assign_map_between_elem_oneface_and_elem_uniqueFace(rank, first, second, map_elem_uniqueFace_to_elem_oneface, map_elem_oneface_to_elem_uniqueFace);
1680  }
1681 }
1682 
1692 template<class T> inline
1693 void complete_map_vertex_to_dof_with_counter_face(hashmap::unordered_map<tuple<T>,
1694  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
1695  hashmap::unordered_map<triple<T>,
1696  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
1697  hashmap::unordered_map<quadruple<T>,
1698  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face,
1699  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
1700  mesh_int_t> & map_vertex_tag_to_dof)
1701 {
1702  MPI_Comm comm = SF_COMM;
1703  int size, rank;
1704  MPI_Comm_size(comm, &size);
1705  MPI_Comm_rank(comm, &rank);
1706 
1707  // look for counter part line_inteface to add to the map key<vertex,tag> -> value<dof>
1708  for(const auto & v : line_face) {
1709 
1710  // if the ranks are equal, then it has already added to the map
1711  if(v.second.first.rank == v.second.second.rank){
1712  continue;
1713  }
1714  // otherwise we look for the counter part with different rank
1715  emi_face<T,tuple<T>> surf_neighbor =
1716  (v.second.first.rank != rank) ? v.second.first : v.second.second;
1717 
1718  {
1719  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1720  Index_tag_old = std::make_pair(surf_neighbor.points.v1,surf_neighbor.tag);
1721  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1722  }
1723 
1724  {
1725  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1726  Index_tag_old = std::make_pair(surf_neighbor.points.v2,surf_neighbor.tag);
1727  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1728  }
1729  }
1730 
1731  // look for counter part line_inteface to add to the map key<vertex,tag> -> value<dof>
1732  for(const auto & v : tri_face) {
1733 
1734  // if the ranks are equal, then it has already added to the map
1735  if(v.second.first.rank == v.second.second.rank){
1736  continue;
1737  }
1738  // otherwise we look for the counter part with different rank
1739  emi_face<T,triple<T>> surf_neighbor =
1740  (v.second.first.rank != rank) ? v.second.first : v.second.second;
1741 
1742  {
1743  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1744  Index_tag_old = std::make_pair(surf_neighbor.points.v1,surf_neighbor.tag);
1745  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1746  }
1747 
1748  {
1749  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1750  Index_tag_old = std::make_pair(surf_neighbor.points.v2,surf_neighbor.tag);
1751  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1752  }
1753 
1754  {
1755  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1756  Index_tag_old = std::make_pair(surf_neighbor.points.v3,surf_neighbor.tag);
1757  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1758  }
1759  }
1760 
1761  // look for counter part line_inteface to add to the map key<vertex,tag> -> value<dof>
1762  for(const auto & v : quad_face) {
1763 
1764  // if the ranks are equal, then it has already added to the map
1765  if(v.second.first.rank == v.second.second.rank){
1766  continue;
1767  }
1768  // otherwise we look for the counter part with different rank
1769  emi_face<T,quadruple<T>> qsurf_neighbor =
1770  (v.second.first.rank != rank) ? v.second.first : v.second.second;
1771 
1772  {
1773  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1774  Index_tag_old = std::make_pair(qsurf_neighbor.points.v1,qsurf_neighbor.tag);
1775  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1776  }
1777 
1778  {
1779  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1780  Index_tag_old = std::make_pair(qsurf_neighbor.points.v2,qsurf_neighbor.tag);
1781  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1782  }
1783 
1784  {
1785  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1786  Index_tag_old = std::make_pair(qsurf_neighbor.points.v3,qsurf_neighbor.tag);
1787  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1788  }
1789 
1790  {
1791  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1792  Index_tag_old = std::make_pair(qsurf_neighbor.points.v4,qsurf_neighbor.tag);
1793  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1794  }
1795  }
1796 
1797  // assign the counter faces which we haven't assigned it yet the new dof
1798  assign_dof_on_counter_face(map_vertex_tag_to_dof,comm);
1799 }
1800 
1830 template<class T, class S> inline
1831 void compute_map_vertex_to_dof(meshdata<T,S> & mesh,
1832  const SF_nbr numbering,
1833  hashmap::unordered_map<T,T> & vertex2ptsdata,
1834  hashmap::unordered_set<int> & extra_tags,
1835  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
1836  mesh_int_t> & map_vertex_tag_to_dof)
1837 {
1838  MPI_Comm comm = mesh.comm;
1839  int size, rank;
1840  MPI_Comm_size(comm, &size);
1841  MPI_Comm_rank(comm, &rank);
1842 
1843  hashmap::unordered_map<mesh_int_t, intersection_data> map_vertex_to_tags_data_ranks;
1844 
1845  const T* con = mesh.con.data();
1846  const T* nbr = mesh.get_numbering(numbering).data();
1847  const SF::vector<T> & rnod = mesh.get_numbering(numbering);
1848 
1849  hashmap::unordered_map<T,T> g2ptsData;
1850 
1851  for(size_t i=0; i<mesh.con.size(); i++)
1852  {
1853  g2ptsData[rnod[con[i]]] = vertex2ptsdata[rnod[con[i]]];
1854  }
1855 
1856  //-----------------------------------------------------------------
1857  // Step 1: Local Information Gathering
1858  //-----------------------------------------------------------------
1859  // - Identify Simple Cases: For any vertex that is not on an interface (i.e., it's completely inside the extracellular space or completely
1860  // inside a unique myocyte), the mapping is simple. The new DOF index is the same as the original vertex index. The function populates the
1861  // output map, map_vertex_tag_to_dof, with these direct mappings, e.g., (vertex {1}, tag_3) -> 1.
1862  // - 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
1863  // MPI ranks that share it, based on the elements the current process knows about. This information is stored in a temporary map called
1864  // map_vertex_to_tags_data_ranks.
1865  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++)
1866  {
1867  T tag = mesh.tag[eidx];
1868  auto pos_extra = extra_tags.find(tag);
1869  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
1870  {
1871  T gIndex = rnod[con[n]];
1872  T data_on_gIndex = g2ptsData[rnod[con[n]]];
1873 
1874  // This handles inner DoFs, either within the extracellular domain or vertices belonging to the intracellular domain.
1875  // and marked as a visited dofs and saved key<gIndex,tag)> -> val<gIndex>
1876  if((pos_extra != extra_tags.end()) || data_on_gIndex==0){
1877  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = std::make_pair(gIndex,tag);
1878  if (map_vertex_tag_to_dof.find(Index_tag_old) == map_vertex_tag_to_dof.end() )
1879  {
1880  map_vertex_tag_to_dof.insert({Index_tag_old,gIndex});
1881  }
1882  }
1883 
1884  // When gIndex is located on an interface — membrane, gap junction, or their intersection.
1885  auto it = map_vertex_to_tags_data_ranks.find(gIndex);
1886  if (it != map_vertex_to_tags_data_ranks.end() )
1887  {
1888  intersection_data& value = it->second;
1889  bool check_tag_rank = true;
1890  for(int i=0; i<MAX_INTERSECTIONS; ++i) {
1891  if(value.tags[i] == tag && value.ranks[i] == rank) {
1892  check_tag_rank = false;
1893  break;
1894  }
1895  }
1896 
1897  // Insert the tag and corresponding rank if they haven’t been inserted yet.
1898  if(check_tag_rank){
1899  for(int i=0; i<MAX_INTERSECTIONS; ++i) {
1900  if(value.tags[i] == -1 && value.ranks[i] == -1) {
1901  value.tags[i] = tag;
1902  value.data[i] = data_on_gIndex;
1903  value.ranks[i] = rank;
1904  break;
1905  }
1906  }
1907  }
1908  }
1909  else
1910  {
1911  // insert gIndex to the map_vertex_to_tags_data_ranks
1912  intersection_data value;
1913  value.tags[0] = tag;
1914  value.data[0] = data_on_gIndex;
1915  value.ranks[0] = rank;
1916  map_vertex_to_tags_data_ranks.insert({gIndex,value});
1917  }
1918  }
1919  }
1920  // until here: map_vertex_to_tags_data_ranks knows about the interfaces for the elements it owns,
1921  // but it doesn't know if a vertex it owns is also part of an interface on another process.
1922 
1923  //-----------------------------------------------------------------
1924  // Step 2: Globalize Information with MPI
1925  //-----------------------------------------------------------------
1926  // - assign_counter_vertices_tuple: This helper function is called to manage the communication. It takes the partial
1927  // map_vertex_to_tags_data_ranks from each process and uses MPI (specifically MPI_Exchange, which is likely a wrapper around routines like
1928  // MPI_Alltoallv) to send and receive data.
1929  // - Build Global Map: After the communication, the partial maps are merged. The result is that the map_vertex_to_tags_data_ranks on every
1930  // process now contains the complete sharing information for every vertex in the entire mesh. Each process now knows the full list of tags
1931  // and ranks associated with any given vertex.
1932  assign_counter_vertices_tuple(map_vertex_to_tags_data_ranks, comm);
1933 
1934  //-----------------------------------------------------------------
1935  // Step 3: Identify and Count DOFs to Be Created
1936  //-----------------------------------------------------------------
1937  // Now that every process has the same global information, they can independently and deterministically decide which new DOFs to create.
1938  // - Analyze Intersections: The code iterates through the now-global map_vertex_to_tags_data_ranks. For each vertex, it inspects the list of
1939  // tags that share it to classify the type of interface.
1940  // - Apply Rules: Based on the interface type, it decides which side needs a new, unique DOF index.
1941  // * Membrane (Myocyte-Extracellular): The extracellular DOF keeps the original vertex index. The myocyte DOF is marked as needing a new, unique index.
1942  // * Gap Junction (Myocyte-Myocyte): To separate the two myocytes, one of them needs a new index. A deterministic rule (e.g., the myocyte
1943  // with the higher tag number gets the new index) is applied to ensure all processes make the same decision.
1944  // * Complex Junctions: For even more complex intersections (e.g., where multiple myocytes and the extracellular space meet), similar
1945  // rules are applied to create the necessary new DOFs.
1946  // - Count Local Contribution: Each process counts how many new DOFs (shift) it is responsible for creating based on these rules.
1947  int shift = 0;
1949  for(const auto & key : map_vertex_to_tags_data_ranks)
1950  {
1951  T gIndex = key.first;
1952  const intersection_data& value = key.second;
1953  T data_on_gIndex = g2ptsData[gIndex];
1954 
1955  // find the first -1 in value
1956  int first_index = 0;
1957  for (int i = 0; i < MAX_INTERSECTIONS; ++i) {
1958  if(value.tags[i] == -1 && value.ranks[i] == -1) {
1959  first_index = i;
1960  break;
1961  }else{
1962  // check if all the intersection_data has the same g2ptsData
1963  if(value.data[i]!=data_on_gIndex)
1964  {
1965  // Error for unhandled cases.
1966  fprintf(stderr,
1967  "WARN: Due to an inconsistency in ptsData computation at gIndex %d with tag %d, "
1968  "the program is stopped while preparing the data structure to introduce DOFs.\n",
1969  gIndex, value.tags[i]);
1970  fprintf(stderr, "\n");
1971  EXIT(1);
1972  }
1973  }
1974  if (i == MAX_INTERSECTIONS -1) first_index = MAX_INTERSECTIONS;
1975  }
1976 
1977  if(data_on_gIndex==1){ // membrane
1978  // check the mumber of tag numbers,
1979  int count_intra_tags = 0;
1980  int count_extra_tags = 0;
1981  for (int i = 0; i < first_index; ++i) {
1982  if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])==extra_tags.end()){
1983  count_intra_tags += 1;
1984  }else if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])!=extra_tags.end()){
1985  count_extra_tags += 1;
1986  }
1987  }
1988 
1989  // Check if the collected data at gIndex is indeed on the membrane.
1990  // A valid membrane vertex must intersect exactly one intracellular tag and one or more extracellular tags.
1991  if(count_intra_tags>1 || count_extra_tags == 0)
1992  {
1993  // Error for unhandled cases.
1994  fprintf(stderr,
1995  "WARN: Due to an issue in defining ptsData for gIndex %d, , on the membrane"
1996  "the program is stopped while preparing the data structure to introduce DOFs.\n",
1997  gIndex);
1998  fprintf(stderr, "\n");
1999  EXIT(1);
2000  }
2001 
2002  int index_intra = -1;
2003  // Find the intracellular tag number. Since data_on_gIndex == 1, there will be only one tag corresponding to the intracellular domain.
2004  for (int i = 0; i < first_index; ++i) {
2005  if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])==extra_tags.end()){
2006  index_intra = i;
2007  break;
2008  }
2009  }
2010  if(index_intra == -1) {
2011  fprintf(stderr, "WARN: unhandled case in membrane while decoupling: gIndex %d", gIndex);
2012  fprintf(stderr, "\n");
2013  EXIT(1);
2014  }
2015  if(index_intra != -1) {
2016  T tag_myocyte = value.tags[index_intra];
2017  T rank_myocyte = value.ranks[index_intra];
2018 
2019  if(rank_myocyte == rank) { // mark intracellular tag with gIndex
2020  std::pair <mesh_int_t,mesh_int_t> Index_tag_next = std::make_pair(gIndex,tag_myocyte);
2021  if (map_mark_new_dofs.find(Index_tag_next) == map_mark_new_dofs.end()) {
2022  map_mark_new_dofs.insert({Index_tag_next,true}); //mark a new index
2023  shift++;
2024  }
2025  }
2026  }
2027  }
2028  else if(data_on_gIndex==2){ // gap junction
2029  // count the number of extra tags and intra tags
2030  int count_intra_tags = 0;
2031  int count_extra_tags = 0;
2032  for (int i = 0; i < first_index; ++i) {
2033  if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])==extra_tags.end()){
2034  count_intra_tags += 1;
2035  }else if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])!=extra_tags.end()){
2036  count_extra_tags += 1;
2037  }
2038  }
2039  if(count_intra_tags!=2 || count_extra_tags != 0){
2040  // Error for unhandled cases.
2041  fprintf(stderr,
2042  "WARN: Due to an issue in defining ptsData for gIndex %d, on the gapjunction"
2043  "the program is stopped while preparing the data structure to introduce DOFs.\n",
2044  gIndex);
2045  fprintf(stderr, "\n");
2046  EXIT(1);
2047  }
2048 
2049  T tag1 = -1;
2050  T tag2 = -1;
2051  T rank1 = -1;
2052  T rank2 = -1;
2053  bool tag1_checked = false;
2054  bool tag2_checked = false;
2055  // find the tag1 and tag2 on the gap junction
2056  for (int i = 0; i < first_index; ++i)
2057  {
2058  auto pos_extra = extra_tags.find(value.tags[i]);
2059  if(pos_extra == extra_tags.end()) {
2060  if(tag1==-1){
2061  tag1 = value.tags[i];
2062  rank1 = value.ranks[i];
2063  tag1_checked = true;
2064  }else if(tag1!=-1){
2065  tag2 = value.tags[i];
2066  rank2 = value.ranks[i];
2067  tag2_checked = true;
2068  break;
2069  }
2070  }
2071  }
2072 
2073  if(!tag1_checked && !tag2_checked || (rank1==rank2 && rank1!=rank) ) {
2074  fprintf(stderr, "WARN: unhandled case in gap junction while decoupling: gIndex %d with tags %d:%d ", gIndex, tag1, tag2);
2075  fprintf(stderr, "\n");
2076  EXIT(1);
2077  }
2078 
2079  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.
2080  T smallest_tag = tag1 < tag2? tag1:tag2;
2081  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = std::make_pair(gIndex,smallest_tag);
2082  if (map_vertex_tag_to_dof.find(Index_tag_old) == map_vertex_tag_to_dof.end() ) {
2083  map_vertex_tag_to_dof.insert({Index_tag_old,gIndex});
2084  }
2085 
2086  T bigger_tag = tag1 < tag2? tag2:tag1;
2087  std::pair <mesh_int_t,mesh_int_t> Index_tag_new = std::make_pair(gIndex,bigger_tag);
2088  if (map_mark_new_dofs.find(Index_tag_new) == map_mark_new_dofs.end() ) {
2089  map_mark_new_dofs.insert({Index_tag_new,true});
2090  shift++;
2091  }
2092  }
2093  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.
2094  if(tag1<tag2) {
2095  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = std::make_pair(gIndex,tag1);
2096  if (map_vertex_tag_to_dof.find(Index_tag_old) == map_vertex_tag_to_dof.end() ) {
2097  map_vertex_tag_to_dof.insert({Index_tag_old,gIndex});
2098  }
2099  } else {
2100  std::pair <mesh_int_t,mesh_int_t> Index_tag_new = std::make_pair(gIndex,tag1);
2101  if (map_mark_new_dofs.find(Index_tag_new) == map_mark_new_dofs.end() ) {
2102  map_mark_new_dofs.insert({Index_tag_new,true});
2103  shift++;
2104  }
2105  }
2106  }
2107  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.
2108  if(tag2<tag1) {
2109  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = std::make_pair(gIndex,tag2);
2110  if (map_vertex_tag_to_dof.find(Index_tag_old) == map_vertex_tag_to_dof.end() ) {
2111  map_vertex_tag_to_dof.insert({Index_tag_old,gIndex});
2112  }
2113  } else {
2114  std::pair <mesh_int_t,mesh_int_t> Index_tag_new = std::make_pair(gIndex,tag2);
2115  if (map_mark_new_dofs.find(Index_tag_new) == map_mark_new_dofs.end() ) {
2116  map_mark_new_dofs.insert({Index_tag_new,true});
2117  shift++;
2118  }
2119  }
2120  }
2121  }
2122  else if(data_on_gIndex==3)
2123  {
2124  // count the number of extra tags and intra tags
2125  int count_intra_tags = 0;
2126  int count_extra_tags = 0;
2127  for (int i = 0; i < first_index; ++i) {
2128  if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])==extra_tags.end()){
2129  count_intra_tags += 1;
2130  }else if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])!=extra_tags.end()){
2131  count_extra_tags += 1;
2132  }
2133  }
2134  if(count_intra_tags!=2 || count_extra_tags == 0){
2135  // Error for unhandled cases.
2136  fprintf(stderr,
2137  "WARN: Due to an issue in defining ptsData for gIndex %d, on intersection between membrane and gapjunction"
2138  "the program is stopped while preparing the data structure to introduce DOFs.\n",
2139  gIndex);
2140  fprintf(stderr, "\n");
2141  EXIT(1);
2142  }
2143 
2144  // All extracellular tags have already been added in the previous loop that iterated over the mesh.
2145  // on the interface between gap juntion and membarne
2146  for (int i = 0; i < first_index; ++i) {
2147  // find the first two tag number which belong to itracellular tags, add for indexing
2148  auto pos_extra = extra_tags.find(value.tags[i]);
2149  if(value.tags[i]!=-1 && value.ranks[i] == rank && pos_extra == extra_tags.end()) { // myocytes, for sure we should have two myocytes
2150  std::pair <mesh_int_t,mesh_int_t> Index_tag_new = std::make_pair(gIndex,value.tags[i]);
2151  if (map_mark_new_dofs.find(Index_tag_new) == map_mark_new_dofs.end() )
2152  {
2153  map_mark_new_dofs.insert({Index_tag_new,true}); //mark a new index
2154  shift++;
2155  }
2156  }
2157  }
2158  }
2159  }
2160 
2161  vector<size_t> dsp_dof(size);
2162  MPI_Allgather(&shift, sizeof(size_t), MPI_BYTE, dsp_dof.data(), sizeof(size_t), MPI_BYTE, comm);
2163 
2164  //-----------------------------------------------------------------
2165  // Step 4: Assign Unique Global Indices in Parallel
2166  //-----------------------------------------------------------------
2167  // The final step is to assign the new indices without any conflicts between processes.
2168  // - Calculate Global Offset (`start`): An MPI_Allgather is used to share the shift counts among all processes. Each process then calculates
2169  // 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
2170  // 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
2171  // of new indices created by different processes will not overlap.
2172  // - Assign New Indices: Finally, the function loops through the vertices it marked for new DOF creation and assigns them a unique global
2173  // index (newIndex = start + count). This final mapping, (original_vertex_index, tag) -> new_unique_dof_index, is inserted into the output
2174  // map map_vertex_tag_to_dof.
2175  T start;
2176  if(rank==0){
2177  start = mesh.g_numpts;
2178  }
2179  else
2180  {
2181  start = mesh.g_numpts;
2182  for (int r = 0; r < rank; ++r)
2183  {
2184  start+= dsp_dof[r];
2185  }
2186  }
2187 
2188  int count = 0;
2189  auto lexicographic_comp_pair = [](const std::pair<mesh_int_t,mesh_int_t> & a, const std::pair<mesh_int_t,mesh_int_t> & b)
2190  {
2191  if (a.second == b.second) return a.first < b.first;
2192  return a.second < b.second;
2193  };
2194 
2195  map_mark_new_dofs.sort(lexicographic_comp_pair);
2196 
2197  for(const auto & entry : map_mark_new_dofs)
2198  {
2199  T newIndex = start+count;
2200  map_vertex_tag_to_dof.insert({entry.first,newIndex});
2201  count++;
2202  }
2203 }
2204 
2211 template<class K, class intersection_indices> inline
2212 void assign_counter_dofs(hashmap::unordered_map<K, intersection_indices>& map, const MPI_Comm comm)
2213 {
2214  size_t dsize = map.size();
2215  vector<K> key_vec(dsize);
2216  vector<intersection_indices> value_vec(dsize);
2217  intersection_indices indices;
2218  size_t idx = 0;
2219  for (const auto& v : map) {
2220  key_vec[idx] = v.first;
2221  value_vec[idx] = v.second;
2222  idx++;
2223  }
2224 
2225  int size, rank;
2226  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
2227 
2228  vector<int> perm, dest(dsize);
2229  for (size_t i = 0; i < dsize; i++)
2230  dest[i] = hashmap::hash_ops<K>::hash(key_vec[i]) % size;
2231 
2232  commgraph<size_t> grph;
2233  grph.configure(dest, comm);
2234  size_t nrecv = sum(grph.rcnt);
2235 
2236  interval(perm, 0, dsize);
2237  binary_sort_copy(dest, perm);
2238 
2239  // fill send buffer and communicate
2240  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
2241  vector<intersection_indices> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
2242  for (size_t i = 0; i < dsize; i++) {
2243  sbuff_key[i] = key_vec[perm[i]];
2244  sbuff_value[i] = value_vec[perm[i]];
2245  }
2246  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
2247  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
2248 
2250  for (size_t i = 0; i < nrecv; i++)
2251  {
2252  auto it = rmap.find(rbuff_key[i]);
2253  if (it != rmap.end())
2254  {
2255  intersection_indices& map_indices = it->second;
2256  intersection_indices& r_indices = rbuff_value[i];
2257 
2258  for (int r_indices : r_indices.indices)
2259  {
2260  if (r_indices == -1) continue;
2261  bool found = false;
2262  for (int m_index : map_indices.indices)
2263  {
2264  if (m_index == r_indices) {
2265  found = true;
2266  break;
2267  }
2268  }
2269  if (!found)
2270  {
2271  for (size_t k = 0; k < MAX_INTERSECTIONS; ++k) {
2272  if (map_indices.indices[k] == -1) {
2273  map_indices.indices[k] = r_indices;
2274  break;
2275  }
2276  }
2277  }
2278  }
2279  }
2280  else
2281  {
2282  rmap.insert({ rbuff_key[i], rbuff_value[i] });
2283  }
2284  }
2285 
2286  for (size_t i = 0; i < nrecv; i++) {
2287  auto it = rmap.find(rbuff_key[i]);
2288  if (it != rmap.end()) rbuff_value[i] = it->second;
2289  }
2290 
2291  // sending the assigned data back to original rank
2292  grph.transpose();
2293  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
2294 
2295  for (size_t i = 0; i < dsize; i++) {
2296  auto it = map.find(sbuff_key[i]);
2297  if (it != map.end()) it->second = ibuff_value[i];
2298  }
2299 }
2310 template<class T, class S> inline
2311 void update_emi_mesh_with_dofs(meshdata<T,S> & mesh,
2312  const SF_nbr numbering,
2313  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
2314  mesh_int_t> & map_vertex_tag_to_dof,
2316 {
2317  // map between old index to set of new indices
2319  MPI_Comm comm = mesh.comm;
2320  mesh.globalize(numbering);
2321  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
2322  T tag = mesh.tag[eidx];
2323  T size_elem = mesh.dsp[eidx+1]-mesh.dsp[eidx];
2324 
2325  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
2326  {
2327  T gIndex = mesh.con[n];
2328 
2329  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2330  Index_tag_old = std::make_pair(gIndex,tag);
2331 
2332  auto it_new = map_vertex_tag_to_dof.find(Index_tag_old);
2333  if (it_new != map_vertex_tag_to_dof.end() )
2334  {
2335  mesh_int_t dof_new = (*it_new).second;
2336  dof2vertex.insert({dof_new, gIndex});
2337  mesh.con[n] = dof_new;
2338 
2339  auto it = vertex2dof.find(gIndex);
2340  if (it != vertex2dof.end())
2341  {
2342  intersection_indices& indices = it->second;
2343  bool found = false;
2344  for(T t : indices.indices) {
2345  if(t == dof_new) {
2346  found = true;
2347  break;
2348  }
2349  }
2350  if(!found) {
2351  // 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.
2352  for(size_t i=0; i<MAX_INTERSECTIONS; ++i) {
2353  if(indices.indices[i] == -1) {
2354  indices.indices[i] = dof_new;
2355  break;
2356  }
2357  }
2358  }
2359  }else{
2360  intersection_indices indices;
2361  indices.indices[0] = dof_new;
2362  vertex2dof.insert({gIndex,indices});
2363  }
2364  }
2365  }
2366  }
2367 
2368  assign_counter_dofs(vertex2dof, comm);
2369 
2370  // Now update dof2vertex using vertex2dof
2371  for (const auto& [old_idx, indices] : vertex2dof) {
2372  for (int new_idx : indices.indices) {
2373  if (new_idx != -1) {
2374  // Avoid overwriting unless necessary
2375  if (dof2vertex.find(new_idx) == dof2vertex.end()) {
2376  dof2vertex.insert({new_idx, old_idx});
2377  }
2378  }
2379  }
2380  }
2381 
2382  mesh.localize(numbering);
2383 }
2384 
2396 template<class T, class S> inline
2397 void update_map_indices_to_petsc(meshdata<T,S> & mesh,
2398  const SF_nbr numbering_ref,const SF_nbr numbering_petsc,
2399  hashmap::unordered_set<int> extra_tags,
2400  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
2401  std::pair<mesh_int_t,mesh_int_t>> & map_vertex_tag_to_dof_petsc,
2403  SF::vector<mesh_int_t> & elemTag_emi_mesh)
2404 {
2405  const SF::vector<SF_int> & ref_nbr = mesh.get_numbering(numbering_ref);
2406  const SF::vector<SF_int> & petsc_nbr = mesh.get_numbering(numbering_petsc);
2407 
2408  elemTag_emi_mesh.resize(mesh.l_numelem);
2409  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
2410  T tag = mesh.tag[eidx];
2411  elemTag_emi_mesh[eidx] = 1; // assigned 1 for for extracellular located on extracellular side
2412  if(extra_tags.find(tag) == extra_tags.end())
2413  elemTag_emi_mesh[eidx] = 2; // assigned 2 for for intracellular located on myocyte
2414  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
2415  {
2416  T l_Indx = mesh.con[n];
2417  T petsc_Idx = petsc_nbr[l_Indx];
2418  T n_Indx = ref_nbr[l_Indx];
2419  T o_Indx = dof2vertex[n_Indx];
2420 
2421  std::pair <mesh_int_t,mesh_int_t> oIndx_tag;
2422  oIndx_tag = std::make_pair(o_Indx,tag);
2423 
2424  auto it_new = map_vertex_tag_to_dof_petsc.find(oIndx_tag);
2425  if (it_new != map_vertex_tag_to_dof_petsc.end() )
2426  {
2427  std::pair <mesh_int_t,mesh_int_t> nIndx_petsc;
2428  nIndx_petsc = std::make_pair(n_Indx,petsc_Idx);
2429 
2430  (*it_new).second = nIndx_petsc;
2431  }
2432  }
2433  }
2434 }
2435 
2448 template<class T, class S>
2449 inline void update_faces_on_surface_mesh_after_decoupling_with_dofs(meshdata<T, S> & mesh,
2450  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
2451  mesh_int_t> & map_vertex_tag_to_dof,
2452  hashmap::unordered_map<tuple<T>,
2453  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
2454  hashmap::unordered_map<triple<T>,
2455  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
2456  hashmap::unordered_map<quadruple<T>,
2457  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face)
2458 {
2459  MPI_Comm comm = SF_COMM;
2460  int size, rank;
2461  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
2462 
2463  SF::vector<mesh_int_t> idxbuff(mesh.con);
2464  binary_sort(idxbuff); unique_resize(idxbuff);
2465 
2468  for(size_t i=0; i<idxbuff.size(); i++){
2469  g2l[idxbuff[i]] = i;
2470  l2g[i] = idxbuff[i];
2471  }
2472 
2473  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
2474  T tag = mesh.tag[eidx];
2475  T size_elem = mesh.dsp[eidx+1]-mesh.dsp[eidx];
2476  std::vector<int> elem_nodes;
2477  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
2478  {
2479  T gIndex = mesh.con[n];
2480  elem_nodes.push_back(gIndex);
2481  }
2482  std::sort(elem_nodes.begin(),elem_nodes.end());
2483  if(elem_nodes.size()==2){
2484  tuple<T> key;
2485  key.v1 = elem_nodes[0];
2486  key.v2 = elem_nodes[1];
2487  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>> value = line_face[key];
2488 
2489  line_face.erase(key);
2490  std::vector<int> new_nodes(2);
2491 
2492  auto tag_key = (value.first.rank == rank) ? value.first.tag : value.second.tag;
2493  {
2494  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2495  Index_tag_old = std::make_pair(elem_nodes[0],tag_key);
2496  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2497  new_nodes[0] = Index_new;
2498  }
2499  {
2500  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2501  Index_tag_old = std::make_pair(elem_nodes[1],tag_key);
2502  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2503  new_nodes[1] = Index_new;
2504  }
2505  std::sort(new_nodes.begin(),new_nodes.end());
2506  // update key
2507  tuple<T> new_key;
2508  new_key.v1 = new_nodes[0];
2509  new_key.v2 = new_nodes[1];
2510 
2511  // first
2512  {
2513  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v1, value.first.tag)];
2514  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v2, value.first.tag)];
2515  value.first.points.v1 = Index_new1;
2516  value.first.points.v2 = Index_new2;
2517  }
2518  // second
2519  {
2520  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v1, value.second.tag)];
2521  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v2, value.second.tag)];
2522  value.second.points.v1 = Index_new1;
2523  value.second.points.v2 = Index_new2;
2524  }
2525  line_face.insert({new_key,value});
2526  }
2527  if(elem_nodes.size()==3){
2528  triple<T> key;
2529  key.v1 = elem_nodes[0];
2530  key.v2 = elem_nodes[1];
2531  key.v3 = elem_nodes[2];
2532  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>> value = tri_face[key];
2533  tri_face.erase(key);
2534  std::vector<int> new_nodes(3);
2535 
2536  auto tag_key = (value.first.rank == rank) ? value.first.tag : value.second.tag;
2537  {
2538  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2539  Index_tag_old = std::make_pair(elem_nodes[0],tag_key);
2540  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2541  new_nodes[0] = Index_new;
2542  }
2543  {
2544  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2545  Index_tag_old = std::make_pair(elem_nodes[1],tag_key);
2546  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2547  new_nodes[1] = Index_new;
2548  }
2549  {
2550  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2551  Index_tag_old = std::make_pair(elem_nodes[2],tag_key);
2552  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2553  new_nodes[2] = Index_new;
2554  }
2555  std::sort(new_nodes.begin(),new_nodes.end());
2556  // update key
2557  triple<T> new_key;
2558  new_key.v1 = new_nodes[0];
2559  new_key.v2 = new_nodes[1];
2560  new_key.v3 = new_nodes[2];
2561 
2562  // first
2563  {
2564  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v1, value.first.tag)];
2565  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v2, value.first.tag)];
2566  mesh_int_t Index_new3 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v3, value.first.tag)];
2567  value.first.points.v1 = Index_new1;
2568  value.first.points.v2 = Index_new2;
2569  value.first.points.v3 = Index_new3;
2570  }
2571  // second
2572  {
2573  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v1, value.second.tag)];
2574  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v2, value.second.tag)];
2575  mesh_int_t Index_new3 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v3, value.second.tag)];
2576  value.second.points.v1 = Index_new1;
2577  value.second.points.v2 = Index_new2;
2578  value.second.points.v3 = Index_new3;
2579  }
2580  tri_face.insert({new_key,value});
2581  }
2582  if(elem_nodes.size()==4){
2583  quadruple<T> key;
2584  key.v1 = elem_nodes[0];
2585  key.v2 = elem_nodes[1];
2586  key.v3 = elem_nodes[2];
2587  key.v4 = elem_nodes[3];
2588  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>> value = quad_face[key];
2589  quad_face.erase(key);
2590  std::vector<int> new_nodes(4);
2591 
2592  auto tag_key = (value.first.rank == rank) ? value.first.tag : value.second.tag;
2593  {
2594  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2595  Index_tag_old = std::make_pair(elem_nodes[0],tag_key);
2596  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2597  new_nodes[0] = Index_new;
2598  }
2599  {
2600  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2601  Index_tag_old = std::make_pair(elem_nodes[1],tag_key);
2602  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2603  new_nodes[1] = Index_new;
2604  }
2605  {
2606  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2607  Index_tag_old = std::make_pair(elem_nodes[2],tag_key);
2608  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2609  new_nodes[2] = Index_new;
2610  }
2611  {
2612  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2613  Index_tag_old = std::make_pair(elem_nodes[3],tag_key);
2614  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2615  new_nodes[3] = Index_new;
2616  }
2617  std::sort(new_nodes.begin(),new_nodes.end());
2618  // update key
2619  quadruple<T> new_key;
2620  new_key.v1 = new_nodes[0];
2621  new_key.v2 = new_nodes[1];
2622  new_key.v3 = new_nodes[2];
2623  new_key.v4 = new_nodes[3];
2624 
2625  // first
2626  {
2627  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v1, value.first.tag)];
2628  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v2, value.first.tag)];
2629  mesh_int_t Index_new3 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v3, value.first.tag)];
2630  mesh_int_t Index_new4 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v4, value.first.tag)];
2631  value.first.points.v1 = Index_new1;
2632  value.first.points.v2 = Index_new2;
2633  value.first.points.v3 = Index_new3;
2634  value.first.points.v4 = Index_new4;
2635  }
2636  // second
2637  {
2638  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v1, value.second.tag)];
2639  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v2, value.second.tag)];
2640  mesh_int_t Index_new3 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v3, value.second.tag)];
2641  mesh_int_t Index_new4 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v4, value.second.tag)];
2642  value.second.points.v1 = Index_new1;
2643  value.second.points.v2 = Index_new2;
2644  value.second.points.v3 = Index_new3;
2645  value.second.points.v4 = Index_new4;
2646  }
2647  quad_face.insert({new_key,value});
2648  }
2649  }
2650 }
2651 
2663 template<class T, class S> inline
2664 void compute_surface_mesh_with_counter_face(meshdata<T, S> & mesh, const SF_nbr numbering,
2665  hashmap::unordered_map<tuple<T>,
2666  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
2667  hashmap::unordered_map<triple<T>,
2668  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
2669  hashmap::unordered_map<quadruple<T>,
2670  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face)
2671 {
2672  mesh.register_numbering(numbering);
2673 
2674  MPI_Comm comm = SF_COMM;
2675  int size, rank;
2676  MPI_Comm_size(comm, &size);
2677  MPI_Comm_rank(comm, &rank);
2678 
2679  // To ensure deterministic order, we must sort the keys of the hashmaps before iterating.
2680  vector<tuple<T>> line_keys;
2681  for(auto const& [key, val] : line_face) line_keys.push_back(key);
2682  std::sort(line_keys.begin(), line_keys.end());
2683 
2684  vector<triple<T>> tri_keys;
2685  for(auto const& [key, val] : tri_face) tri_keys.push_back(key);
2686  std::sort(tri_keys.begin(), tri_keys.end());
2687 
2688  vector<quadruple<T>> quad_keys;
2689  for(auto const& [key, val] : quad_face) quad_keys.push_back(key);
2690  std::sort(quad_keys.begin(), quad_keys.end());
2691 
2692  vector<T> cnt(mesh.l_numelem);
2693  size_t idx = 0, cidx = 0;
2694 
2695  // Add faces to the surface mesh used in ionic computation.
2696  // For each {line,tri,quad}_face, add faces only if one of the pair matches the current rank,
2697  // or if both faces have the same rank, add both.
2698  for(const auto & key : line_keys) {
2699  const auto & v = line_face.at(key);
2700  bool both_faces = (v.first.rank == v.second.rank) ? true: false;
2701 
2702  cnt[idx] = 2;
2703 
2704  emi_face<T,tuple<T>> surf_first;
2705  emi_face<T,tuple<T>> surf_second;
2706  if(v.first.rank == rank)
2707  {
2708  surf_first = v.first;
2709  surf_second = v.second;
2710  }
2711  else
2712  {
2713  surf_first = v.second;
2714  surf_second = v.first;
2715  }
2716 
2717  mesh.type[idx] = Line;
2718  mesh.tag[idx] = surf_first.tag;
2719  mesh.con[cidx + 0] = surf_first.points.v1;
2720  mesh.con[cidx + 1] = surf_first.points.v2;
2721  idx += 1;
2722  cidx += 2;
2723 
2724  if(both_faces){
2725 
2726  cnt[idx] = 2;
2727  mesh.type[idx] = Line;
2728  mesh.tag[idx] = surf_second.tag;
2729  mesh.con[cidx + 0] = surf_second.points.v1;
2730  mesh.con[cidx + 1] = surf_second.points.v2;
2731 
2732  idx += 1;
2733  cidx += 2;
2734  }
2735  }
2736 
2737  for(const auto & key : tri_keys) {
2738  const auto & v = tri_face.at(key);
2739  bool both_faces = (v.first.rank == v.second.rank) ? true: false;
2740 
2741  cnt[idx] = 3;
2742 
2743  emi_face<T,triple<T>> surf_first;
2744  emi_face<T,triple<T>> surf_second;
2745  if(v.first.rank == rank)
2746  {
2747  surf_first = v.first;
2748  surf_second = v.second;
2749  }
2750  else
2751  {
2752  surf_first = v.second;
2753  surf_second = v.first;
2754  }
2755 
2756  mesh.type[idx] = Tri;
2757  mesh.tag[idx] = surf_first.tag;
2758  mesh.con[cidx + 0] = surf_first.points.v1;
2759  mesh.con[cidx + 1] = surf_first.points.v2;
2760  mesh.con[cidx + 2] = surf_first.points.v3;
2761 
2762  idx += 1;
2763  cidx += 3;
2764 
2765  if(both_faces){
2766 
2767  cnt[idx] = 3;
2768 
2769  mesh.type[idx] = Tri;
2770  mesh.tag[idx] = surf_second.tag;
2771  mesh.con[cidx + 0] = surf_second.points.v1;
2772  mesh.con[cidx + 1] = surf_second.points.v2;
2773  mesh.con[cidx + 2] = surf_second.points.v3;
2774 
2775  idx += 1;
2776  cidx += 3;
2777  }
2778  }
2779 
2780  for(const auto & key : quad_keys) {
2781  const auto & v = quad_face.at(key);
2782  bool both_faces = (v.first.rank == v.second.rank) ? true: false;
2783  cnt[idx] = 4;
2784 
2785  emi_face<T,quadruple<T>> surf_first;
2786  emi_face<T,quadruple<T>> surf_second;
2787  if(v.first.rank == rank)
2788  {
2789  surf_first = v.first;
2790  surf_second = v.second;
2791  }
2792  else
2793  {
2794  surf_first = v.second;
2795  surf_second = v.first;
2796  }
2797 
2798  mesh.type[idx] = Quad;
2799  mesh.tag[idx] = surf_first.tag;
2800  mesh.con[cidx + 0] = surf_first.points.v1;
2801  mesh.con[cidx + 1] = surf_first.points.v2;
2802  mesh.con[cidx + 2] = surf_first.points.v3;
2803  mesh.con[cidx + 3] = surf_first.points.v4;
2804 
2805  idx += 1;
2806  cidx += 4;
2807 
2808  if(both_faces){
2809 
2810  cnt[idx] = 4;
2811 
2812  mesh.type[idx] = Quad;
2813  mesh.tag[idx] = surf_second.tag;
2814  mesh.con[cidx + 0] = surf_second.points.v1;
2815  mesh.con[cidx + 1] = surf_second.points.v2;
2816  mesh.con[cidx + 2] = surf_second.points.v3;
2817  mesh.con[cidx + 3] = surf_second.points.v4;
2818 
2819  idx += 1;
2820  cidx += 4;
2821  }
2822  }
2823  dsp_from_cnt(cnt, mesh.dsp);
2824 }
2825 
2856 template<class T, class S> inline
2857 void compute_surface_mesh_with_unique_face(meshdata<T, S> & mesh, const SF_nbr numbering,
2858  hashmap::unordered_map<tuple<T>,
2859  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
2860  hashmap::unordered_map<triple<T>,
2861  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
2862  hashmap::unordered_map<quadruple<T>,
2863  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face,
2864  hashmap::unordered_map<mesh_int_t, std::pair<mesh_int_t, mesh_int_t>> & map_elem_uniqueFace_to_tags)
2865 {
2866  mesh.register_numbering(numbering);
2867 
2868  MPI_Comm comm = SF_COMM;
2869  int size, rank;
2870  MPI_Comm_size(comm, &size);
2871  MPI_Comm_rank(comm, &rank);
2872 
2873  // To ensure deterministic order, we must sort the keys of the hashmaps before iterating.
2874  vector<tuple<T>> line_keys;
2875  for(auto const& [key, val] : line_face) line_keys.push_back(key);
2876  std::sort(line_keys.begin(), line_keys.end());
2877 
2878  vector<triple<T>> tri_keys;
2879  for(auto const& [key, val] : tri_face) tri_keys.push_back(key);
2880  std::sort(tri_keys.begin(), tri_keys.end());
2881 
2882  vector<quadruple<T>> quad_keys;
2883  for(auto const& [key, val] : quad_face) quad_keys.push_back(key);
2884  std::sort(quad_keys.begin(), quad_keys.end());
2885 
2886  vector<T> cnt(mesh.l_numelem);
2887  size_t idx = 0, cidx = 0;
2888 
2889  // Add faces to the surface mesh used in ionic computation.
2890  // For each {line,tri,quad}_face, add faces only if one of the pair matches the current rank,
2891  // or if both faces have the same rank, add both.
2892  for(const auto & key : line_keys) {
2893  const auto & v = line_face.at(key);
2894 
2895  bool to_take_face = false;
2896 
2897  emi_face<T,tuple<T>> surf_take;
2898  std::pair <mesh_int_t,mesh_int_t> tag_pairs_take;
2899  {
2900  if(v.first.rank == rank && v.first.mark_to_take == true)
2901  {
2902  surf_take = v.first;
2903  to_take_face = true;
2904  tag_pairs_take = std::make_pair(v.first.tag,v.second.tag);
2905  }
2906  else if(v.second.rank == rank && v.second.mark_to_take == true)
2907  {
2908  surf_take = v.second;
2909  to_take_face = true;
2910  tag_pairs_take = std::make_pair(v.second.tag,v.first.tag);
2911  }
2912 
2913  if(to_take_face)
2914  {
2915  map_elem_uniqueFace_to_tags.insert({idx,tag_pairs_take});
2916  cnt[idx] = 2;
2917  mesh.type[idx] = Line;
2918  mesh.tag[idx] = surf_take.tag;
2919  mesh.con[cidx + 0] = surf_take.points.v1;
2920  mesh.con[cidx + 1] = surf_take.points.v2;
2921  idx += 1;
2922  cidx += 2;
2923  }
2924  }
2925  }
2926 
2927  for(const auto & key : tri_keys) {
2928  const auto & v = tri_face.at(key);
2929 
2930  bool to_take_face = false;
2931 
2932  emi_face<T,triple<T>> surf_take;
2933  std::pair <mesh_int_t,mesh_int_t> tag_pairs_take;
2934  {
2935  if(v.first.rank == rank && v.first.mark_to_take == true)
2936  {
2937  surf_take = v.first;
2938  to_take_face = true;
2939  tag_pairs_take = std::make_pair(v.first.tag,v.second.tag);
2940  }
2941  else if(v.second.rank == rank && v.second.mark_to_take == true)
2942  {
2943  surf_take = v.second;
2944  to_take_face = true;
2945  tag_pairs_take = std::make_pair(v.second.tag,v.first.tag);
2946  }
2947 
2948  if(to_take_face)
2949  {
2950  map_elem_uniqueFace_to_tags.insert({idx,tag_pairs_take});
2951  cnt[idx] = 3;
2952  mesh.type[idx] = Tri;
2953  mesh.tag[idx] = surf_take.tag;
2954  mesh.con[cidx + 0] = surf_take.points.v1;
2955  mesh.con[cidx + 1] = surf_take.points.v2;
2956  mesh.con[cidx + 2] = surf_take.points.v3;
2957 
2958  idx += 1;
2959  cidx += 3;
2960  }
2961  }
2962  }
2963 
2964  for(const auto & key : quad_keys) {
2965  const auto & v = quad_face.at(key);
2966  cnt[idx] = 4;
2967 
2968  bool to_take_face = false;
2969 
2970  emi_face<T,quadruple<T>> surf_take;
2971  std::pair <mesh_int_t,mesh_int_t> tag_pairs_take;
2972  {
2973  if(v.first.rank == rank && v.first.mark_to_take == true)
2974  {
2975  surf_take = v.first;
2976  to_take_face = true;
2977  tag_pairs_take = std::make_pair(v.first.tag,v.second.tag);
2978  }
2979  else if(v.second.rank == rank && v.second.mark_to_take == true)
2980  {
2981  surf_take = v.second;
2982  to_take_face = true;
2983  tag_pairs_take = std::make_pair(v.second.tag,v.first.tag);
2984  }
2985 
2986  if(to_take_face)
2987  {
2988  map_elem_uniqueFace_to_tags.insert({idx,tag_pairs_take});
2989  cnt[idx] = 4;
2990  mesh.type[idx] = Quad;
2991  mesh.tag[idx] = surf_take.tag;
2992  mesh.con[cidx + 0] = surf_take.points.v1;
2993  mesh.con[cidx + 1] = surf_take.points.v2;
2994  mesh.con[cidx + 2] = surf_take.points.v3;
2995  mesh.con[cidx + 3] = surf_take.points.v4;
2996 
2997  idx += 1;
2998  cidx += 4;
2999  }
3000  }
3001  }
3002  dsp_from_cnt(cnt, mesh.dsp);
3003 }
3004 
3024 template<class T> inline
3025 void create_reverse_elem_mapping_between_surface_meshes(
3026  hashmap::unordered_map<tuple<T>, std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>>& line_face,
3027  hashmap::unordered_map<triple<T>, std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>>& tri_face,
3028  hashmap::unordered_map<quadruple<T>, std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>>& quad_face,
3029  SF::vector<T>& vec_one_to_both_face,
3030  MPI_Comm comm)
3031 {
3032  int rank;
3033  MPI_Comm_rank(comm, &rank);
3034 
3035  // To ensure deterministic order, we must sort the keys of the hashmaps before iterating.
3036  vector<tuple<T>> line_keys;
3037  for(auto const& [key, val] : line_face) line_keys.push_back(key);
3038  std::sort(line_keys.begin(), line_keys.end());
3039 
3040  vector<triple<T>> tri_keys;
3041  for(auto const& [key, val] : tri_face) tri_keys.push_back(key);
3042  std::sort(tri_keys.begin(), tri_keys.end());
3043 
3044  vector<quadruple<T>> quad_keys;
3045  for(auto const& [key, val] : quad_face) quad_keys.push_back(key);
3046  std::sort(quad_keys.begin(), quad_keys.end());
3047 
3048  size_t surf_elem_idx = 0;
3049  size_t w_counter_elem_idx = 0;
3050 
3051  size_t lsize_line = 0;
3052  size_t lsize_tri = 0;
3053  size_t lsize_quad = 0;
3054 
3055  // Count how many both-face elements are local on this rank.
3056  // If both faces are on the same rank, the both-face mesh has two entries.
3057  for(const auto& key : line_keys) {
3058  const auto& v = line_face.at(key);
3059  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3060  if (!local_involved) {
3061  continue;
3062  }
3063  if(v.first.rank == v.second.rank) lsize_line+=2;
3064  else lsize_line+=1;
3065  }
3066  for(const auto& key : tri_keys) {
3067  const auto& v = tri_face.at(key);
3068  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3069  if (!local_involved) {
3070  continue;
3071  }
3072  if(v.first.rank == v.second.rank) lsize_tri+=2;
3073  else lsize_tri+=1;
3074  }
3075  for(const auto& key : quad_keys) {
3076  const auto& v = quad_face.at(key);
3077  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3078  if (!local_involved) {
3079  continue;
3080  }
3081  if(v.first.rank == v.second.rank) lsize_quad+=2;
3082  else lsize_quad+=1;
3083  }
3084 
3085  vec_one_to_both_face.resize(lsize_line + lsize_tri + lsize_quad);
3086 
3087  // Fill mapping in deterministic order:
3088  // each both-face element gets the index of its corresponding one-face element.
3089  for (const auto& key : line_keys) {
3090  const auto& v = line_face.at(key);
3091  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3092  if (!local_involved) {
3093  continue;
3094  }
3095  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3096  if (v.first.rank == v.second.rank) {
3097  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3098  }
3099  surf_elem_idx++;
3100  }
3101 
3102  for (const auto& key : tri_keys) {
3103  const auto& v = tri_face.at(key);
3104  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3105  if (!local_involved) {
3106  continue;
3107  }
3108  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3109  if (v.first.rank == v.second.rank) {
3110  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3111  }
3112  surf_elem_idx++;
3113  }
3114 
3115  for (const auto& key : quad_keys) {
3116  const auto& v = quad_face.at(key);
3117  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3118  if (!local_involved) {
3119  continue;
3120  }
3121  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3122  if (v.first.rank == v.second.rank) {
3123  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3124  }
3125  surf_elem_idx++;
3126  }
3127 }
3128 
3137 template<class T> inline
3138 void added_counter_faces_to_map(hashmap::unordered_map<tuple<T>,
3139  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
3140  hashmap::unordered_map<triple<T>,
3141  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
3142  hashmap::unordered_map<quadruple<T>,
3143  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face)
3144 {
3145  for(const auto & v : line_face) {
3146  SF::emi_face<mesh_int_t,SF::tuple<mesh_int_t>> surf_first = v.second.first;
3147  SF::emi_face<mesh_int_t,SF::tuple<mesh_int_t>> surf_second = v.second.second;
3148 
3149  std::vector<mesh_int_t> elem_nodes_first; elem_nodes_first.resize(2);
3150  std::vector<mesh_int_t> elem_nodes_second; elem_nodes_second.resize(2);
3151 
3152  SF::tuple<mesh_int_t> key_first;
3153  elem_nodes_first[0] = surf_first.points.v1;
3154  elem_nodes_first[1] = surf_first.points.v2;
3155  std::sort(elem_nodes_first.begin(),elem_nodes_first.end());
3156  key_first.v1 = elem_nodes_first[0];
3157  key_first.v2 = elem_nodes_first[1];
3158 
3159  SF::tuple<mesh_int_t> key_second;
3160  elem_nodes_second[0] = surf_second.points.v1;
3161  elem_nodes_second[1] = surf_second.points.v2;
3162  std::sort(elem_nodes_second.begin(),elem_nodes_second.end());
3163  key_second.v1 = elem_nodes_second[0];
3164  key_second.v2 = elem_nodes_second[1];
3165 
3166  if(v.first.v1!=key_first.v1 and v.first.v2!=key_first.v2){
3167  line_face.insert({key_first,v.second});
3168  }
3169  if(v.first.v1!=key_second.v1 and v.first.v2!=key_second.v2){
3170  line_face.insert({key_second,v.second});
3171 
3172  }
3173  }
3174 
3175  for(const auto & v : tri_face) {
3176  SF::emi_face<mesh_int_t,SF::triple<mesh_int_t>> surf_first = v.second.first;
3177  SF::emi_face<mesh_int_t,SF::triple<mesh_int_t>> surf_second = v.second.second;
3178 
3179  std::vector<mesh_int_t> elem_nodes_first; elem_nodes_first.resize(3);
3180  std::vector<mesh_int_t> elem_nodes_second; elem_nodes_second.resize(3);
3181 
3182  SF::triple<mesh_int_t> key_first;
3183  elem_nodes_first[0] = surf_first.points.v1;
3184  elem_nodes_first[1] = surf_first.points.v2;
3185  elem_nodes_first[2] = surf_first.points.v3;
3186  std::sort(elem_nodes_first.begin(),elem_nodes_first.end());
3187  key_first.v1 = elem_nodes_first[0];
3188  key_first.v2 = elem_nodes_first[1];
3189  key_first.v3 = elem_nodes_first[2];
3190 
3191  SF::triple<mesh_int_t> key_second;
3192  elem_nodes_second[0] = surf_second.points.v1;
3193  elem_nodes_second[1] = surf_second.points.v2;
3194  elem_nodes_second[2] = surf_second.points.v3;
3195  std::sort(elem_nodes_second.begin(),elem_nodes_second.end());
3196  key_second.v1 = elem_nodes_second[0];
3197  key_second.v2 = elem_nodes_second[1];
3198  key_second.v3 = elem_nodes_second[2];
3199 
3200  if(v.first.v1!=key_first.v1 and v.first.v2!=key_first.v2 and v.first.v3!=key_first.v3){
3201  tri_face.insert({key_first,v.second});
3202  }
3203  if(v.first.v1!=key_second.v1 and v.first.v2!=key_second.v2 and v.first.v3!=key_second.v3){
3204  tri_face.insert({key_second,v.second});
3205  }
3206  }
3207 
3208  for(const auto & v : quad_face) {
3209  SF::emi_face<mesh_int_t,SF::quadruple<mesh_int_t>> surf_first = v.second.first;
3210  SF::emi_face<mesh_int_t,SF::quadruple<mesh_int_t>> surf_second = v.second.second;
3211 
3212  std::vector<mesh_int_t> elem_nodes_first; elem_nodes_first.resize(4);
3213  std::vector<mesh_int_t> elem_nodes_second; elem_nodes_second.resize(4);
3214 
3215  SF::quadruple<mesh_int_t> key_first;
3216  elem_nodes_first[0] = surf_first.points.v1;
3217  elem_nodes_first[1] = surf_first.points.v2;
3218  elem_nodes_first[2] = surf_first.points.v3;
3219  elem_nodes_first[3] = surf_first.points.v4;
3220  std::sort(elem_nodes_first.begin(),elem_nodes_first.end());
3221  key_first.v1 = elem_nodes_first[0];
3222  key_first.v2 = elem_nodes_first[1];
3223  key_first.v3 = elem_nodes_first[2];
3224  key_first.v4 = elem_nodes_first[3];
3225 
3226  SF::quadruple<mesh_int_t> key_second;
3227  elem_nodes_second[0] = surf_second.points.v1;
3228  elem_nodes_second[1] = surf_second.points.v2;
3229  elem_nodes_second[2] = surf_second.points.v3;
3230  elem_nodes_second[3] = surf_second.points.v4;
3231  std::sort(elem_nodes_second.begin(),elem_nodes_second.end());
3232  key_second.v1 = elem_nodes_second[0];
3233  key_second.v2 = elem_nodes_second[1];
3234  key_second.v3 = elem_nodes_second[2];
3235  key_second.v4 = elem_nodes_second[3];
3236 
3237  if(v.first.v1!=key_first.v1 and v.first.v2!=key_first.v2 and v.first.v3!=key_first.v3 and v.first.v4!=key_first.v4){
3238  quad_face.insert({key_first,v.second});
3239  }
3240  if(v.first.v1!=key_second.v1 and v.first.v2!=key_second.v2 and v.first.v3!=key_second.v3 and v.first.v4!=key_second.v4){
3241  quad_face.insert({key_second,v.second});
3242  }
3243  }
3244 }
3245 
3252 template<class K> inline
3253 void assign_counter_tags(hashmap::unordered_map<K, intersection_tags>& map, const MPI_Comm comm)
3254 {
3255  size_t dsize = map.size();
3256  vector<K> key_vec(dsize);
3257  vector<intersection_tags> value_vec(dsize);
3258 
3259  size_t idx = 0;
3260  for (const auto& v : map) {
3261  key_vec[idx] = v.first;
3262  value_vec[idx] = v.second;
3263  idx++;
3264  }
3265 
3266  int size, rank;
3267  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
3268 
3269  vector<int> perm, dest(dsize);
3270  for (size_t i = 0; i < dsize; i++)
3271  dest[i] = hashmap::hash_ops<K>::hash(key_vec[i]) % size;
3272 
3273  commgraph<size_t> grph;
3274  grph.configure(dest, comm);
3275  size_t nrecv = sum(grph.rcnt);
3276 
3277  interval(perm, 0, dsize);
3278  binary_sort_copy(dest, perm);
3279 
3280  // fill send buffer and communicate
3281  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
3282  vector<intersection_tags> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
3283  for (size_t i = 0; i < dsize; i++) {
3284  sbuff_key[i] = key_vec[perm[i]];
3285  sbuff_value[i] = value_vec[perm[i]];
3286  }
3287  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
3288  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
3289 
3291  for (size_t i = 0; i < nrecv; i++)
3292  {
3293  auto it = rmap.find(rbuff_key[i]);
3294  if (it != rmap.end())
3295  {
3296  intersection_tags& map_tags = it->second;
3297  intersection_tags& r_tags = rbuff_value[i];
3298 
3299  for (int r_tag : r_tags.tags)
3300  {
3301  if (r_tag == -1) continue;
3302  bool found = false;
3303  for (int m_tag : map_tags.tags)
3304  {
3305  if (m_tag == r_tag) {
3306  found = true;
3307  break;
3308  }
3309  }
3310  if (!found)
3311  {
3312  for (size_t k = 0; k < MAX_INTERSECTIONS; ++k) {
3313  if (map_tags.tags[k] == -1) {
3314  map_tags.tags[k] = r_tag;
3315  break;
3316  }
3317  }
3318  }
3319  }
3320  }
3321  else
3322  {
3323  rmap.insert({ rbuff_key[i], rbuff_value[i] });
3324  }
3325  }
3326 
3327  for (size_t i = 0; i < nrecv; i++) {
3328  auto it = rmap.find(rbuff_key[i]);
3329  if (it != rmap.end()) rbuff_value[i] = it->second;
3330  }
3331 
3332  // sending the assigned data back to original rank
3333  grph.transpose();
3334  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
3335 
3336  for (size_t i = 0; i < dsize; i++) {
3337  auto it = map.find(sbuff_key[i]);
3338  if (it != map.end()) it->second = ibuff_value[i];
3339  }
3340 }
3341 
3363 template<class T, class S> inline
3364 void compute_ptsdata_from_original_mesh(meshdata<T,S> & mesh,
3365  const SF_nbr numbering,
3366  hashmap::unordered_map<T,T> & vertex2ptsdata,
3367  hashmap::unordered_set<int> extra_tags,
3368  hashmap::unordered_set<int> intra_tags)
3369 {
3370  MPI_Comm comm = mesh.comm;
3371  int size, rank;
3372  MPI_Comm_size(comm, &size);
3373  MPI_Comm_rank(comm, &rank);
3374 
3375  const T* con = mesh.con.data();
3376  const SF::vector<T> & rnod = mesh.get_numbering(numbering);
3377 
3378  // Step 1: Gather all unique region tags associated with each vertex.
3379  // This map is local to each process and will store a list of tags for each vertex
3380  // based on the elements owned by the current process.
3382  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++)
3383  {
3384  T tag = mesh.tag[eidx];
3385  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
3386  {
3387  T gIndex = rnod[con[n]];
3388 
3389  auto it_new = map_index_to_tags.find(gIndex);
3390  if (it_new != map_index_to_tags.end())
3391  {
3392  // If vertex is already in the map, add the new tag if it's not already present.
3393  intersection_tags& tags = it_new->second;
3394  bool found = false;
3395  for(T t : tags.tags) {
3396  if(t == tag) {
3397  found = true;
3398  break;
3399  }
3400  }
3401  if(!found) {
3402  // This is an intracellular tag. Find an empty slot and add it.
3403  for(int i=0; i <MAX_INTERSECTIONS; ++i) {
3404  if (tags.tags[i] == -1) {
3405  tags.tags[i] = tag;
3406  break;
3407  }
3408  }
3409  }
3410  }
3411  else{
3412  // If vertex is not in the map, add it with the current tag.
3413  intersection_tags tags;
3414  tags.tags[0] = tag;
3415  map_index_to_tags.insert({gIndex,tags});
3416  }
3417  }
3418  }
3419 
3420  // Step 2: Globalize the tag information.
3421  // After this call, `map_index_to_tags` on every process will contain the
3422  // complete list of unique tags for every vertex in the entire mesh.
3423  assign_counter_tags(map_index_to_tags, comm);
3424 
3425  // Step 3: Classify each vertex based on its global list of tags.
3426  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++)
3427  {
3428  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
3429  {
3430  T gIndex = rnod[con[n]];
3431  intersection_tags& vec_tags = map_index_to_tags[gIndex];
3432  int count_intra_tags = 0;
3433  int count_extra_tags = 0;
3434 
3435  // Count the number of unique intracellular and extracellular tags for the vertex.
3436  int count_tags = 0;
3437  for(T t : vec_tags.tags) {
3438  if(intra_tags.count(t)) count_intra_tags++;
3439  if(extra_tags.count(t)) count_extra_tags++;
3440  }
3441 
3442  // Apply classification rules based on the tag counts.
3443  if(count_extra_tags>=1 && count_intra_tags==0){
3444  vertex2ptsdata[gIndex] = 0; // Interior to an extracellular region
3445  }else if(count_extra_tags==0 && count_intra_tags==1){
3446  vertex2ptsdata[gIndex] = 0; // Interior to a myocyte
3447  }else if(count_extra_tags>=1 && count_intra_tags==1){
3448  vertex2ptsdata[gIndex] = 1; // On a membrane (1 myocyte, 1 extracellular)
3449  }else if(count_extra_tags==0 && count_intra_tags==2){
3450  vertex2ptsdata[gIndex] = 2; // On a gap junction (2 myocytes)
3451  }else if(count_extra_tags>=1 && count_intra_tags==2){
3452  vertex2ptsdata[gIndex] = 3; // On a complex junction (2 myocytes, 1 extracellular)
3453  }
3454  else if(count_intra_tags>2){
3455  // Error for unhandled cases.
3456  fprintf(stderr, "More than two intracellular are connected %d. Tags:", gIndex);
3457  for(T tag : vec_tags.tags) if(tag != -1) fprintf(stderr, " %d", tag);
3458  fprintf(stderr, "\n");
3459  EXIT(1);
3460  }
3461  else{
3462  // Error for unhandled cases.
3463  fprintf(stderr, "WARN: unhandled case in ptsData computation for gIndex %d. Tags:", gIndex);
3464  for(T tag : vec_tags.tags) if(tag != -1) fprintf(stderr, " %d", tag);
3465  fprintf(stderr, "\n");
3466  EXIT(1);
3467  }
3468  }
3469  }
3470 }
3471 
3472 
3480 template <class T, class S, class V, class emi_index_rank>
3481 inline void assemble_map_both_to_unique(SF::abstract_matrix<T, S>& op,
3483  std::pair<emi_index_rank, emi_index_rank>>& map,
3484  const SF::meshdata<T, V>& unique_mesh,
3485  const SF::meshdata<T, V>& both_mesh)
3486 {
3487  int rank;
3488  MPI_Comm_rank(both_mesh.comm, &rank);
3489 
3490  SF::vector<long int> layout_both;
3491  SF::layout_from_count<long int>(both_mesh.l_numelem, layout_both, both_mesh.comm);
3492 
3493  SF::vector<long int> layout_unique;
3494  SF::layout_from_count<long int>(unique_mesh.l_numelem, layout_unique, unique_mesh.comm);
3495 
3496  SF::vector<SF_int> row_idx(1), col_idx(1);
3497  SF::dmat<SF_real> ebuff(1, 1);
3498 
3499  // Assemble: For each unique face that this rank owns, average the two both-face values
3500  for (const auto& [local_unique_idx, both_pair] : map) {
3501  const auto& first_both = both_pair.first;
3502  const auto& second_both = both_pair.second;
3503 
3504  // Get global row index for this unique face
3505  T global_unique_row = layout_unique[rank] + local_unique_idx;
3506  row_idx[0] = global_unique_row;
3507 
3508  // Determine valid global columns
3509  bool first_valid = (first_both.index >= 0 && first_both.rank >= 0);
3510  bool second_valid = (second_both.index >= 0 && second_both.rank >= 0);
3511 
3512  T global_both_col_first = -1;
3513  T global_both_col_second = -1;
3514  if (first_valid) {
3515  global_both_col_first = layout_both[first_both.rank] + first_both.index;
3516  }
3517  if (second_valid) {
3518  global_both_col_second = layout_both[second_both.rank] + second_both.index;
3519  }
3520 
3521  // If both map to the same column, treat as single-sided
3522  if (first_valid && second_valid && global_both_col_first == global_both_col_second) {
3523  second_valid = false;
3524  }
3525 
3526  // Count how many both-face entries we have (1 or 2)
3527  int count = 0;
3528  if (first_valid) count++;
3529  if (second_valid) count++;
3530 
3531  if (count == 0) continue; // No data to average
3532 
3533  double weight = 0.5; // Average the two sides
3534  if (count == 1) weight = 1.0; // Only one side available
3535 
3536  ebuff.assign(1, 1, weight);
3537 
3538  // Set first both-face entry
3539  if (first_valid) {
3540  col_idx[0] = global_both_col_first;
3541  op.set_values(row_idx, col_idx, ebuff.data(), false);
3542  }
3543 
3544  // Set second both-face entry
3545  if (second_valid) {
3546  col_idx[0] = global_both_col_second;
3547  op.set_values(row_idx, col_idx, ebuff.data(), false);
3548  }
3549  }
3550 
3551  op.finish_assembly();
3552 }
3553 
3570 template<class T>
3571 inline void restrict_to_membrane(vector<T> & v, hashmap::unordered_map<T,T> & dof2ptsData, const meshdata<int, float> & mesh)
3572 {
3573  const SF::vector<mesh_int_t> & rnod = mesh.get_numbering(SF::NBR_REF);
3576  for(size_t i=0; i<rnod.size(); i++){
3577  g2l[rnod[i]] = i;
3578  l2g[i] = rnod[i];
3579  }
3580 
3581  size_t widx = 0;
3582 
3583  for(size_t i=0; i<v.size(); i++){
3584  T g = l2g[v[i]];
3585 
3586  if (dof2ptsData[g] > 0) { // dof2ptsData requires global index
3587  v[widx++] = v[i];
3588  }
3589  }
3590 
3591  v.resize(widx);
3592 }
3593 
3594 }
3595 
3596 #endif
3597 #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:593
void sort(Compare comp=Compare())
Sort data entries.
Definition: hashmap.hpp:659
void insert(InputIterator first, InputIterator last)
Insert Iterator range.
Definition: hashmap.hpp:539
T & at(const K &key)
Data access by key.
Definition: hashmap.hpp:613
size_t size() const
Definition: hashmap.hpp:687
iterator find(const K &key)
Definition: hashmap.hpp:1006
hm_int count(const K &key) const
Definition: hashmap.hpp:992
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
@ Line
Definition: SF_container.h:61
@ 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
Base hashing class.
Definition: hashmap.hpp:88
static hm_uint hash(const T &a)
Definition: hashmap.hpp:92