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 
1287 template<class T, class S> inline
1288 void extract_face_based_tags(meshdata<T,S> & mesh,
1289  const SF_nbr numbering,
1290  hashmap::unordered_map<T,T> & vertex2ptsdata,
1291  hashmap::unordered_map<tuple<T> ,
1292  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
1293  hashmap::unordered_map<triple<T>,
1294  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
1295  hashmap::unordered_map<quadruple<T>,
1296  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face,
1297  hashmap::unordered_set<int> & extra_tags,
1298  hashmap::unordered_set<int> & intra_tags,
1299  meshdata<T,S> & surfmesh,
1300  meshdata<T,S> & surfmesh_w_counter,
1301  meshdata<T,S> & surfmesh_unique_face,
1302  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,
1303  hashmap::unordered_map<mesh_int_t, emi_index_rank<mesh_int_t>> & map_elem_oneface_to_elem_uniqueFace)
1304 {
1305  // ============================================================================
1306  // STEP 1: Initialize surface meshes and extract interface faces
1307  // ============================================================================
1308  surfmesh.register_numbering(SF::NBR_REF);
1309  surfmesh_w_counter.register_numbering(SF::NBR_REF);
1310  surfmesh_unique_face.register_numbering(SF::NBR_REF);
1311 
1312  // Extract faces with counter faces on the interfaces to build surface mesh
1313 
1314  compute_surface_with_tags(mesh, numbering, vertex2ptsdata, extra_tags, line_face, tri_face, quad_face);
1315 
1316  int size, rank;
1317  MPI_Comm_size(mesh.comm, &size);
1318  MPI_Comm_rank(mesh.comm, &rank);
1319 
1320  long int g_num_line = line_face.size();
1321  long int g_num_tri = tri_face.size();
1322  long int g_num_quad = quad_face.size();
1323  MPI_Allreduce(MPI_IN_PLACE, &g_num_line, 1, MPI_LONG, MPI_SUM, mesh.comm);
1324  MPI_Allreduce(MPI_IN_PLACE, &g_num_tri, 1, MPI_LONG, MPI_SUM, mesh.comm);
1325  MPI_Allreduce(MPI_IN_PLACE, &g_num_quad, 1, MPI_LONG, MPI_SUM, mesh.comm);
1326 
1327  // Warn if any extra-extra faces remain (should be none after filtering)
1328  {
1329  long int local_extra_extra = 0;
1330  for (const auto& [key, v] : line_face) {
1331  if (extra_tags.find(v.first.tag) != extra_tags.end() &&
1332  extra_tags.find(v.second.tag) != extra_tags.end()) {
1333  local_extra_extra++;
1334  }
1335  }
1336  for (const auto& [key, v] : tri_face) {
1337  if (extra_tags.find(v.first.tag) != extra_tags.end() &&
1338  extra_tags.find(v.second.tag) != extra_tags.end()) {
1339  local_extra_extra++;
1340  }
1341  }
1342  for (const auto& [key, v] : quad_face) {
1343  if (extra_tags.find(v.first.tag) != extra_tags.end() &&
1344  extra_tags.find(v.second.tag) != extra_tags.end()) {
1345  local_extra_extra++;
1346  }
1347  }
1348  long int global_extra_extra = 0;
1349  MPI_Allreduce(&local_extra_extra, &global_extra_extra, 1, MPI_LONG, MPI_SUM, mesh.comm);
1350  if (global_extra_extra > 0 && rank == 0) {
1351  fprintf(stderr, "WARN: extra-extra faces detected before filtering: %ld\n", global_extra_extra);
1352  }
1353  }
1354 
1355  surfmesh.g_numelem = g_num_line + g_num_tri + g_num_quad;
1356  surfmesh.l_numelem = line_face.size() + tri_face.size() + quad_face.size();
1357 
1358  vector<T> cnt(surfmesh.l_numelem, 0);
1359  surfmesh.tag.resize(surfmesh.l_numelem);
1360  surfmesh.type.resize(surfmesh.l_numelem);
1361  surfmesh.con.resize(line_face.size() * 2 + tri_face.size() * 3 + quad_face.size() * 4);
1362 
1363  // ============================================================================
1364  // STEP 2: Sort face keys for deterministic ordering across MPI ranks
1365  // ============================================================================
1366  // To ensure deterministic order, we must sort the keys of the hashmaps before iterating.
1367  // This is crucial for MPI consistency - all ranks must process faces in the same order.
1368  vector<tuple<T>> line_keys;
1369  for(auto const& [key, val] : line_face) line_keys.push_back(key);
1370  std::sort(line_keys.begin(), line_keys.end());
1371 
1372  vector<triple<T>> tri_keys;
1373  for(auto const& [key, val] : tri_face) tri_keys.push_back(key);
1374  std::sort(tri_keys.begin(), tri_keys.end());
1375 
1376  vector<quadruple<T>> quad_keys;
1377  for(auto const& [key, val] : quad_face) quad_keys.push_back(key);
1378  std::sort(quad_keys.begin(), quad_keys.end());
1379 
1380  // Temporary hashmaps to track unique-face assignments before MPI communication
1381  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;
1382  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;
1383  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;
1384 
1385  // Counters for local element indices:
1386  // - idx: index into surfmesh (one-side mesh)
1387  // - lsize_*: count of elements in surfmesh_w_counter
1388  // - lsize_unique_*: count of elements in surfmesh_unique_face
1389  mesh_int_t idx = 0, cidx = 0;
1390  mesh_int_t lsize_line = 0;
1391  mesh_int_t lsize_tri = 0;
1392  mesh_int_t lsize_quad = 0;
1393 
1394  mesh_int_t lsize_unique_line = 0;
1395  mesh_int_t lsize_unique_tri = 0;
1396  mesh_int_t lsize_unique_quad = 0;
1397 
1398  // ============================================================================
1399  // STEP 3: Process each face and determine unique-face ownership
1400  // ============================================================================
1401  // For each interface face, we:
1402  // 1. Add it to surfmesh
1403  for(auto const& key : line_keys) {
1404  auto& v = line_face.at(key);
1405  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
1406  if (!local_involved) {
1407  continue;
1408  }
1409  cnt[idx] = 2;
1410 
1411  // Mark gap-junction faces explicitly; membrane is the default.
1412  if ((intra_tags.find(v.first.tag) != intra_tags.end()) && (intra_tags.find(v.second.tag) != intra_tags.end())) {
1413  v.first.mem = 2; // to set the face located on gap-junction
1414  v.second.mem = 2; // to set the face located on gap-junction
1415  }
1416 
1417  // changing the strategy to select the face with the one with the same rank
1418  emi_face<T,tuple<T>> surf_neighbor =
1419  (v.first.rank == rank) ? v.first : v.second;
1420 
1421  surfmesh.type[idx] = Line;
1422  surfmesh.tag[idx] = surf_neighbor.tag;
1423  surfmesh.con[cidx + 0] = surf_neighbor.points.v1;
1424  surfmesh.con[cidx + 1] = surf_neighbor.points.v2;
1425 
1426  // If the ranks are equal, both faces will be added later.
1427  if(v.first.rank == v.second.rank)
1428  lsize_line+=2;
1429  else
1430  lsize_line+=1;
1431 
1432  assign_ownership_rank_on_faces<T, tuple<T>, tuple<T>>( key, v, rank,
1433  idx,
1434  lsize_unique_line,
1435  extra_tags, intra_tags,
1436  line_unique_face_to_elements,
1437  map_elem_oneface_to_elem_uniqueFace);
1438  idx += 1;
1439  cidx += 2;
1440  }
1441 
1442  for(auto const& key : tri_keys) {
1443  auto& v = tri_face.at(key);
1444  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
1445  if (!local_involved) {
1446  continue;
1447  }
1448  cnt[idx] = 3;
1449 
1450  if ((intra_tags.find(v.first.tag) != intra_tags.end()) && (intra_tags.find(v.second.tag) != intra_tags.end())) {
1451  v.first.mem = 2; // to set the face located on gap-junction
1452  v.second.mem = 2; // to set the face located on gap-junction
1453  }
1454 
1455  // changing the strategy to select the face with the one with the same rank
1456  emi_face<T,triple<T>> surf_neighbor =
1457  (v.first.rank == rank) ? v.first : v.second;
1458 
1459  surfmesh.type[idx] = Tri;
1460  surfmesh.tag[idx] = surf_neighbor.tag;
1461  surfmesh.con[cidx + 0] = surf_neighbor.points.v1;
1462  surfmesh.con[cidx + 1] = surf_neighbor.points.v2;
1463  surfmesh.con[cidx + 2] = surf_neighbor.points.v3;
1464 
1465  // If the ranks are equal, both faces will be added later.
1466  if(v.first.rank == v.second.rank)
1467  lsize_tri+=2;
1468  else
1469  lsize_tri+=1;
1470 
1471  assign_ownership_rank_on_faces<T, triple<T>, triple<T>>( key, v, rank,
1472  idx,
1473  lsize_unique_tri,
1474  extra_tags, intra_tags,
1475  tri_unique_face_to_elements,
1476  map_elem_oneface_to_elem_uniqueFace);
1477  idx += 1;
1478  cidx += 3;
1479  }
1480 
1481  for(auto const& key : quad_keys) {
1482  auto& v = quad_face.at(key);
1483  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
1484  if (!local_involved) {
1485  continue;
1486  }
1487  cnt[idx] = 4;
1488 
1489  if ((intra_tags.find(v.first.tag) != intra_tags.end()) && (intra_tags.find(v.second.tag) != intra_tags.end())) {
1490  v.first.mem = 2; // to set the face located on gap-junction
1491  v.second.mem = 2; // to set the face located on gap-junction
1492  }
1493 
1494  // changing the strategy to select the face with the one with the same rank
1495  emi_face<T,quadruple<T>> qsurf_neighbor =
1496  (v.first.rank == rank) ? v.first : v.second;
1497 
1498  surfmesh.type[idx] = Quad;
1499  surfmesh.tag[idx] = qsurf_neighbor.tag;
1500  surfmesh.con[cidx + 0] = qsurf_neighbor.points.v1;
1501  surfmesh.con[cidx + 1] = qsurf_neighbor.points.v2;
1502  surfmesh.con[cidx + 2] = qsurf_neighbor.points.v3;
1503  surfmesh.con[cidx + 3] = qsurf_neighbor.points.v4;
1504 
1505  // If the ranks are equal, both faces will be added later.
1506  if(v.first.rank == v.second.rank)
1507  lsize_quad+=2;
1508  else
1509  lsize_quad+=1;
1510 
1511  assign_ownership_rank_on_faces<T, quadruple<T>, quadruple<T>>( key, v, rank,
1512  idx,
1513  lsize_unique_quad,
1514  extra_tags, intra_tags,
1515  quad_unique_face_to_elements,
1516  map_elem_oneface_to_elem_uniqueFace);
1517  idx += 1;
1518  cidx += 4;
1519  }
1520 
1521  surfmesh.l_numelem = idx;
1522  surfmesh.tag.resize(idx);
1523  surfmesh.type.resize(idx);
1524  cnt.resize(idx);
1525  dsp_from_cnt(cnt, surfmesh.dsp);
1526  surfmesh.con.resize(cidx);
1527 
1528  long int g_num_surf = surfmesh.l_numelem;
1529  MPI_Allreduce(MPI_IN_PLACE, &g_num_surf, 1, MPI_LONG, MPI_SUM, mesh.comm);
1530  surfmesh.g_numelem = g_num_surf;
1531 
1532  if(rank ==0){
1533  std::cout << "surfmesh.g_numelem: " << surfmesh.g_numelem <<std::endl;
1534  std::cout << "surfmesh.l_numelem: " << surfmesh.l_numelem <<std::endl;
1535  }
1536 
1537  {
1538  long int g_num_w_counter_line = lsize_line;
1539  long int g_num_w_counter_tri = lsize_tri;
1540  long int g_num_w_counter_quad = lsize_quad;
1541 
1542  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_counter_line, 1, MPI_LONG, MPI_SUM, SF_COMM);
1543  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_counter_tri, 1, MPI_LONG, MPI_SUM, SF_COMM);
1544  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_counter_quad, 1, MPI_LONG, MPI_SUM, SF_COMM);
1545 
1546  // Assign global and local elements for the surface mesh used in ionic computation.
1547  surfmesh_w_counter.g_numelem = g_num_w_counter_line+g_num_w_counter_tri+g_num_w_counter_quad;
1548  surfmesh_w_counter.l_numelem = lsize_line + lsize_tri + lsize_quad;
1549 
1550  surfmesh_w_counter.tag.resize(surfmesh_w_counter.l_numelem);
1551  surfmesh_w_counter.type.resize(surfmesh_w_counter.l_numelem);
1552  surfmesh_w_counter.con.resize(lsize_line * 2 + lsize_tri * 3 + lsize_quad * 4);
1553 
1554  if(rank ==0){
1555  std::cout << "surfmesh_w_counter.g_numelem: " << surfmesh_w_counter.g_numelem <<std::endl;
1556  std::cout << "surfmesh_w_counter.l_numelem: " << surfmesh_w_counter.l_numelem <<std::endl;
1557  }
1558  }
1559 
1560  {
1561  long int g_num_w_unique_line = lsize_unique_line;
1562  long int g_num_w_unique_tri = lsize_unique_tri;
1563  long int g_num_w_unique_quad = lsize_unique_quad;
1564 
1565  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_unique_line, 1, MPI_LONG, MPI_SUM, SF_COMM);
1566  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_unique_tri, 1, MPI_LONG, MPI_SUM, SF_COMM);
1567  MPI_Allreduce(MPI_IN_PLACE, &g_num_w_unique_quad, 1, MPI_LONG, MPI_SUM, SF_COMM);
1568 
1569  // Assign global and local elements for the surface mesh used in ionic computation.
1570  surfmesh_unique_face.g_numelem = g_num_w_unique_line+g_num_w_unique_tri+g_num_w_unique_quad;
1571  surfmesh_unique_face.l_numelem = lsize_unique_line + lsize_unique_tri + lsize_unique_quad;
1572 
1573  surfmesh_unique_face.tag.resize(surfmesh_unique_face.l_numelem);
1574  surfmesh_unique_face.type.resize(surfmesh_unique_face.l_numelem);
1575  surfmesh_unique_face.con.resize(lsize_unique_line * 2 + lsize_unique_tri * 3 + lsize_unique_quad * 4);
1576 
1577  if(rank ==0){
1578  std::cout << "surfmesh_unique_face.g_numelem: " << surfmesh_unique_face.g_numelem <<std::endl;
1579  std::cout << "surfmesh_unique_face.l_numelem: " << surfmesh_unique_face.l_numelem <<std::endl;
1580  }
1581  }
1582 
1583  // ============================================================================
1584  // STEP 4: MPI communication to share unique-face ownership between ranks
1585  // ============================================================================
1586  // When faces span multiple ranks, only one rank "owns" the unique-face.
1587  // assign_unique_first_face communicates the owner's index_unique to non-owner ranks
1588  // so they can set up their map_elem_oneface_to_elem_uniqueFace mappings correctly.
1589 
1590  // DEBUG: Before exchange
1591  #ifdef EMI_DEBUG_MESH
1592  {
1593  fprintf(stderr, "RANK %d BEFORE exchange: line_hash=%zu, tri_hash=%zu, quad_hash=%zu (expected unique total=%ld)\n",
1594  rank, line_unique_face_to_elements.size(), tri_unique_face_to_elements.size(),
1595  quad_unique_face_to_elements.size(), lsize_unique_line + lsize_unique_tri + lsize_unique_quad);
1596  fflush(stderr);
1597  }
1598  #endif
1599 
1600  if(size>1)
1601  {
1602  // Find the counter faces corresponding to the selected potential faces.
1603  assign_unique_first_face(line_unique_face_to_elements, mesh.comm);
1604  assign_unique_first_face(tri_unique_face_to_elements, mesh.comm);
1605  assign_unique_first_face(quad_unique_face_to_elements, mesh.comm);
1606  }
1607 
1608  // DEBUG: After exchange
1609  #ifdef EMI_DEBUG_MESH
1610  {
1611  int valid_line = 0, valid_tri = 0, valid_quad = 0;
1612  for (const auto& [key, val] : line_unique_face_to_elements) {
1613  if (val.first.index_unique >= 0 || val.second.index_unique >= 0) valid_line++;
1614  }
1615  for (const auto& [key, val] : tri_unique_face_to_elements) {
1616  if (val.first.index_unique >= 0 || val.second.index_unique >= 0) valid_tri++;
1617  }
1618  for (const auto& [key, val] : quad_unique_face_to_elements) {
1619  if (val.first.index_unique >= 0 || val.second.index_unique >= 0) valid_quad++;
1620  }
1621  fprintf(stderr, "RANK %d AFTER exchange: line_hash=%zu, tri_hash=%zu, quad_hash=%zu\n",
1622  rank, line_unique_face_to_elements.size(), tri_unique_face_to_elements.size(),
1623  quad_unique_face_to_elements.size());
1624  fprintf(stderr, "RANK %d VALID entries (index_unique>=0): line=%d, tri=%d, quad=%d\n",
1625  rank, valid_line, valid_tri, valid_quad);
1626  fflush(stderr);
1627  }
1628  #endif
1629 
1630  // ============================================================================
1631  // STEP 5: Finalize mappings after MPI communication
1632  // ============================================================================
1633  for(auto it = line_unique_face_to_elements.begin(); it != line_unique_face_to_elements.end(); ++it) {
1634 
1635  auto& first = it->second.first;
1636  auto& second = it->second.second;
1637 
1638  assign_map_between_elem_oneface_and_elem_uniqueFace(rank, first, second, map_elem_uniqueFace_to_elem_oneface, map_elem_oneface_to_elem_uniqueFace);
1639  }
1640 
1641  for(auto it = tri_unique_face_to_elements.begin(); it != tri_unique_face_to_elements.end(); ++it) {
1642  auto& first = it->second.first;
1643  auto& second = it->second.second;
1644 
1645  assign_map_between_elem_oneface_and_elem_uniqueFace(rank, first, second, map_elem_uniqueFace_to_elem_oneface, map_elem_oneface_to_elem_uniqueFace);
1646  }
1647 
1648  for(auto it = quad_unique_face_to_elements.begin(); it != quad_unique_face_to_elements.end(); ++it) {
1649  auto& first = it->second.first;
1650  auto& second = it->second.second;
1651 
1652  assign_map_between_elem_oneface_and_elem_uniqueFace(rank, first, second, map_elem_uniqueFace_to_elem_oneface, map_elem_oneface_to_elem_uniqueFace);
1653  }
1654 }
1655 
1665 template<class T> inline
1666 void complete_map_vertex_to_dof_with_counter_face(hashmap::unordered_map<tuple<T>,
1667  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
1668  hashmap::unordered_map<triple<T>,
1669  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
1670  hashmap::unordered_map<quadruple<T>,
1671  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face,
1672  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
1673  mesh_int_t> & map_vertex_tag_to_dof)
1674 {
1675  MPI_Comm comm = SF_COMM;
1676  int size, rank;
1677  MPI_Comm_size(comm, &size);
1678  MPI_Comm_rank(comm, &rank);
1679 
1680  // look for counter part line_inteface to add to the map key<vertex,tag> -> value<dof>
1681  for(const auto & v : line_face) {
1682 
1683  // if the ranks are equal, then it has already added to the map
1684  if(v.second.first.rank == v.second.second.rank){
1685  continue;
1686  }
1687  // otherwise we look for the counter part with different rank
1688  emi_face<T,tuple<T>> surf_neighbor =
1689  (v.second.first.rank != rank) ? v.second.first : v.second.second;
1690 
1691  {
1692  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1693  Index_tag_old = std::make_pair(surf_neighbor.points.v1,surf_neighbor.tag);
1694  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1695  }
1696 
1697  {
1698  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1699  Index_tag_old = std::make_pair(surf_neighbor.points.v2,surf_neighbor.tag);
1700  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1701  }
1702  }
1703 
1704  // look for counter part line_inteface to add to the map key<vertex,tag> -> value<dof>
1705  for(const auto & v : tri_face) {
1706 
1707  // if the ranks are equal, then it has already added to the map
1708  if(v.second.first.rank == v.second.second.rank){
1709  continue;
1710  }
1711  // otherwise we look for the counter part with different rank
1712  emi_face<T,triple<T>> surf_neighbor =
1713  (v.second.first.rank != rank) ? v.second.first : v.second.second;
1714 
1715  {
1716  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1717  Index_tag_old = std::make_pair(surf_neighbor.points.v1,surf_neighbor.tag);
1718  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1719  }
1720 
1721  {
1722  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1723  Index_tag_old = std::make_pair(surf_neighbor.points.v2,surf_neighbor.tag);
1724  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1725  }
1726 
1727  {
1728  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1729  Index_tag_old = std::make_pair(surf_neighbor.points.v3,surf_neighbor.tag);
1730  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1731  }
1732  }
1733 
1734  // look for counter part line_inteface to add to the map key<vertex,tag> -> value<dof>
1735  for(const auto & v : quad_face) {
1736 
1737  // if the ranks are equal, then it has already added to the map
1738  if(v.second.first.rank == v.second.second.rank){
1739  continue;
1740  }
1741  // otherwise we look for the counter part with different rank
1742  emi_face<T,quadruple<T>> qsurf_neighbor =
1743  (v.second.first.rank != rank) ? v.second.first : v.second.second;
1744 
1745  {
1746  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1747  Index_tag_old = std::make_pair(qsurf_neighbor.points.v1,qsurf_neighbor.tag);
1748  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1749  }
1750 
1751  {
1752  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1753  Index_tag_old = std::make_pair(qsurf_neighbor.points.v2,qsurf_neighbor.tag);
1754  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1755  }
1756 
1757  {
1758  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1759  Index_tag_old = std::make_pair(qsurf_neighbor.points.v3,qsurf_neighbor.tag);
1760  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1761  }
1762 
1763  {
1764  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
1765  Index_tag_old = std::make_pair(qsurf_neighbor.points.v4,qsurf_neighbor.tag);
1766  map_vertex_tag_to_dof.insert({Index_tag_old,-1});
1767  }
1768  }
1769 
1770  // assign the counter faces which we haven't assigned it yet the new dof
1771  assign_dof_on_counter_face(map_vertex_tag_to_dof,comm);
1772 }
1773 
1803 template<class T, class S> inline
1804 void compute_map_vertex_to_dof(meshdata<T,S> & mesh,
1805  const SF_nbr numbering,
1806  hashmap::unordered_map<T,T> & vertex2ptsdata,
1807  hashmap::unordered_set<int> & extra_tags,
1808  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
1809  mesh_int_t> & map_vertex_tag_to_dof)
1810 {
1811  MPI_Comm comm = mesh.comm;
1812  int size, rank;
1813  MPI_Comm_size(comm, &size);
1814  MPI_Comm_rank(comm, &rank);
1815 
1816  hashmap::unordered_map<mesh_int_t, intersection_data> map_vertex_to_tags_data_ranks;
1817 
1818  const T* con = mesh.con.data();
1819  const T* nbr = mesh.get_numbering(numbering).data();
1820  const SF::vector<T> & rnod = mesh.get_numbering(numbering);
1821 
1822  hashmap::unordered_map<T,T> g2ptsData;
1823 
1824  for(size_t i=0; i<mesh.con.size(); i++)
1825  {
1826  g2ptsData[rnod[con[i]]] = vertex2ptsdata[rnod[con[i]]];
1827  }
1828 
1829  //-----------------------------------------------------------------
1830  // Step 1: Local Information Gathering
1831  //-----------------------------------------------------------------
1832  // - Identify Simple Cases: For any vertex that is not on an interface (i.e., it's completely inside the extracellular space or completely
1833  // inside a unique myocyte), the mapping is simple. The new DOF index is the same as the original vertex index. The function populates the
1834  // output map, map_vertex_tag_to_dof, with these direct mappings, e.g., (vertex {1}, tag_3) -> 1.
1835  // - 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
1836  // MPI ranks that share it, based on the elements the current process knows about. This information is stored in a temporary map called
1837  // map_vertex_to_tags_data_ranks.
1838  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++)
1839  {
1840  T tag = mesh.tag[eidx];
1841  auto pos_extra = extra_tags.find(tag);
1842  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
1843  {
1844  T gIndex = rnod[con[n]];
1845  T data_on_gIndex = g2ptsData[rnod[con[n]]];
1846 
1847  // This handles inner DoFs, either within the extracellular domain or vertices belonging to the intracellular domain.
1848  // and marked as a visited dofs and saved key<gIndex,tag)> -> val<gIndex>
1849  if((pos_extra != extra_tags.end()) || data_on_gIndex==0){
1850  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = std::make_pair(gIndex,tag);
1851  if (map_vertex_tag_to_dof.find(Index_tag_old) == map_vertex_tag_to_dof.end() )
1852  {
1853  map_vertex_tag_to_dof.insert({Index_tag_old,gIndex});
1854  }
1855  }
1856 
1857  // When gIndex is located on an interface — membrane, gap junction, or their intersection.
1858  auto it = map_vertex_to_tags_data_ranks.find(gIndex);
1859  if (it != map_vertex_to_tags_data_ranks.end() )
1860  {
1861  intersection_data& value = it->second;
1862  bool check_tag_rank = true;
1863  for(int i=0; i<MAX_INTERSECTIONS; ++i) {
1864  if(value.tags[i] == tag && value.ranks[i] == rank) {
1865  check_tag_rank = false;
1866  break;
1867  }
1868  }
1869 
1870  // Insert the tag and corresponding rank if they haven’t been inserted yet.
1871  if(check_tag_rank){
1872  for(int i=0; i<MAX_INTERSECTIONS; ++i) {
1873  if(value.tags[i] == -1 && value.ranks[i] == -1) {
1874  value.tags[i] = tag;
1875  value.data[i] = data_on_gIndex;
1876  value.ranks[i] = rank;
1877  break;
1878  }
1879  }
1880  }
1881  }
1882  else
1883  {
1884  // insert gIndex to the map_vertex_to_tags_data_ranks
1885  intersection_data value;
1886  value.tags[0] = tag;
1887  value.data[0] = data_on_gIndex;
1888  value.ranks[0] = rank;
1889  map_vertex_to_tags_data_ranks.insert({gIndex,value});
1890  }
1891  }
1892  }
1893  // until here: map_vertex_to_tags_data_ranks knows about the interfaces for the elements it owns,
1894  // but it doesn't know if a vertex it owns is also part of an interface on another process.
1895 
1896  //-----------------------------------------------------------------
1897  // Step 2: Globalize Information with MPI
1898  //-----------------------------------------------------------------
1899  // - assign_counter_vertices_tuple: This helper function is called to manage the communication. It takes the partial
1900  // map_vertex_to_tags_data_ranks from each process and uses MPI (specifically MPI_Exchange, which is likely a wrapper around routines like
1901  // MPI_Alltoallv) to send and receive data.
1902  // - Build Global Map: After the communication, the partial maps are merged. The result is that the map_vertex_to_tags_data_ranks on every
1903  // process now contains the complete sharing information for every vertex in the entire mesh. Each process now knows the full list of tags
1904  // and ranks associated with any given vertex.
1905  assign_counter_vertices_tuple(map_vertex_to_tags_data_ranks, comm);
1906 
1907  //-----------------------------------------------------------------
1908  // Step 3: Identify and Count DOFs to Be Created
1909  //-----------------------------------------------------------------
1910  // Now that every process has the same global information, they can independently and deterministically decide which new DOFs to create.
1911  // - Analyze Intersections: The code iterates through the now-global map_vertex_to_tags_data_ranks. For each vertex, it inspects the list of
1912  // tags that share it to classify the type of interface.
1913  // - Apply Rules: Based on the interface type, it decides which side needs a new, unique DOF index.
1914  // * Membrane (Myocyte-Extracellular): The extracellular DOF keeps the original vertex index. The myocyte DOF is marked as needing a new, unique index.
1915  // * Gap Junction (Myocyte-Myocyte): To separate the two myocytes, one of them needs a new index. A deterministic rule (e.g., the myocyte
1916  // with the higher tag number gets the new index) is applied to ensure all processes make the same decision.
1917  // * Complex Junctions: For even more complex intersections (e.g., where multiple myocytes and the extracellular space meet), similar
1918  // rules are applied to create the necessary new DOFs.
1919  // - Count Local Contribution: Each process counts how many new DOFs (shift) it is responsible for creating based on these rules.
1920  int shift = 0;
1922  for(const auto & key : map_vertex_to_tags_data_ranks)
1923  {
1924  T gIndex = key.first;
1925  const intersection_data& value = key.second;
1926  T data_on_gIndex = g2ptsData[gIndex];
1927 
1928  // find the first -1 in value
1929  int first_index = 0;
1930  for (int i = 0; i < MAX_INTERSECTIONS; ++i) {
1931  if(value.tags[i] == -1 && value.ranks[i] == -1) {
1932  first_index = i;
1933  break;
1934  }else{
1935  // check if all the intersection_data has the same g2ptsData
1936  if(value.data[i]!=data_on_gIndex)
1937  {
1938  // Error for unhandled cases.
1939  fprintf(stderr,
1940  "WARN: Due to an inconsistency in ptsData computation at gIndex %d with tag %d, "
1941  "the program is stopped while preparing the data structure to introduce DOFs.\n",
1942  gIndex, value.tags[i]);
1943  fprintf(stderr, "\n");
1944  EXIT(1);
1945  }
1946  }
1947  if (i == MAX_INTERSECTIONS -1) first_index = MAX_INTERSECTIONS;
1948  }
1949 
1950  if(data_on_gIndex==1){ // membrane
1951  // check the mumber of tag numbers,
1952  int count_intra_tags = 0;
1953  int count_extra_tags = 0;
1954  for (int i = 0; i < first_index; ++i) {
1955  if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])==extra_tags.end()){
1956  count_intra_tags += 1;
1957  }else if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])!=extra_tags.end()){
1958  count_extra_tags += 1;
1959  }
1960  }
1961 
1962  // Check if the collected data at gIndex is indeed on the membrane.
1963  // A valid membrane vertex must intersect exactly one intracellular tag and one or more extracellular tags.
1964  if(count_intra_tags>1 || count_extra_tags == 0)
1965  {
1966  // Error for unhandled cases.
1967  fprintf(stderr,
1968  "WARN: Due to an issue in defining ptsData for gIndex %d, , on the membrane"
1969  "the program is stopped while preparing the data structure to introduce DOFs.\n",
1970  gIndex);
1971  fprintf(stderr, "\n");
1972  EXIT(1);
1973  }
1974 
1975  int index_intra = -1;
1976  // Find the intracellular tag number. Since data_on_gIndex == 1, there will be only one tag corresponding to the intracellular domain.
1977  for (int i = 0; i < first_index; ++i) {
1978  if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])==extra_tags.end()){
1979  index_intra = i;
1980  break;
1981  }
1982  }
1983  if(index_intra == -1) {
1984  fprintf(stderr, "WARN: unhandled case in membrane while decoupling: gIndex %d", gIndex);
1985  fprintf(stderr, "\n");
1986  EXIT(1);
1987  }
1988  if(index_intra != -1) {
1989  T tag_myocyte = value.tags[index_intra];
1990  T rank_myocyte = value.ranks[index_intra];
1991 
1992  if(rank_myocyte == rank) { // mark intracellular tag with gIndex
1993  std::pair <mesh_int_t,mesh_int_t> Index_tag_next = std::make_pair(gIndex,tag_myocyte);
1994  if (map_mark_new_dofs.find(Index_tag_next) == map_mark_new_dofs.end()) {
1995  map_mark_new_dofs.insert({Index_tag_next,true}); //mark a new index
1996  shift++;
1997  }
1998  }
1999  }
2000  }
2001  else if(data_on_gIndex==2){ // gap junction
2002  // count the number of extra tags and intra tags
2003  int count_intra_tags = 0;
2004  int count_extra_tags = 0;
2005  for (int i = 0; i < first_index; ++i) {
2006  if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])==extra_tags.end()){
2007  count_intra_tags += 1;
2008  }else if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])!=extra_tags.end()){
2009  count_extra_tags += 1;
2010  }
2011  }
2012  if(count_intra_tags!=2 || count_extra_tags != 0){
2013  // Error for unhandled cases.
2014  fprintf(stderr,
2015  "WARN: Due to an issue in defining ptsData for gIndex %d, on the gapjunction"
2016  "the program is stopped while preparing the data structure to introduce DOFs.\n",
2017  gIndex);
2018  fprintf(stderr, "\n");
2019  EXIT(1);
2020  }
2021 
2022  T tag1 = -1;
2023  T tag2 = -1;
2024  T rank1 = -1;
2025  T rank2 = -1;
2026  bool tag1_checked = false;
2027  bool tag2_checked = false;
2028  // find the tag1 and tag2 on the gap junction
2029  for (int i = 0; i < first_index; ++i)
2030  {
2031  auto pos_extra = extra_tags.find(value.tags[i]);
2032  if(pos_extra == extra_tags.end()) {
2033  if(tag1==-1){
2034  tag1 = value.tags[i];
2035  rank1 = value.ranks[i];
2036  tag1_checked = true;
2037  }else if(tag1!=-1){
2038  tag2 = value.tags[i];
2039  rank2 = value.ranks[i];
2040  tag2_checked = true;
2041  break;
2042  }
2043  }
2044  }
2045 
2046  if(!tag1_checked && !tag2_checked || (rank1==rank2 && rank1!=rank) ) {
2047  fprintf(stderr, "WARN: unhandled case in gap junction while decoupling: gIndex %d with tags %d:%d ", gIndex, tag1, tag2);
2048  fprintf(stderr, "\n");
2049  EXIT(1);
2050  }
2051 
2052  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.
2053  T smallest_tag = tag1 < tag2? tag1:tag2;
2054  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = std::make_pair(gIndex,smallest_tag);
2055  if (map_vertex_tag_to_dof.find(Index_tag_old) == map_vertex_tag_to_dof.end() ) {
2056  map_vertex_tag_to_dof.insert({Index_tag_old,gIndex});
2057  }
2058 
2059  T bigger_tag = tag1 < tag2? tag2:tag1;
2060  std::pair <mesh_int_t,mesh_int_t> Index_tag_new = std::make_pair(gIndex,bigger_tag);
2061  if (map_mark_new_dofs.find(Index_tag_new) == map_mark_new_dofs.end() ) {
2062  map_mark_new_dofs.insert({Index_tag_new,true});
2063  shift++;
2064  }
2065  }
2066  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.
2067  if(tag1<tag2) {
2068  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = std::make_pair(gIndex,tag1);
2069  if (map_vertex_tag_to_dof.find(Index_tag_old) == map_vertex_tag_to_dof.end() ) {
2070  map_vertex_tag_to_dof.insert({Index_tag_old,gIndex});
2071  }
2072  } else {
2073  std::pair <mesh_int_t,mesh_int_t> Index_tag_new = std::make_pair(gIndex,tag1);
2074  if (map_mark_new_dofs.find(Index_tag_new) == map_mark_new_dofs.end() ) {
2075  map_mark_new_dofs.insert({Index_tag_new,true});
2076  shift++;
2077  }
2078  }
2079  }
2080  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.
2081  if(tag2<tag1) {
2082  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = std::make_pair(gIndex,tag2);
2083  if (map_vertex_tag_to_dof.find(Index_tag_old) == map_vertex_tag_to_dof.end() ) {
2084  map_vertex_tag_to_dof.insert({Index_tag_old,gIndex});
2085  }
2086  } else {
2087  std::pair <mesh_int_t,mesh_int_t> Index_tag_new = std::make_pair(gIndex,tag2);
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  }
2094  }
2095  else if(data_on_gIndex==3)
2096  {
2097  // count the number of extra tags and intra tags
2098  int count_intra_tags = 0;
2099  int count_extra_tags = 0;
2100  for (int i = 0; i < first_index; ++i) {
2101  if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])==extra_tags.end()){
2102  count_intra_tags += 1;
2103  }else if(value.tags[i]!=-1 && extra_tags.find(value.tags[i])!=extra_tags.end()){
2104  count_extra_tags += 1;
2105  }
2106  }
2107  if(count_intra_tags!=2 || count_extra_tags == 0){
2108  // Error for unhandled cases.
2109  fprintf(stderr,
2110  "WARN: Due to an issue in defining ptsData for gIndex %d, on intersection between membrane and gapjunction"
2111  "the program is stopped while preparing the data structure to introduce DOFs.\n",
2112  gIndex);
2113  fprintf(stderr, "\n");
2114  EXIT(1);
2115  }
2116 
2117  // All extracellular tags have already been added in the previous loop that iterated over the mesh.
2118  // on the interface between gap juntion and membarne
2119  for (int i = 0; i < first_index; ++i) {
2120  // find the first two tag number which belong to itracellular tags, add for indexing
2121  auto pos_extra = extra_tags.find(value.tags[i]);
2122  if(value.tags[i]!=-1 && value.ranks[i] == rank && pos_extra == extra_tags.end()) { // myocytes, for sure we should have two myocytes
2123  std::pair <mesh_int_t,mesh_int_t> Index_tag_new = std::make_pair(gIndex,value.tags[i]);
2124  if (map_mark_new_dofs.find(Index_tag_new) == map_mark_new_dofs.end() )
2125  {
2126  map_mark_new_dofs.insert({Index_tag_new,true}); //mark a new index
2127  shift++;
2128  }
2129  }
2130  }
2131  }
2132  }
2133 
2134  vector<size_t> dsp_dof(size);
2135  MPI_Allgather(&shift, sizeof(size_t), MPI_BYTE, dsp_dof.data(), sizeof(size_t), MPI_BYTE, comm);
2136 
2137  //-----------------------------------------------------------------
2138  // Step 4: Assign Unique Global Indices in Parallel
2139  //-----------------------------------------------------------------
2140  // The final step is to assign the new indices without any conflicts between processes.
2141  // - Calculate Global Offset (`start`): An MPI_Allgather is used to share the shift counts among all processes. Each process then calculates
2142  // 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
2143  // 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
2144  // of new indices created by different processes will not overlap.
2145  // - Assign New Indices: Finally, the function loops through the vertices it marked for new DOF creation and assigns them a unique global
2146  // index (newIndex = start + count). This final mapping, (original_vertex_index, tag) -> new_unique_dof_index, is inserted into the output
2147  // map map_vertex_tag_to_dof.
2148  T start;
2149  if(rank==0){
2150  start = mesh.g_numpts;
2151  }
2152  else
2153  {
2154  start = mesh.g_numpts;
2155  for (int r = 0; r < rank; ++r)
2156  {
2157  start+= dsp_dof[r];
2158  }
2159  }
2160 
2161  int count = 0;
2162  auto lexicographic_comp_pair = [](const std::pair<mesh_int_t,mesh_int_t> & a, const std::pair<mesh_int_t,mesh_int_t> & b)
2163  {
2164  if (a.second == b.second) return a.first < b.first;
2165  return a.second < b.second;
2166  };
2167 
2168  map_mark_new_dofs.sort(lexicographic_comp_pair);
2169 
2170  for(const auto & entry : map_mark_new_dofs)
2171  {
2172  T newIndex = start+count;
2173  map_vertex_tag_to_dof.insert({entry.first,newIndex});
2174  count++;
2175  }
2176 }
2177 
2184 template<class K, class intersection_indices> inline
2185 void assign_counter_dofs(hashmap::unordered_map<K, intersection_indices>& map, const MPI_Comm comm)
2186 {
2187  size_t dsize = map.size();
2188  vector<K> key_vec(dsize);
2189  vector<intersection_indices> value_vec(dsize);
2190  intersection_indices indices;
2191  size_t idx = 0;
2192  for (const auto& v : map) {
2193  key_vec[idx] = v.first;
2194  value_vec[idx] = v.second;
2195  idx++;
2196  }
2197 
2198  int size, rank;
2199  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
2200 
2201  vector<int> perm, dest(dsize);
2202  for (size_t i = 0; i < dsize; i++)
2203  dest[i] = hashmap::hash_ops<K>::hash(key_vec[i]) % size;
2204 
2205  commgraph<size_t> grph;
2206  grph.configure(dest, comm);
2207  size_t nrecv = sum(grph.rcnt);
2208 
2209  interval(perm, 0, dsize);
2210  binary_sort_copy(dest, perm);
2211 
2212  // fill send buffer and communicate
2213  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
2214  vector<intersection_indices> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
2215  for (size_t i = 0; i < dsize; i++) {
2216  sbuff_key[i] = key_vec[perm[i]];
2217  sbuff_value[i] = value_vec[perm[i]];
2218  }
2219  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
2220  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
2221 
2223  for (size_t i = 0; i < nrecv; i++)
2224  {
2225  auto it = rmap.find(rbuff_key[i]);
2226  if (it != rmap.end())
2227  {
2228  intersection_indices& map_indices = it->second;
2229  intersection_indices& r_indices = rbuff_value[i];
2230 
2231  for (int r_indices : r_indices.indices)
2232  {
2233  if (r_indices == -1) continue;
2234  bool found = false;
2235  for (int m_index : map_indices.indices)
2236  {
2237  if (m_index == r_indices) {
2238  found = true;
2239  break;
2240  }
2241  }
2242  if (!found)
2243  {
2244  for (size_t k = 0; k < MAX_INTERSECTIONS; ++k) {
2245  if (map_indices.indices[k] == -1) {
2246  map_indices.indices[k] = r_indices;
2247  break;
2248  }
2249  }
2250  }
2251  }
2252  }
2253  else
2254  {
2255  rmap.insert({ rbuff_key[i], rbuff_value[i] });
2256  }
2257  }
2258 
2259  for (size_t i = 0; i < nrecv; i++) {
2260  auto it = rmap.find(rbuff_key[i]);
2261  if (it != rmap.end()) rbuff_value[i] = it->second;
2262  }
2263 
2264  // sending the assigned data back to original rank
2265  grph.transpose();
2266  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
2267 
2268  for (size_t i = 0; i < dsize; i++) {
2269  auto it = map.find(sbuff_key[i]);
2270  if (it != map.end()) it->second = ibuff_value[i];
2271  }
2272 }
2283 template<class T, class S> inline
2284 void update_emi_mesh_with_dofs(meshdata<T,S> & mesh,
2285  const SF_nbr numbering,
2286  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
2287  mesh_int_t> & map_vertex_tag_to_dof,
2289 {
2290  // map between old index to set of new indices
2292  MPI_Comm comm = mesh.comm;
2293  mesh.globalize(numbering);
2294  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
2295  T tag = mesh.tag[eidx];
2296  T size_elem = mesh.dsp[eidx+1]-mesh.dsp[eidx];
2297 
2298  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
2299  {
2300  T gIndex = mesh.con[n];
2301 
2302  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2303  Index_tag_old = std::make_pair(gIndex,tag);
2304 
2305  auto it_new = map_vertex_tag_to_dof.find(Index_tag_old);
2306  if (it_new != map_vertex_tag_to_dof.end() )
2307  {
2308  mesh_int_t dof_new = (*it_new).second;
2309  dof2vertex.insert({dof_new, gIndex});
2310  mesh.con[n] = dof_new;
2311 
2312  auto it = vertex2dof.find(gIndex);
2313  if (it != vertex2dof.end())
2314  {
2315  intersection_indices& indices = it->second;
2316  bool found = false;
2317  for(T t : indices.indices) {
2318  if(t == dof_new) {
2319  found = true;
2320  break;
2321  }
2322  }
2323  if(!found) {
2324  // 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.
2325  for(size_t i=0; i<MAX_INTERSECTIONS; ++i) {
2326  if(indices.indices[i] == -1) {
2327  indices.indices[i] = dof_new;
2328  break;
2329  }
2330  }
2331  }
2332  }else{
2333  intersection_indices indices;
2334  indices.indices[0] = dof_new;
2335  vertex2dof.insert({gIndex,indices});
2336  }
2337  }
2338  }
2339  }
2340 
2341  assign_counter_dofs(vertex2dof, comm);
2342 
2343  // Now update dof2vertex using vertex2dof
2344  for (const auto& [old_idx, indices] : vertex2dof) {
2345  for (int new_idx : indices.indices) {
2346  if (new_idx != -1) {
2347  // Avoid overwriting unless necessary
2348  if (dof2vertex.find(new_idx) == dof2vertex.end()) {
2349  dof2vertex.insert({new_idx, old_idx});
2350  }
2351  }
2352  }
2353  }
2354 
2355  mesh.localize(numbering);
2356 }
2357 
2369 template<class T, class S> inline
2370 void update_map_indices_to_petsc(meshdata<T,S> & mesh,
2371  const SF_nbr numbering_ref,const SF_nbr numbering_petsc,
2372  hashmap::unordered_set<int> extra_tags,
2373  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
2374  std::pair<mesh_int_t,mesh_int_t>> & map_vertex_tag_to_dof_petsc,
2376  SF::vector<mesh_int_t> & elemTag_emi_mesh)
2377 {
2378  const SF::vector<SF_int> & ref_nbr = mesh.get_numbering(numbering_ref);
2379  const SF::vector<SF_int> & petsc_nbr = mesh.get_numbering(numbering_petsc);
2380 
2381  elemTag_emi_mesh.resize(mesh.l_numelem);
2382  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
2383  T tag = mesh.tag[eidx];
2384  elemTag_emi_mesh[eidx] = 1; // assigned 1 for for extracellular located on extracellular side
2385  if(extra_tags.find(tag) == extra_tags.end())
2386  elemTag_emi_mesh[eidx] = 2; // assigned 2 for for intracellular located on myocyte
2387  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
2388  {
2389  T l_Indx = mesh.con[n];
2390  T petsc_Idx = petsc_nbr[l_Indx];
2391  T n_Indx = ref_nbr[l_Indx];
2392  T o_Indx = dof2vertex[n_Indx];
2393 
2394  std::pair <mesh_int_t,mesh_int_t> oIndx_tag;
2395  oIndx_tag = std::make_pair(o_Indx,tag);
2396 
2397  auto it_new = map_vertex_tag_to_dof_petsc.find(oIndx_tag);
2398  if (it_new != map_vertex_tag_to_dof_petsc.end() )
2399  {
2400  std::pair <mesh_int_t,mesh_int_t> nIndx_petsc;
2401  nIndx_petsc = std::make_pair(n_Indx,petsc_Idx);
2402 
2403  (*it_new).second = nIndx_petsc;
2404  }
2405  }
2406  }
2407 }
2408 
2421 template<class T, class S>
2422 inline void update_faces_on_surface_mesh_after_decoupling_with_dofs(meshdata<T, S> & mesh,
2423  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>,
2424  mesh_int_t> & map_vertex_tag_to_dof,
2425  hashmap::unordered_map<tuple<T>,
2426  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
2427  hashmap::unordered_map<triple<T>,
2428  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
2429  hashmap::unordered_map<quadruple<T>,
2430  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face)
2431 {
2432  MPI_Comm comm = SF_COMM;
2433  int size, rank;
2434  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
2435 
2436  SF::vector<mesh_int_t> idxbuff(mesh.con);
2437  binary_sort(idxbuff); unique_resize(idxbuff);
2438 
2441  for(size_t i=0; i<idxbuff.size(); i++){
2442  g2l[idxbuff[i]] = i;
2443  l2g[i] = idxbuff[i];
2444  }
2445 
2446  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++) {
2447  T tag = mesh.tag[eidx];
2448  T size_elem = mesh.dsp[eidx+1]-mesh.dsp[eidx];
2449  std::vector<int> elem_nodes;
2450  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
2451  {
2452  T gIndex = mesh.con[n];
2453  elem_nodes.push_back(gIndex);
2454  }
2455  std::sort(elem_nodes.begin(),elem_nodes.end());
2456  if(elem_nodes.size()==2){
2457  tuple<T> key;
2458  key.v1 = elem_nodes[0];
2459  key.v2 = elem_nodes[1];
2460  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>> value = line_face[key];
2461 
2462  line_face.erase(key);
2463  std::vector<int> new_nodes(2);
2464 
2465  auto tag_key = (value.first.rank == rank) ? value.first.tag : value.second.tag;
2466  {
2467  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2468  Index_tag_old = std::make_pair(elem_nodes[0],tag_key);
2469  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2470  new_nodes[0] = Index_new;
2471  }
2472  {
2473  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2474  Index_tag_old = std::make_pair(elem_nodes[1],tag_key);
2475  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2476  new_nodes[1] = Index_new;
2477  }
2478  std::sort(new_nodes.begin(),new_nodes.end());
2479  // update key
2480  tuple<T> new_key;
2481  new_key.v1 = new_nodes[0];
2482  new_key.v2 = new_nodes[1];
2483 
2484  // first
2485  {
2486  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v1, value.first.tag)];
2487  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v2, value.first.tag)];
2488  value.first.points.v1 = Index_new1;
2489  value.first.points.v2 = Index_new2;
2490  }
2491  // second
2492  {
2493  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v1, value.second.tag)];
2494  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v2, value.second.tag)];
2495  value.second.points.v1 = Index_new1;
2496  value.second.points.v2 = Index_new2;
2497  }
2498  line_face.insert({new_key,value});
2499  }
2500  if(elem_nodes.size()==3){
2501  triple<T> key;
2502  key.v1 = elem_nodes[0];
2503  key.v2 = elem_nodes[1];
2504  key.v3 = elem_nodes[2];
2505  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>> value = tri_face[key];
2506  tri_face.erase(key);
2507  std::vector<int> new_nodes(3);
2508 
2509  auto tag_key = (value.first.rank == rank) ? value.first.tag : value.second.tag;
2510  {
2511  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2512  Index_tag_old = std::make_pair(elem_nodes[0],tag_key);
2513  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2514  new_nodes[0] = Index_new;
2515  }
2516  {
2517  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2518  Index_tag_old = std::make_pair(elem_nodes[1],tag_key);
2519  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2520  new_nodes[1] = Index_new;
2521  }
2522  {
2523  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2524  Index_tag_old = std::make_pair(elem_nodes[2],tag_key);
2525  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2526  new_nodes[2] = Index_new;
2527  }
2528  std::sort(new_nodes.begin(),new_nodes.end());
2529  // update key
2530  triple<T> new_key;
2531  new_key.v1 = new_nodes[0];
2532  new_key.v2 = new_nodes[1];
2533  new_key.v3 = new_nodes[2];
2534 
2535  // first
2536  {
2537  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v1, value.first.tag)];
2538  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v2, value.first.tag)];
2539  mesh_int_t Index_new3 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v3, value.first.tag)];
2540  value.first.points.v1 = Index_new1;
2541  value.first.points.v2 = Index_new2;
2542  value.first.points.v3 = Index_new3;
2543  }
2544  // second
2545  {
2546  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v1, value.second.tag)];
2547  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v2, value.second.tag)];
2548  mesh_int_t Index_new3 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v3, value.second.tag)];
2549  value.second.points.v1 = Index_new1;
2550  value.second.points.v2 = Index_new2;
2551  value.second.points.v3 = Index_new3;
2552  }
2553  tri_face.insert({new_key,value});
2554  }
2555  if(elem_nodes.size()==4){
2556  quadruple<T> key;
2557  key.v1 = elem_nodes[0];
2558  key.v2 = elem_nodes[1];
2559  key.v3 = elem_nodes[2];
2560  key.v4 = elem_nodes[3];
2561  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>> value = quad_face[key];
2562  quad_face.erase(key);
2563  std::vector<int> new_nodes(4);
2564 
2565  auto tag_key = (value.first.rank == rank) ? value.first.tag : value.second.tag;
2566  {
2567  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2568  Index_tag_old = std::make_pair(elem_nodes[0],tag_key);
2569  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2570  new_nodes[0] = Index_new;
2571  }
2572  {
2573  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2574  Index_tag_old = std::make_pair(elem_nodes[1],tag_key);
2575  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2576  new_nodes[1] = Index_new;
2577  }
2578  {
2579  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2580  Index_tag_old = std::make_pair(elem_nodes[2],tag_key);
2581  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2582  new_nodes[2] = Index_new;
2583  }
2584  {
2585  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
2586  Index_tag_old = std::make_pair(elem_nodes[3],tag_key);
2587  mesh_int_t Index_new = map_vertex_tag_to_dof[Index_tag_old];
2588  new_nodes[3] = Index_new;
2589  }
2590  std::sort(new_nodes.begin(),new_nodes.end());
2591  // update key
2592  quadruple<T> new_key;
2593  new_key.v1 = new_nodes[0];
2594  new_key.v2 = new_nodes[1];
2595  new_key.v3 = new_nodes[2];
2596  new_key.v4 = new_nodes[3];
2597 
2598  // first
2599  {
2600  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v1, value.first.tag)];
2601  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v2, value.first.tag)];
2602  mesh_int_t Index_new3 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v3, value.first.tag)];
2603  mesh_int_t Index_new4 = map_vertex_tag_to_dof[std::make_pair(value.first.points.v4, value.first.tag)];
2604  value.first.points.v1 = Index_new1;
2605  value.first.points.v2 = Index_new2;
2606  value.first.points.v3 = Index_new3;
2607  value.first.points.v4 = Index_new4;
2608  }
2609  // second
2610  {
2611  mesh_int_t Index_new1 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v1, value.second.tag)];
2612  mesh_int_t Index_new2 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v2, value.second.tag)];
2613  mesh_int_t Index_new3 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v3, value.second.tag)];
2614  mesh_int_t Index_new4 = map_vertex_tag_to_dof[std::make_pair(value.second.points.v4, value.second.tag)];
2615  value.second.points.v1 = Index_new1;
2616  value.second.points.v2 = Index_new2;
2617  value.second.points.v3 = Index_new3;
2618  value.second.points.v4 = Index_new4;
2619  }
2620  quad_face.insert({new_key,value});
2621  }
2622  }
2623 }
2624 
2636 template<class T, class S> inline
2637 void compute_surface_mesh_with_counter_face(meshdata<T, S> & mesh, const SF_nbr numbering,
2638  hashmap::unordered_map<tuple<T>,
2639  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
2640  hashmap::unordered_map<triple<T>,
2641  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
2642  hashmap::unordered_map<quadruple<T>,
2643  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face)
2644 {
2645  mesh.register_numbering(numbering);
2646 
2647  MPI_Comm comm = SF_COMM;
2648  int size, rank;
2649  MPI_Comm_size(comm, &size);
2650  MPI_Comm_rank(comm, &rank);
2651 
2652  // To ensure deterministic order, we must sort the keys of the hashmaps before iterating.
2653  vector<tuple<T>> line_keys;
2654  for(auto const& [key, val] : line_face) line_keys.push_back(key);
2655  std::sort(line_keys.begin(), line_keys.end());
2656 
2657  vector<triple<T>> tri_keys;
2658  for(auto const& [key, val] : tri_face) tri_keys.push_back(key);
2659  std::sort(tri_keys.begin(), tri_keys.end());
2660 
2661  vector<quadruple<T>> quad_keys;
2662  for(auto const& [key, val] : quad_face) quad_keys.push_back(key);
2663  std::sort(quad_keys.begin(), quad_keys.end());
2664 
2665  vector<T> cnt(mesh.l_numelem);
2666  size_t idx = 0, cidx = 0;
2667 
2668  // Add faces to the surface mesh used in ionic computation.
2669  // For each {line,tri,quad}_face, add faces only if one of the pair matches the current rank,
2670  // or if both faces have the same rank, add both.
2671  for(const auto & key : line_keys) {
2672  const auto & v = line_face.at(key);
2673  bool both_faces = (v.first.rank == v.second.rank) ? true: false;
2674 
2675  cnt[idx] = 2;
2676 
2677  emi_face<T,tuple<T>> surf_first;
2678  emi_face<T,tuple<T>> surf_second;
2679  if(v.first.rank == rank)
2680  {
2681  surf_first = v.first;
2682  surf_second = v.second;
2683  }
2684  else
2685  {
2686  surf_first = v.second;
2687  surf_second = v.first;
2688  }
2689 
2690  mesh.type[idx] = Line;
2691  mesh.tag[idx] = surf_first.tag;
2692  mesh.con[cidx + 0] = surf_first.points.v1;
2693  mesh.con[cidx + 1] = surf_first.points.v2;
2694  idx += 1;
2695  cidx += 2;
2696 
2697  if(both_faces){
2698 
2699  cnt[idx] = 2;
2700  mesh.type[idx] = Line;
2701  mesh.tag[idx] = surf_second.tag;
2702  mesh.con[cidx + 0] = surf_second.points.v1;
2703  mesh.con[cidx + 1] = surf_second.points.v2;
2704 
2705  idx += 1;
2706  cidx += 2;
2707  }
2708  }
2709 
2710  for(const auto & key : tri_keys) {
2711  const auto & v = tri_face.at(key);
2712  bool both_faces = (v.first.rank == v.second.rank) ? true: false;
2713 
2714  cnt[idx] = 3;
2715 
2716  emi_face<T,triple<T>> surf_first;
2717  emi_face<T,triple<T>> surf_second;
2718  if(v.first.rank == rank)
2719  {
2720  surf_first = v.first;
2721  surf_second = v.second;
2722  }
2723  else
2724  {
2725  surf_first = v.second;
2726  surf_second = v.first;
2727  }
2728 
2729  mesh.type[idx] = Tri;
2730  mesh.tag[idx] = surf_first.tag;
2731  mesh.con[cidx + 0] = surf_first.points.v1;
2732  mesh.con[cidx + 1] = surf_first.points.v2;
2733  mesh.con[cidx + 2] = surf_first.points.v3;
2734 
2735  idx += 1;
2736  cidx += 3;
2737 
2738  if(both_faces){
2739 
2740  cnt[idx] = 3;
2741 
2742  mesh.type[idx] = Tri;
2743  mesh.tag[idx] = surf_second.tag;
2744  mesh.con[cidx + 0] = surf_second.points.v1;
2745  mesh.con[cidx + 1] = surf_second.points.v2;
2746  mesh.con[cidx + 2] = surf_second.points.v3;
2747 
2748  idx += 1;
2749  cidx += 3;
2750  }
2751  }
2752 
2753  for(const auto & key : quad_keys) {
2754  const auto & v = quad_face.at(key);
2755  bool both_faces = (v.first.rank == v.second.rank) ? true: false;
2756  cnt[idx] = 4;
2757 
2758  emi_face<T,quadruple<T>> surf_first;
2759  emi_face<T,quadruple<T>> surf_second;
2760  if(v.first.rank == rank)
2761  {
2762  surf_first = v.first;
2763  surf_second = v.second;
2764  }
2765  else
2766  {
2767  surf_first = v.second;
2768  surf_second = v.first;
2769  }
2770 
2771  mesh.type[idx] = Quad;
2772  mesh.tag[idx] = surf_first.tag;
2773  mesh.con[cidx + 0] = surf_first.points.v1;
2774  mesh.con[cidx + 1] = surf_first.points.v2;
2775  mesh.con[cidx + 2] = surf_first.points.v3;
2776  mesh.con[cidx + 3] = surf_first.points.v4;
2777 
2778  idx += 1;
2779  cidx += 4;
2780 
2781  if(both_faces){
2782 
2783  cnt[idx] = 4;
2784 
2785  mesh.type[idx] = Quad;
2786  mesh.tag[idx] = surf_second.tag;
2787  mesh.con[cidx + 0] = surf_second.points.v1;
2788  mesh.con[cidx + 1] = surf_second.points.v2;
2789  mesh.con[cidx + 2] = surf_second.points.v3;
2790  mesh.con[cidx + 3] = surf_second.points.v4;
2791 
2792  idx += 1;
2793  cidx += 4;
2794  }
2795  }
2796  dsp_from_cnt(cnt, mesh.dsp);
2797 }
2798 
2829 template<class T, class S> inline
2830 void compute_surface_mesh_with_unique_face(meshdata<T, S> & mesh, const SF_nbr numbering,
2831  hashmap::unordered_map<tuple<T>,
2832  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
2833  hashmap::unordered_map<triple<T>,
2834  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
2835  hashmap::unordered_map<quadruple<T>,
2836  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face,
2837  hashmap::unordered_map<mesh_int_t, std::pair<mesh_int_t, mesh_int_t>> & map_elem_uniqueFace_to_tags)
2838 {
2839  mesh.register_numbering(numbering);
2840 
2841  MPI_Comm comm = SF_COMM;
2842  int size, rank;
2843  MPI_Comm_size(comm, &size);
2844  MPI_Comm_rank(comm, &rank);
2845 
2846  // To ensure deterministic order, we must sort the keys of the hashmaps before iterating.
2847  vector<tuple<T>> line_keys;
2848  for(auto const& [key, val] : line_face) line_keys.push_back(key);
2849  std::sort(line_keys.begin(), line_keys.end());
2850 
2851  vector<triple<T>> tri_keys;
2852  for(auto const& [key, val] : tri_face) tri_keys.push_back(key);
2853  std::sort(tri_keys.begin(), tri_keys.end());
2854 
2855  vector<quadruple<T>> quad_keys;
2856  for(auto const& [key, val] : quad_face) quad_keys.push_back(key);
2857  std::sort(quad_keys.begin(), quad_keys.end());
2858 
2859  vector<T> cnt(mesh.l_numelem);
2860  size_t idx = 0, cidx = 0;
2861 
2862  // Add faces to the surface mesh used in ionic computation.
2863  // For each {line,tri,quad}_face, add faces only if one of the pair matches the current rank,
2864  // or if both faces have the same rank, add both.
2865  for(const auto & key : line_keys) {
2866  const auto & v = line_face.at(key);
2867 
2868  bool to_take_face = false;
2869 
2870  emi_face<T,tuple<T>> surf_take;
2871  std::pair <mesh_int_t,mesh_int_t> tag_pairs_take;
2872  {
2873  if(v.first.rank == rank && v.first.mark_to_take == true)
2874  {
2875  surf_take = v.first;
2876  to_take_face = true;
2877  tag_pairs_take = std::make_pair(v.first.tag,v.second.tag);
2878  }
2879  else if(v.second.rank == rank && v.second.mark_to_take == true)
2880  {
2881  surf_take = v.second;
2882  to_take_face = true;
2883  tag_pairs_take = std::make_pair(v.second.tag,v.first.tag);
2884  }
2885 
2886  if(to_take_face)
2887  {
2888  map_elem_uniqueFace_to_tags.insert({idx,tag_pairs_take});
2889  cnt[idx] = 2;
2890  mesh.type[idx] = Line;
2891  mesh.tag[idx] = surf_take.tag;
2892  mesh.con[cidx + 0] = surf_take.points.v1;
2893  mesh.con[cidx + 1] = surf_take.points.v2;
2894  idx += 1;
2895  cidx += 2;
2896  }
2897  }
2898  }
2899 
2900  for(const auto & key : tri_keys) {
2901  const auto & v = tri_face.at(key);
2902 
2903  bool to_take_face = false;
2904 
2905  emi_face<T,triple<T>> surf_take;
2906  std::pair <mesh_int_t,mesh_int_t> tag_pairs_take;
2907  {
2908  if(v.first.rank == rank && v.first.mark_to_take == true)
2909  {
2910  surf_take = v.first;
2911  to_take_face = true;
2912  tag_pairs_take = std::make_pair(v.first.tag,v.second.tag);
2913  }
2914  else if(v.second.rank == rank && v.second.mark_to_take == true)
2915  {
2916  surf_take = v.second;
2917  to_take_face = true;
2918  tag_pairs_take = std::make_pair(v.second.tag,v.first.tag);
2919  }
2920 
2921  if(to_take_face)
2922  {
2923  map_elem_uniqueFace_to_tags.insert({idx,tag_pairs_take});
2924  cnt[idx] = 3;
2925  mesh.type[idx] = Tri;
2926  mesh.tag[idx] = surf_take.tag;
2927  mesh.con[cidx + 0] = surf_take.points.v1;
2928  mesh.con[cidx + 1] = surf_take.points.v2;
2929  mesh.con[cidx + 2] = surf_take.points.v3;
2930 
2931  idx += 1;
2932  cidx += 3;
2933  }
2934  }
2935  }
2936 
2937  for(const auto & key : quad_keys) {
2938  const auto & v = quad_face.at(key);
2939  cnt[idx] = 4;
2940 
2941  bool to_take_face = false;
2942 
2943  emi_face<T,quadruple<T>> surf_take;
2944  std::pair <mesh_int_t,mesh_int_t> tag_pairs_take;
2945  {
2946  if(v.first.rank == rank && v.first.mark_to_take == true)
2947  {
2948  surf_take = v.first;
2949  to_take_face = true;
2950  tag_pairs_take = std::make_pair(v.first.tag,v.second.tag);
2951  }
2952  else if(v.second.rank == rank && v.second.mark_to_take == true)
2953  {
2954  surf_take = v.second;
2955  to_take_face = true;
2956  tag_pairs_take = std::make_pair(v.second.tag,v.first.tag);
2957  }
2958 
2959  if(to_take_face)
2960  {
2961  map_elem_uniqueFace_to_tags.insert({idx,tag_pairs_take});
2962  cnt[idx] = 4;
2963  mesh.type[idx] = Quad;
2964  mesh.tag[idx] = surf_take.tag;
2965  mesh.con[cidx + 0] = surf_take.points.v1;
2966  mesh.con[cidx + 1] = surf_take.points.v2;
2967  mesh.con[cidx + 2] = surf_take.points.v3;
2968  mesh.con[cidx + 3] = surf_take.points.v4;
2969 
2970  idx += 1;
2971  cidx += 4;
2972  }
2973  }
2974  }
2975  dsp_from_cnt(cnt, mesh.dsp);
2976 }
2977 
2997 template<class T> inline
2998 void create_reverse_elem_mapping_between_surface_meshes(
2999  hashmap::unordered_map<tuple<T>, std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>>& line_face,
3000  hashmap::unordered_map<triple<T>, std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>>& tri_face,
3001  hashmap::unordered_map<quadruple<T>, std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>>& quad_face,
3002  SF::vector<T>& vec_one_to_both_face,
3003  MPI_Comm comm)
3004 {
3005  int rank;
3006  MPI_Comm_rank(comm, &rank);
3007 
3008  // To ensure deterministic order, we must sort the keys of the hashmaps before iterating.
3009  vector<tuple<T>> line_keys;
3010  for(auto const& [key, val] : line_face) line_keys.push_back(key);
3011  std::sort(line_keys.begin(), line_keys.end());
3012 
3013  vector<triple<T>> tri_keys;
3014  for(auto const& [key, val] : tri_face) tri_keys.push_back(key);
3015  std::sort(tri_keys.begin(), tri_keys.end());
3016 
3017  vector<quadruple<T>> quad_keys;
3018  for(auto const& [key, val] : quad_face) quad_keys.push_back(key);
3019  std::sort(quad_keys.begin(), quad_keys.end());
3020 
3021  size_t surf_elem_idx = 0;
3022  size_t w_counter_elem_idx = 0;
3023 
3024  size_t lsize_line = 0;
3025  size_t lsize_tri = 0;
3026  size_t lsize_quad = 0;
3027 
3028  // Count how many both-face elements are local on this rank.
3029  // If both faces are on the same rank, the both-face mesh has two entries.
3030  for(const auto& key : line_keys) {
3031  const auto& v = line_face.at(key);
3032  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3033  if (!local_involved) {
3034  continue;
3035  }
3036  if(v.first.rank == v.second.rank) lsize_line+=2;
3037  else lsize_line+=1;
3038  }
3039  for(const auto& key : tri_keys) {
3040  const auto& v = tri_face.at(key);
3041  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3042  if (!local_involved) {
3043  continue;
3044  }
3045  if(v.first.rank == v.second.rank) lsize_tri+=2;
3046  else lsize_tri+=1;
3047  }
3048  for(const auto& key : quad_keys) {
3049  const auto& v = quad_face.at(key);
3050  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3051  if (!local_involved) {
3052  continue;
3053  }
3054  if(v.first.rank == v.second.rank) lsize_quad+=2;
3055  else lsize_quad+=1;
3056  }
3057 
3058  vec_one_to_both_face.resize(lsize_line + lsize_tri + lsize_quad);
3059 
3060  // Fill mapping in deterministic order:
3061  // each both-face element gets the index of its corresponding one-face element.
3062  for (const auto& key : line_keys) {
3063  const auto& v = line_face.at(key);
3064  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3065  if (!local_involved) {
3066  continue;
3067  }
3068  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3069  if (v.first.rank == v.second.rank) {
3070  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3071  }
3072  surf_elem_idx++;
3073  }
3074 
3075  for (const auto& key : tri_keys) {
3076  const auto& v = tri_face.at(key);
3077  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3078  if (!local_involved) {
3079  continue;
3080  }
3081  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3082  if (v.first.rank == v.second.rank) {
3083  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3084  }
3085  surf_elem_idx++;
3086  }
3087 
3088  for (const auto& key : quad_keys) {
3089  const auto& v = quad_face.at(key);
3090  const bool local_involved = (v.first.rank == rank) || (v.second.rank == rank);
3091  if (!local_involved) {
3092  continue;
3093  }
3094  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3095  if (v.first.rank == v.second.rank) {
3096  vec_one_to_both_face[w_counter_elem_idx++] = surf_elem_idx;
3097  }
3098  surf_elem_idx++;
3099  }
3100 }
3101 
3110 template<class T> inline
3111 void added_counter_faces_to_map(hashmap::unordered_map<tuple<T>,
3112  std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>> & line_face,
3113  hashmap::unordered_map<triple<T>,
3114  std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>> & tri_face,
3115  hashmap::unordered_map<quadruple<T>,
3116  std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>> & quad_face)
3117 {
3118  // Buffer inserts not to mutate an unordered_map while iterating over it and later add to unordered_map
3119  std::vector<std::pair<tuple<T>, std::pair<emi_face<T,tuple<T>>, emi_face<T,tuple<T>>>>> line_inserts;
3120  for(const auto & v : line_face) {
3121  SF::emi_face<mesh_int_t,SF::tuple<mesh_int_t>> surf_first = v.second.first;
3122  SF::emi_face<mesh_int_t,SF::tuple<mesh_int_t>> surf_second = v.second.second;
3123 
3124  std::vector<mesh_int_t> elem_nodes_first; elem_nodes_first.resize(2);
3125  std::vector<mesh_int_t> elem_nodes_second; elem_nodes_second.resize(2);
3126 
3127  SF::tuple<mesh_int_t> key_first;
3128  elem_nodes_first[0] = surf_first.points.v1;
3129  elem_nodes_first[1] = surf_first.points.v2;
3130  std::sort(elem_nodes_first.begin(),elem_nodes_first.end());
3131  key_first.v1 = elem_nodes_first[0];
3132  key_first.v2 = elem_nodes_first[1];
3133 
3134  SF::tuple<mesh_int_t> key_second;
3135  elem_nodes_second[0] = surf_second.points.v1;
3136  elem_nodes_second[1] = surf_second.points.v2;
3137  std::sort(elem_nodes_second.begin(),elem_nodes_second.end());
3138  key_second.v1 = elem_nodes_second[0];
3139  key_second.v2 = elem_nodes_second[1];
3140 
3141  const bool same_key_first = (v.first.v1 == key_first.v1 && v.first.v2 == key_first.v2);
3142  if (!same_key_first && line_face.find(key_first) == line_face.end()) {
3143  line_inserts.push_back({key_first, v.second});
3144  }
3145  const bool same_key_second = (v.first.v1 == key_second.v1 && v.first.v2 == key_second.v2);
3146  if (!same_key_second && line_face.find(key_second) == line_face.end()) {
3147  line_inserts.push_back({key_second, v.second});
3148  }
3149  }
3150  for (const auto & entry : line_inserts) line_face.insert(entry);
3151 
3152  std::vector<std::pair<triple<T>, std::pair<emi_face<T,triple<T>>, emi_face<T,triple<T>>>>> tri_inserts;
3153  for(const auto & v : tri_face) {
3154  SF::emi_face<mesh_int_t,SF::triple<mesh_int_t>> surf_first = v.second.first;
3155  SF::emi_face<mesh_int_t,SF::triple<mesh_int_t>> surf_second = v.second.second;
3156 
3157  std::vector<mesh_int_t> elem_nodes_first; elem_nodes_first.resize(3);
3158  std::vector<mesh_int_t> elem_nodes_second; elem_nodes_second.resize(3);
3159 
3160  SF::triple<mesh_int_t> key_first;
3161  elem_nodes_first[0] = surf_first.points.v1;
3162  elem_nodes_first[1] = surf_first.points.v2;
3163  elem_nodes_first[2] = surf_first.points.v3;
3164  std::sort(elem_nodes_first.begin(),elem_nodes_first.end());
3165  key_first.v1 = elem_nodes_first[0];
3166  key_first.v2 = elem_nodes_first[1];
3167  key_first.v3 = elem_nodes_first[2];
3168 
3169  SF::triple<mesh_int_t> key_second;
3170  elem_nodes_second[0] = surf_second.points.v1;
3171  elem_nodes_second[1] = surf_second.points.v2;
3172  elem_nodes_second[2] = surf_second.points.v3;
3173  std::sort(elem_nodes_second.begin(),elem_nodes_second.end());
3174  key_second.v1 = elem_nodes_second[0];
3175  key_second.v2 = elem_nodes_second[1];
3176  key_second.v3 = elem_nodes_second[2];
3177 
3178  const bool same_key_first = (v.first.v1 == key_first.v1 &&
3179  v.first.v2 == key_first.v2 &&
3180  v.first.v3 == key_first.v3);
3181  if (!same_key_first && tri_face.find(key_first) == tri_face.end()) {
3182  tri_inserts.push_back({key_first, v.second});
3183  }
3184  const bool same_key_second = (v.first.v1 == key_second.v1 &&
3185  v.first.v2 == key_second.v2 &&
3186  v.first.v3 == key_second.v3);
3187  if (!same_key_second && tri_face.find(key_second) == tri_face.end()) {
3188  tri_inserts.push_back({key_second, v.second});
3189  }
3190  }
3191  for (const auto & entry : tri_inserts) tri_face.insert(entry);
3192 
3193  std::vector<std::pair<quadruple<T>, std::pair<emi_face<T,quadruple<T>>, emi_face<T,quadruple<T>>>>> quad_inserts;
3194  for(const auto & v : quad_face) {
3195  SF::emi_face<mesh_int_t,SF::quadruple<mesh_int_t>> surf_first = v.second.first;
3196  SF::emi_face<mesh_int_t,SF::quadruple<mesh_int_t>> surf_second = v.second.second;
3197 
3198  std::vector<mesh_int_t> elem_nodes_first; elem_nodes_first.resize(4);
3199  std::vector<mesh_int_t> elem_nodes_second; elem_nodes_second.resize(4);
3200 
3201  SF::quadruple<mesh_int_t> key_first;
3202  elem_nodes_first[0] = surf_first.points.v1;
3203  elem_nodes_first[1] = surf_first.points.v2;
3204  elem_nodes_first[2] = surf_first.points.v3;
3205  elem_nodes_first[3] = surf_first.points.v4;
3206  std::sort(elem_nodes_first.begin(),elem_nodes_first.end());
3207  key_first.v1 = elem_nodes_first[0];
3208  key_first.v2 = elem_nodes_first[1];
3209  key_first.v3 = elem_nodes_first[2];
3210  key_first.v4 = elem_nodes_first[3];
3211 
3212  SF::quadruple<mesh_int_t> key_second;
3213  elem_nodes_second[0] = surf_second.points.v1;
3214  elem_nodes_second[1] = surf_second.points.v2;
3215  elem_nodes_second[2] = surf_second.points.v3;
3216  elem_nodes_second[3] = surf_second.points.v4;
3217  std::sort(elem_nodes_second.begin(),elem_nodes_second.end());
3218  key_second.v1 = elem_nodes_second[0];
3219  key_second.v2 = elem_nodes_second[1];
3220  key_second.v3 = elem_nodes_second[2];
3221  key_second.v4 = elem_nodes_second[3];
3222 
3223  const bool same_key_first = (v.first.v1 == key_first.v1 &&
3224  v.first.v2 == key_first.v2 &&
3225  v.first.v3 == key_first.v3 &&
3226  v.first.v4 == key_first.v4);
3227  if (!same_key_first && quad_face.find(key_first) == quad_face.end()) {
3228  quad_inserts.push_back({key_first, v.second});
3229  }
3230  const bool same_key_second = (v.first.v1 == key_second.v1 &&
3231  v.first.v2 == key_second.v2 &&
3232  v.first.v3 == key_second.v3 &&
3233  v.first.v4 == key_second.v4);
3234  if (!same_key_second && quad_face.find(key_second) == quad_face.end()) {
3235  quad_inserts.push_back({key_second, v.second});
3236  }
3237  }
3238  for (const auto & entry : quad_inserts) quad_face.insert(entry);
3239 }
3240 
3247 template<class K> inline
3248 void assign_counter_tags(hashmap::unordered_map<K, intersection_tags>& map, const MPI_Comm comm)
3249 {
3250  size_t dsize = map.size();
3251  vector<K> key_vec(dsize);
3252  vector<intersection_tags> value_vec(dsize);
3253 
3254  size_t idx = 0;
3255  for (const auto& v : map) {
3256  key_vec[idx] = v.first;
3257  value_vec[idx] = v.second;
3258  idx++;
3259  }
3260 
3261  int size, rank;
3262  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
3263 
3264  vector<int> perm, dest(dsize);
3265  for (size_t i = 0; i < dsize; i++)
3266  dest[i] = hashmap::hash_ops<K>::hash(key_vec[i]) % size;
3267 
3268  commgraph<size_t> grph;
3269  grph.configure(dest, comm);
3270  size_t nrecv = sum(grph.rcnt);
3271 
3272  interval(perm, 0, dsize);
3273  binary_sort_copy(dest, perm);
3274 
3275  // fill send buffer and communicate
3276  vector<K> sbuff_key(dsize), rbuff_key(nrecv);
3277  vector<intersection_tags> sbuff_value(dsize), ibuff_value(dsize), rbuff_value(nrecv);
3278  for (size_t i = 0; i < dsize; i++) {
3279  sbuff_key[i] = key_vec[perm[i]];
3280  sbuff_value[i] = value_vec[perm[i]];
3281  }
3282  MPI_Exchange(grph, sbuff_key, rbuff_key, comm);
3283  MPI_Exchange(grph, sbuff_value, rbuff_value, comm);
3284 
3286  for (size_t i = 0; i < nrecv; i++)
3287  {
3288  auto it = rmap.find(rbuff_key[i]);
3289  if (it != rmap.end())
3290  {
3291  intersection_tags& map_tags = it->second;
3292  intersection_tags& r_tags = rbuff_value[i];
3293 
3294  for (int r_tag : r_tags.tags)
3295  {
3296  if (r_tag == -1) continue;
3297  bool found = false;
3298  for (int m_tag : map_tags.tags)
3299  {
3300  if (m_tag == r_tag) {
3301  found = true;
3302  break;
3303  }
3304  }
3305  if (!found)
3306  {
3307  for (size_t k = 0; k < MAX_INTERSECTIONS; ++k) {
3308  if (map_tags.tags[k] == -1) {
3309  map_tags.tags[k] = r_tag;
3310  break;
3311  }
3312  }
3313  }
3314  }
3315  }
3316  else
3317  {
3318  rmap.insert({ rbuff_key[i], rbuff_value[i] });
3319  }
3320  }
3321 
3322  for (size_t i = 0; i < nrecv; i++) {
3323  auto it = rmap.find(rbuff_key[i]);
3324  if (it != rmap.end()) rbuff_value[i] = it->second;
3325  }
3326 
3327  // sending the assigned data back to original rank
3328  grph.transpose();
3329  MPI_Exchange(grph, rbuff_value, ibuff_value, comm);
3330 
3331  for (size_t i = 0; i < dsize; i++) {
3332  auto it = map.find(sbuff_key[i]);
3333  if (it != map.end()) it->second = ibuff_value[i];
3334  }
3335 }
3336 
3358 template<class T, class S> inline
3359 void compute_ptsdata_from_original_mesh(meshdata<T,S> & mesh,
3360  const SF_nbr numbering,
3361  hashmap::unordered_map<T,T> & vertex2ptsdata,
3362  hashmap::unordered_set<int> extra_tags,
3363  hashmap::unordered_set<int> intra_tags)
3364 {
3365  MPI_Comm comm = mesh.comm;
3366  int size, rank;
3367  MPI_Comm_size(comm, &size);
3368  MPI_Comm_rank(comm, &rank);
3369 
3370  const T* con = mesh.con.data();
3371  const SF::vector<T> & rnod = mesh.get_numbering(numbering);
3372 
3373  // Step 1: Gather all unique region tags associated with each vertex.
3374  // This map is local to each process and will store a list of tags for each vertex
3375  // based on the elements owned by the current process.
3377  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++)
3378  {
3379  T tag = mesh.tag[eidx];
3380  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
3381  {
3382  T gIndex = rnod[con[n]];
3383 
3384  auto it_new = map_index_to_tags.find(gIndex);
3385  if (it_new != map_index_to_tags.end())
3386  {
3387  // If vertex is already in the map, add the new tag if it's not already present.
3388  intersection_tags& tags = it_new->second;
3389  bool found = false;
3390  for(T t : tags.tags) {
3391  if(t == tag) {
3392  found = true;
3393  break;
3394  }
3395  }
3396  if(!found) {
3397  // This is an intracellular tag. Find an empty slot and add it.
3398  for(int i=0; i <MAX_INTERSECTIONS; ++i) {
3399  if (tags.tags[i] == -1) {
3400  tags.tags[i] = tag;
3401  break;
3402  }
3403  }
3404  }
3405  }
3406  else{
3407  // If vertex is not in the map, add it with the current tag.
3408  intersection_tags tags;
3409  tags.tags[0] = tag;
3410  map_index_to_tags.insert({gIndex,tags});
3411  }
3412  }
3413  }
3414 
3415  // Step 2: Globalize the tag information.
3416  // After this call, `map_index_to_tags` on every process will contain the
3417  // complete list of unique tags for every vertex in the entire mesh.
3418  assign_counter_tags(map_index_to_tags, comm);
3419 
3420  // Step 3: Classify each vertex based on its global list of tags.
3421  for(size_t eidx = 0; eidx < mesh.l_numelem; eidx++)
3422  {
3423  for (int n = mesh.dsp[eidx]; n < mesh.dsp[eidx+1];n++)
3424  {
3425  T gIndex = rnod[con[n]];
3426  intersection_tags& vec_tags = map_index_to_tags[gIndex];
3427  int count_intra_tags = 0;
3428  int count_extra_tags = 0;
3429 
3430  // Count the number of unique intracellular and extracellular tags for the vertex.
3431  int count_tags = 0;
3432  for(T t : vec_tags.tags) {
3433  if(intra_tags.count(t)) count_intra_tags++;
3434  if(extra_tags.count(t)) count_extra_tags++;
3435  }
3436 
3437  // Apply classification rules based on the tag counts.
3438  if(count_extra_tags>=1 && count_intra_tags==0){
3439  vertex2ptsdata[gIndex] = 0; // Interior to an extracellular region
3440  }else if(count_extra_tags==0 && count_intra_tags==1){
3441  vertex2ptsdata[gIndex] = 0; // Interior to a myocyte
3442  }else if(count_extra_tags>=1 && count_intra_tags==1){
3443  vertex2ptsdata[gIndex] = 1; // On a membrane (1 myocyte, 1 extracellular)
3444  }else if(count_extra_tags==0 && count_intra_tags==2){
3445  vertex2ptsdata[gIndex] = 2; // On a gap junction (2 myocytes)
3446  }else if(count_extra_tags>=1 && count_intra_tags==2){
3447  vertex2ptsdata[gIndex] = 3; // On a complex junction (2 myocytes, 1 extracellular)
3448  }
3449  else if(count_intra_tags>2){
3450  // Error for unhandled cases.
3451  fprintf(stderr, "More than two intracellular are connected %d. Tags:", gIndex);
3452  for(T tag : vec_tags.tags) if(tag != -1) fprintf(stderr, " %d", tag);
3453  fprintf(stderr, "\n");
3454  EXIT(1);
3455  }
3456  else{
3457  // Error for unhandled cases.
3458  fprintf(stderr, "WARN: unhandled case in ptsData computation for gIndex %d. Tags:", gIndex);
3459  for(T tag : vec_tags.tags) if(tag != -1) fprintf(stderr, " %d", tag);
3460  fprintf(stderr, "\n");
3461  EXIT(1);
3462  }
3463  }
3464  }
3465 }
3466 
3467 
3483 template <class T, class S, class V, class emi_index_rank>
3484 inline void assemble_map_both_to_unique(SF::abstract_matrix<T, S>& op,
3486  std::pair<emi_index_rank, emi_index_rank>>& map,
3487  const SF::meshdata<T, V>& unique_mesh,
3488  const SF::meshdata<T, V>& both_mesh)
3489 {
3490  // Recover global row/column ownership so locally assembled entries can be
3491  // inserted directly with global matrix indices.
3492  int rank;
3493  MPI_Comm_rank(both_mesh.comm, &rank);
3494 
3495  SF::vector<long int> layout_both;
3496  SF::layout_from_count<long int>(both_mesh.l_numelem, layout_both, both_mesh.comm);
3497 
3498  SF::vector<long int> layout_unique;
3499  SF::layout_from_count<long int>(unique_mesh.l_numelem, layout_unique, unique_mesh.comm);
3500 
3501  SF::vector<SF_int> row_idx(1), col_idx(1);
3502  SF::dmat<SF_real> ebuff(1, 1);
3503 
3504  // Assemble one restriction row per locally owned unique face.
3505  for (const auto& [local_unique_idx, both_pair] : map) {
3506  const auto& first_both = both_pair.first;
3507  const auto& second_both = both_pair.second;
3508 
3509  // The unique-face row is owned by this rank, so its global row index comes
3510  // from the unique-face layout plus the local unique-face element id.
3511  T global_unique_row = layout_unique[rank] + local_unique_idx;
3512  row_idx[0] = global_unique_row;
3513 
3514  // Resolve the candidate both-face columns and discard invalid sides.
3515  bool first_valid = (first_both.index >= 0 && first_both.rank >= 0);
3516  bool second_valid = (second_both.index >= 0 && second_both.rank >= 0);
3517 
3518  T global_both_col_first = -1;
3519  T global_both_col_second = -1;
3520  if (first_valid) {
3521  global_both_col_first = layout_both[first_both.rank] + first_both.index;
3522  }
3523  if (second_valid) {
3524  global_both_col_second = layout_both[second_both.rank] + second_both.index;
3525  }
3526 
3527  // If both sides collapse to the same global both-face column, do not insert
3528  // the same contribution twice; treat it as a single-sided restriction row.
3529  if (first_valid && second_valid && global_both_col_first == global_both_col_second) {
3530  second_valid = false;
3531  }
3532 
3533  // Count how many distinct both-face columns contribute to this unique face.
3534  int count = 0;
3535  if (first_valid) count++;
3536  if (second_valid) count++;
3537 
3538  if (count == 0) continue; // No data to average
3539 
3540  // Use arithmetic averaging when both sides exist; otherwise preserve the
3541  // single available value without halving it.
3542  double weight = 0.5;
3543  if (count == 1) weight = 1.0;
3544 
3545  ebuff.assign(1, 1, weight);
3546 
3547  // Insert the first contributing both-face column.
3548  if (first_valid) {
3549  col_idx[0] = global_both_col_first;
3550  op.set_values(row_idx, col_idx, ebuff.data(), false);
3551  }
3552 
3553  // Insert the second contributing both-face column when present.
3554  if (second_valid) {
3555  col_idx[0] = global_both_col_second;
3556  op.set_values(row_idx, col_idx, ebuff.data(), false);
3557  }
3558  }
3559 
3560  op.finish_assembly();
3561 }
3562 
3579 template<class T>
3580 inline void restrict_to_membrane(vector<T> & v, hashmap::unordered_map<T,T> & dof2ptsData, const meshdata<int, float> & mesh)
3581 {
3582  const SF::vector<mesh_int_t> & rnod = mesh.get_numbering(SF::NBR_REF);
3585  for(size_t i=0; i<rnod.size(); i++){
3586  g2l[rnod[i]] = i;
3587  l2g[i] = rnod[i];
3588  }
3589 
3590  size_t widx = 0;
3591 
3592  for(size_t i=0; i<v.size(); i++){
3593  T g = l2g[v[i]];
3594 
3595  if (dof2ptsData[g] > 0) { // dof2ptsData requires global index
3596  v[widx++] = v[i];
3597  }
3598  }
3599 
3600  v.resize(widx);
3601 }
3602 
3603 }
3604 
3605 #endif
3606 #endif
int mesh_int_t
Definition: SF_container.h:46
#define SF_COMM
the default SlimFem MPI communicator
Definition: SF_globals.h:27
Functions related to EMI mesh IO.
Functions handling a distributed mesh.
Definition: mesher.cc:237
virtual void finish_assembly()=0
virtual void set_values(const vector< T > &row_idx, const vector< T > &col_idx, const vector< S > &vals, bool add)=0
Dense matrix class.
Definition: dense_mat.hpp:43
The mesh storage class. It contains both element and vertex data.
Definition: SF_container.h:396
size_t l_numelem
local number of elements
Definition: SF_container.h:399
MPI_Comm comm
the parallel mesh is defined on a MPI world
Definition: SF_container.h:404
A vector storing arbitrary data.
Definition: SF_vector.h:43
size_t size() const
The current size of the vector.
Definition: SF_vector.h:104
void resize(size_t n)
Resize a vector.
Definition: SF_vector.h:209
iterator find(const K &key)
Search for key. Return iterator.
Definition: hashmap.hpp:637
void sort(Compare comp=Compare())
Sort data entries.
Definition: hashmap.hpp:703
void insert(InputIterator first, InputIterator last)
Insert Iterator range.
Definition: hashmap.hpp:583
T & at(const K &key)
Data access by key.
Definition: hashmap.hpp:657
size_t size() const
Definition: hashmap.hpp:731
iterator find(const K &key)
Definition: hashmap.hpp:1092
hm_int count(const K &key) const
Definition: hashmap.hpp:1078
@ Line
Definition: filament.h:16
@ Quad
Definition: filament.h:16
@ Tri
Definition: filament.h:16
Definition: dense_mat.hpp:34
void dsp_from_cnt(const vector< T > &cnt, vector< T > &dsp)
Compute displacements from counts.
Definition: SF_vector.h:310
void interval(vector< T > &vec, size_t start, size_t end)
Create an integer interval between start and end.
Definition: SF_vector.h:350
void binary_sort_copy(vector< T > &_V, vector< S > &_W)
Definition: SF_sort.h:296
T sum(const vector< T > &vec)
Compute sum of a vector's entries.
Definition: SF_vector.h:340
void unique_resize(vector< T > &_P)
Definition: SF_sort.h:348
void count(const vector< T > &data, vector< S > &cnt)
Count number of occurrences of indices.
Definition: SF_vector.h:332
void MPI_Exchange(commgraph< T > &grph, vector< S > &send, vector< S > &recv, MPI_Comm comm)
Exchange data in parallel over MPI.
Definition: SF_network.h:47
void binary_sort(vector< T > &_V)
Definition: SF_sort.h:284
@ Tri
Definition: SF_container.h:60
@ Prism
Definition: SF_container.h:58
@ Pyramid
Definition: SF_container.h:57
@ Tetra
Definition: SF_container.h:54
@ Quad
Definition: SF_container.h:59
@ Hexa
Definition: SF_container.h:55
SF_nbr
Enumeration encoding the different supported numberings.
Definition: SF_container.h:200
@ NBR_ELEM_REF
The element numbering of the reference mesh (the one stored on HD).
Definition: SF_container.h:204
@ NBR_REF
The nodal numbering of the reference mesh (the one stored on HD).
Definition: SF_container.h:201
void dsp_from_cnt(const std::vector< T > &cnt, std::vector< T > &dsp)
Compute displacements from counts.
Definition: kdpart.hpp:140
Base hashing class.
Definition: hashmap.hpp:88
static hm_uint hash(const T &a)
Definition: hashmap.hpp:92