openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
emi.cc
Go to the documentation of this file.
1 //
2 // Copyright (C) 2020 openCARP project
3 //
4 // This program is licensed under the openCARP Academic Public License (APL)
5 // v1.0: You can use and redistribute it and/or modify it in non-commercial
6 // academic environments under the terms of APL as published by the openCARP
7 // project v1.0, or (at your option) any later version. Commercial use requires
8 // a commercial license (info@opencarp.org).
9 //
10 // This program is distributed without any warranty; see the openCARP APL for
11 // more details.
12 //
13 // You should have received a copy of the openCARP APL along with this program
14 // and can find it online: http://www.opencarp.org/license
15 // ----------------------------------------------------------------------------
16 
24 #if WITH_EMI_MODEL
25 
26 #include "petsc_utils.h"
27 #include "physics_types.h"
28 #include "timers.h"
29 #include "stimulate.h"
30 #include "electric_integrators.h"
31 #include "SF_init.h" // for SF::init_xxx()
32 #include "emi.h"
33 #include <array>
34 #include <cmath>
35 
36 #include "SF_petsc_matrix.h"
37 
38 #ifdef WITH_CALIPER
39 #include "caliper/cali.h"
40 #else
41 #include "caliper_hooks.h"
42 #endif
43 
44 namespace opencarp {
45 
46 void log_mesh_local_element_ranges(const sf_mesh& emi_mesh,
47  const sf_mesh& emi_surfmesh_w_counter_face,
48  const sf_mesh& emi_surfmesh_unique_face)
49 {
50  int rank = 0;
51  int comm_size = 0;
52  MPI_Comm_rank(emi_surfmesh_w_counter_face.comm, &rank);
53  MPI_Comm_size(emi_surfmesh_w_counter_face.comm, &comm_size);
54 
55  const unsigned long long local_emi_elems =
56  static_cast<unsigned long long>(emi_mesh.l_numelem);
57  const unsigned long long local_both_face_elems =
58  static_cast<unsigned long long>(emi_surfmesh_w_counter_face.l_numelem);
59  const unsigned long long local_unique_face_elems =
60  static_cast<unsigned long long>(emi_surfmesh_unique_face.l_numelem);
61 
62  std::vector<unsigned long long> all_emi_elems;
63  std::vector<unsigned long long> all_both_face_elems;
64  std::vector<unsigned long long> all_unique_face_elems;
65  if (rank == 0) {
66  all_emi_elems.resize(comm_size, 0);
67  all_both_face_elems.resize(comm_size, 0);
68  all_unique_face_elems.resize(comm_size, 0);
69  }
70 
71  MPI_Gather(&local_emi_elems, 1, MPI_UNSIGNED_LONG_LONG,
72  rank == 0 ? all_emi_elems.data() : nullptr, 1, MPI_UNSIGNED_LONG_LONG,
73  0, emi_surfmesh_w_counter_face.comm);
74  MPI_Gather(&local_both_face_elems, 1, MPI_UNSIGNED_LONG_LONG,
75  rank == 0 ? all_both_face_elems.data() : nullptr, 1, MPI_UNSIGNED_LONG_LONG,
76  0, emi_surfmesh_w_counter_face.comm);
77  MPI_Gather(&local_unique_face_elems, 1, MPI_UNSIGNED_LONG_LONG,
78  rank == 0 ? all_unique_face_elems.data() : nullptr, 1, MPI_UNSIGNED_LONG_LONG,
79  0, emi_surfmesh_w_counter_face.comm);
80 
81  if (rank != 0) return;
82 
83  const auto print_min_max = [](const char* label, const std::vector<unsigned long long>& counts) {
84  if (counts.empty()) return;
85 
86  unsigned long long min_val = counts[0];
87  unsigned long long max_val = counts[0];
88  int min_rank = 0;
89  int max_rank = 0;
90 
91  for (int r = 1; r < static_cast<int>(counts.size()); ++r) {
92  if (counts[r] < min_val) {
93  min_val = counts[r];
94  min_rank = r;
95  }
96  if (counts[r] > max_val) {
97  max_val = counts[r];
98  max_rank = r;
99  }
100  }
101 
102  log_msg(NULL, 0, 0, " %s: \n\t\t min=%llu on rank=%d, \n\t\t max=%llu on rank=%d\n",
103  label, min_val, min_rank, max_val, max_rank);
104  };
105  log_msg(NULL, 0, 0, "\n**********************************");
106  log_msg(NULL, 0, 0, "min/max number of local-element ranges:");
107  print_min_max("emi_mesh", all_emi_elems);
108  print_min_max("emi_surfmesh_w_counter_face", all_both_face_elems);
109  print_min_max("emi_surfmesh_unique_face", all_unique_face_elems);
110  log_msg(NULL, 0, 0, "**********************************");
111 }
112 
113 #ifdef EMI_DEBUG_MESH
114 void log_lhs_positive_definite_probe(SF::abstract_matrix<SF_int, SF_real>* mat,
115  FILE_SPEC logger,
116  const char* stage,
117  int num_trials = 5)
118 {
119  (void)logger;
120  if (param_globals::flavor != std::string("petsc")) return;
121 
122  auto* petsc_mat = dynamic_cast<SF::petsc_matrix*>(mat);
123  if (petsc_mat == nullptr) return;
124 
125  Vec x = NULL, y = NULL;
126  MatCreateVecs(petsc_mat->data, &x, &y);
127 
128  PetscRandom rnd = NULL;
129  PetscRandomCreate(PETSC_COMM_WORLD, &rnd);
130  PetscRandomSetFromOptions(rnd);
131 
132  int rank = 0;
133  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
134 
135  PetscReal min_q = PETSC_MAX_REAL;
136  PetscReal max_q = -PETSC_MAX_REAL;
137  PetscInt nonpos_count = 0;
138 
139  for (int i = 0; i < num_trials; ++i) {
140  VecSetRandom(x, rnd);
141  MatMult(petsc_mat->data, x, y);
142 
143  PetscScalar q_scalar = 0.0;
144  VecDot(x, y, &q_scalar);
145 
146  const PetscReal q = PetscRealPart(q_scalar);
147  min_q = std::min(min_q, q);
148  max_q = std::max(max_q, q);
149  if (q <= 0.0) nonpos_count++;
150  }
151 
152  if (rank == 0) {
153  PetscPrintf(PETSC_COMM_SELF,
154  "%s: SPD probe with %d random vectors: min(x^T A x)=%g, max(x^T A x)=%g, nonpositive=%d\n",
155  stage, num_trials, double(min_q), double(max_q), int(nonpos_count));
156  }
157 
158  PetscRandomDestroy(&rnd);
159  VecDestroy(&x);
160  VecDestroy(&y);
161 }
162 
163 PetscScalar emi_probe_vector_entry(const PetscInt gid, const int probe_id)
164 {
165  const double x = static_cast<double>(gid + 1);
166 
167  switch (probe_id) {
168  case 0:
169  return std::sin(1.0e-3 * x) + 0.5 * std::cos(3.0e-3 * x);
170  case 1:
171  return std::cos(7.0e-4 * x) - 0.35 * std::sin(2.0e-3 * x);
172  default:
173  return 0.75 * std::sin(1.3e-3 * x) + 0.25 * std::cos(4.0e-3 * x);
174  }
175 }
176 
177 void log_lhs_operator_probe(SF::abstract_matrix<SF_int, SF_real>* mat,
178  FILE_SPEC logger,
179  const char* stage,
180  int num_probes = 3)
181 {
182  (void)logger;
183  if (param_globals::flavor != std::string("petsc")) return;
184 
185  auto* petsc_mat = dynamic_cast<SF::petsc_matrix*>(mat);
186  if (petsc_mat == nullptr) return;
187 
188  Vec x = NULL, y = NULL;
189  MatCreateVecs(petsc_mat->data, &x, &y);
190 
191  PetscInt i_start = 0, i_end = 0;
192  VecGetOwnershipRange(x, &i_start, &i_end);
193 
194  PetscScalar* x_arr = NULL;
195  int rank = 0;
196  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
197 
198  for (int probe_id = 0; probe_id < num_probes; ++probe_id) {
199  VecGetArray(x, &x_arr);
200  for (PetscInt i = i_start; i < i_end; ++i) {
201  x_arr[i - i_start] = emi_probe_vector_entry(i, probe_id);
202  }
203  VecRestoreArray(x, &x_arr);
204 
205  MatMult(petsc_mat->data, x, y);
206 
207  PetscReal y_norm2 = 0.0;
208  PetscReal y_norminf = 0.0;
209  PetscScalar y_sum = 0.0;
210  PetscScalar xAy = 0.0;
211  VecNorm(y, NORM_2, &y_norm2);
212  VecNorm(y, NORM_INFINITY, &y_norminf);
213  VecSum(y, &y_sum);
214  VecDot(x, y, &xAy);
215 
216  PetscScalar weighted_checksum_local = 0.0;
217  const PetscScalar* y_arr = NULL;
218  VecGetArrayRead(y, &y_arr);
219  for (PetscInt i = i_start; i < i_end; ++i) {
220  const PetscScalar weight = static_cast<PetscScalar>(i + 1);
221  weighted_checksum_local += weight * y_arr[i - i_start];
222  }
223  VecRestoreArrayRead(y, &y_arr);
224 
225  PetscScalar weighted_checksum = 0.0;
226  MPI_Allreduce(&weighted_checksum_local, &weighted_checksum, 1, MPIU_SCALAR, MPI_SUM, PETSC_COMM_WORLD);
227 
228  if (rank == 0) {
229  PetscPrintf(PETSC_COMM_SELF,
230  "%s: operator probe %d ||Ax||_2=%g, ||Ax||_inf=%g, sum(Ax)=%g, x^T A x=%g, weighted_checksum=%g\n",
231  stage, probe_id + 1, double(y_norm2), double(y_norminf), double(PetscRealPart(y_sum)),
232  double(PetscRealPart(xAy)), double(PetscRealPart(weighted_checksum)));
233  }
234  }
235 
236  VecDestroy(&x);
237  VecDestroy(&y);
238 }
239 
240 void log_rhs_probe(SF::abstract_vector<SF_int, SF_real>* vec,
241  const char* stage)
242 {
243  if (param_globals::flavor != std::string("petsc")) return;
244 
245  auto* petsc_vec = dynamic_cast<SF::petsc_vector*>(vec);
246  if (petsc_vec == nullptr) return;
247 
248  PetscReal norm2 = 0.0;
249  PetscReal norminf = 0.0;
250  PetscScalar sum = 0.0;
251  VecNorm(petsc_vec->data, NORM_2, &norm2);
252  VecNorm(petsc_vec->data, NORM_INFINITY, &norminf);
253  VecSum(petsc_vec->data, &sum);
254 
255  PetscInt i_start = 0, i_end = 0;
256  VecGetOwnershipRange(petsc_vec->data, &i_start, &i_end);
257 
258  PetscScalar weighted_checksum_local = 0.0;
259  const PetscScalar* arr = NULL;
260  VecGetArrayRead(petsc_vec->data, &arr);
261  for (PetscInt i = i_start; i < i_end; ++i) {
262  const PetscScalar weight = static_cast<PetscScalar>(i + 1);
263  weighted_checksum_local += weight * arr[i - i_start];
264  }
265  VecRestoreArrayRead(petsc_vec->data, &arr);
266 
267  PetscScalar weighted_checksum = 0.0;
268  MPI_Allreduce(&weighted_checksum_local, &weighted_checksum, 1, MPIU_SCALAR, MPI_SUM, PETSC_COMM_WORLD);
269 
270  int rank = 0;
271  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
272  if (rank == 0) {
273  PetscPrintf(PETSC_COMM_SELF,
274  "%s: rhs probe ||b||_2=%g, ||b||_inf=%g, sum(b)=%g, weighted_checksum=%g\n",
275  stage, double(norm2), double(norminf), double(PetscRealPart(sum)),
276  double(PetscRealPart(weighted_checksum)));
277  }
278 }
279 
280 void log_linear_system_probe(SF::abstract_matrix<SF_int, SF_real>* mat,
282  const char* stage,
283  int num_probes = 3)
284 {
285  if (param_globals::flavor != std::string("petsc")) return;
286 
287  auto* petsc_mat = dynamic_cast<SF::petsc_matrix*>(mat);
288  auto* petsc_vec = dynamic_cast<SF::petsc_vector*>(vec);
289  if (petsc_mat == nullptr || petsc_vec == nullptr) return;
290 
291  Vec x = NULL, ax = NULL, residual = NULL;
292  MatCreateVecs(petsc_mat->data, &x, &ax);
293  VecDuplicate(ax, &residual);
294 
295  PetscInt i_start = 0, i_end = 0;
296  VecGetOwnershipRange(x, &i_start, &i_end);
297 
298  int rank = 0;
299  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
300 
301  for (int probe_id = 0; probe_id < num_probes; ++probe_id) {
302  PetscScalar* x_arr = NULL;
303  VecGetArray(x, &x_arr);
304  for (PetscInt i = i_start; i < i_end; ++i) {
305  x_arr[i - i_start] = emi_probe_vector_entry(i, probe_id);
306  }
307  VecRestoreArray(x, &x_arr);
308 
309  MatMult(petsc_mat->data, x, ax);
310  VecWAXPY(residual, -1.0, petsc_vec->data, ax);
311 
312  PetscReal residual_norm2 = 0.0;
313  PetscReal residual_norminf = 0.0;
314  PetscScalar residual_sum = 0.0;
315  PetscScalar xTResidual = 0.0;
316  VecNorm(residual, NORM_2, &residual_norm2);
317  VecNorm(residual, NORM_INFINITY, &residual_norminf);
318  VecSum(residual, &residual_sum);
319  VecDot(x, residual, &xTResidual);
320 
321  PetscScalar weighted_checksum_local = 0.0;
322  const PetscScalar* residual_arr = NULL;
323  VecGetArrayRead(residual, &residual_arr);
324  for (PetscInt i = i_start; i < i_end; ++i) {
325  const PetscScalar weight = static_cast<PetscScalar>(i + 1);
326  weighted_checksum_local += weight * residual_arr[i - i_start];
327  }
328  VecRestoreArrayRead(residual, &residual_arr);
329 
330  PetscScalar weighted_checksum = 0.0;
331  MPI_Allreduce(&weighted_checksum_local, &weighted_checksum, 1, MPIU_SCALAR, MPI_SUM, PETSC_COMM_WORLD);
332 
333  if (rank == 0) {
334  PetscPrintf(PETSC_COMM_SELF,
335  "%s: linear-system probe %d ||Ax-b||_2=%g, ||Ax-b||_inf=%g, sum(Ax-b)=%g, x^T(Ax-b)=%g, weighted_checksum=%g\n",
336  stage, probe_id + 1, double(residual_norm2), double(residual_norminf),
337  double(PetscRealPart(residual_sum)), double(PetscRealPart(xTResidual)),
338  double(PetscRealPart(weighted_checksum)));
339  }
340  }
341 
342  VecDestroy(&x);
343  VecDestroy(&ax);
344  VecDestroy(&residual);
345 }
346 #endif
347 
348 void set_elec_tissue_properties_emi_volume(MaterialType* mtype, hashmap::unordered_set<int> &extra_tags, hashmap::unordered_set<int> &intra_tags, FILE_SPEC logger)
349 {
351  MaterialType *m = mtype;
352 
353  // initialize random conductivity fluctuation structure with PrM values
354  m->regions.resize(param_globals::num_gregions);
355 
356  const char* grid_name = "emi_grid_domain";
357  log_msg(logger, 0, 0, "Setting up %s tissue poperties for %d regions ..", grid_name,
358  param_globals::num_gregions);
359 
360  char buf[64];
361  RegionSpecs* reg = m->regions.data();
362 
363  // default tags for extra and intra domains
364  hashmap::unordered_set<int> extra_tags_default = extra_tags;
365  hashmap::unordered_set<int> intra_tags_default = intra_tags;
366 
367  for (size_t i=0; i<m->regions.size(); i++)
368  {
369  for (int j=0;j<param_globals::gregion[i].num_IDs;j++)
370  {
371  int tag = param_globals::gregion[i].ID[j];
372 
373  // removed tags explicitly defined by user in param_globals::gregion[i>1]
374  if(extra_tags_default.find(tag) != extra_tags_default.end())
375  extra_tags_default.erase(tag);
376 
377  if(intra_tags_default.find(tag) != intra_tags_default.end())
378  intra_tags_default.erase(tag);
379  }
380  }
381 
382  for (size_t i=0; i<m->regions.size(); i++, reg++)
383  {
384  if(!strcmp(param_globals::gregion[i].name, "")) {
385  snprintf(buf, sizeof buf, ", gregion_%d", int(i));
386  param_globals::gregion[i].name = dupstr(buf);
387  }
388 
389  // copy metadata into RegionSpecs
390  reg->regname = strdup(param_globals::gregion[i].name);
391  reg->regID = i;
392 
393  if(i==0) // default Extracellular region
394  reg->nsubregs = extra_tags_default.size();
395  if(i==1) // default intracellular region
396  reg->nsubregs = intra_tags_default.size();
397  if(i>1) // optional: rest of other param_globals::gregion[i>1] defined by user
398  reg->nsubregs = param_globals::gregion[i].num_IDs;
399 
400  if(!reg->nsubregs)
401  reg->subregtags = NULL;
402  else
403  {
404  reg->subregtags = new int[reg->nsubregs];
405 
406  if(i==0){
407  int j = 0;
408  for (int tag : extra_tags_default) {
409  reg->subregtags[j] = tag;
410  j++;
411  }
412  }
413  else if(i==1){
414  int j = 0;
415  for (int tag : intra_tags_default) {
416  reg->subregtags[j] = tag;
417  j++;
418  }
419  }
420  else{
421  for (int j=0;j<reg->nsubregs;j++)
422  reg->subregtags[j] = param_globals::gregion[i].ID[j]; // explicit tags defined by user
423  }
424  }
425 
426  // describe material in given region
427  elecMaterial *emat = new elecMaterial();
428  emat->material_type = ElecMat;
429 
430  // Isotropic conductivity is considered in EMI model.
431  emat->InVal[0] = param_globals::gregion[i].g_bath;
432  emat->InVal[1] = param_globals::gregion[i].g_bath;
433  emat->InVal[2] = param_globals::gregion[i].g_bath;
434 
435  emat->ExVal[0] = param_globals::gregion[i].g_bath;
436  emat->ExVal[1] = param_globals::gregion[i].g_bath;
437  emat->ExVal[2] = param_globals::gregion[i].g_bath;
438 
439  emat->BathVal[0] = param_globals::gregion[i].g_bath;
440  emat->BathVal[1] = param_globals::gregion[i].g_bath;
441  emat->BathVal[2] = param_globals::gregion[i].g_bath;
442 
443  // convert units from S/m -> mS/um
444  for (int j=0; j<3; j++) {
445  emat->InVal[j] *= 1e-3 * param_globals::gregion[i].g_mult;
446  emat->ExVal[j] *= 1e-3 * param_globals::gregion[i].g_mult;
447  emat->BathVal[j] *= 1e-3 * param_globals::gregion[i].g_mult;
448  }
449  reg->material = emat;
450  }
451 
452  if(strlen(param_globals::gi_scale_vec))
453  {
454  mesh_t mt = emi_msh;
455  sf_mesh & mesh = get_mesh(mt);
456  const char* file = param_globals::gi_scale_vec;
457 
458  size_t num_file_entries = SF::root_count_ascii_lines(file, mesh.comm);
459 
460  if(num_file_entries != mesh.g_numelem)
461  log_msg(0,4,0, "%s warning: number of %s conductivity scaling entries does not match number of elements!",
462  __func__, get_mesh_type_name(mt));
463 
464  // set up parallel element vector and read data
465  sf_vec *escale;
466  SF::init_vector(&escale, get_mesh(mt), 1, sf_vec::elemwise);
467  escale->read_ascii(param_globals::gi_scale_vec);
468 
469  if(get_size() > 1) {
470  // set up element vector permutation and permute
471  if(get_permutation(mt, ELEM_PETSC_TO_CANONICAL, 1) == NULL) {
473  }
475  sc(*escale, false);
476  }
477 
478  // copy data into SF::vector
479  SF_real* p = escale->ptr();
480  m->el_scale.assign(p, p + escale->lsize());
481  escale->release_ptr(p);
482  }
483 }
484 
485 void parabolic_solver_emi::init()
486 {
488  double t0, t1, dur;
489  get_time(t0);
490  stats.init_logger("par_stats.dat");
491 
492  // Create/initialise linear solver object: PETSc gets constructed with default settings
493  SF::init_solver(&lin_solver);
494 
495  // EMI mesh: After decoupling the interfaces, we obtain a volumetric mesh with new dofs for all vertices.
496  sf_mesh & emi_mesh = get_mesh(emi_msh);
497  // Surface mesh called emi_surface_counter_msh:
498  // This mesh contains only the interface faces, where the indices correspond
499  // to the vertices of the original input mesh.
500  // Each face has a corresponding counter-face, and they are identified by different tag numbers in each rank.
501  sf_mesh & emi_surfmesh_w_counter_face = get_mesh(emi_surface_counter_msh);
502  sf_mesh & emi_surfmesh_one_side = get_mesh(emi_surface_msh);
503  sf_mesh & emi_surfmesh_unique_face = get_mesh(emi_surface_unique_face_msh);
504 
505  int max_row_entries_emi = max_nodal_edgecount(emi_mesh);
506  int max_row_entries_emi_surfmesh = max_nodal_edgecount(emi_surfmesh_w_counter_face);
507 
508  int rank;
509  MPI_Comm_rank(emi_surfmesh_w_counter_face.comm, &rank);
510 
511  sf_vec::ltype alg_type = sf_vec::algebraic;
512  sf_vec::ltype alg_surface_type = sf_vec::elemwise;
513 
514  int dpn = 1;
515 
516  //-----------------------------------------------------------------
517  // setup vectors
518  //-----------------------------------------------------------------
519  SF::init_vector(&ui, emi_mesh, dpn, alg_type);
520  SF::init_vector(&dui, emi_mesh, dpn, alg_type);
521  SF::init_vector(&ui_pre, emi_mesh, dpn, alg_type);
522  SF::init_vector(&Irhs, emi_mesh, dpn, alg_type);
523  SF::init_vector(&Iij_stim, emi_mesh, dpn, alg_type);
524  SF::init_vector(&Iij_temp, emi_mesh, dpn, alg_type);
525  SF::init_vector(&vb_both_face, emi_surfmesh_w_counter_face, dpn, alg_surface_type);
526  SF::init_vector(&vb_unique_face, emi_surfmesh_unique_face, dpn, alg_surface_type);
527  SF::init_vector(&Ib_both_face, emi_surfmesh_w_counter_face, dpn, alg_surface_type);
528  SF::init_vector(&Ib_unique_face, emi_surfmesh_unique_face, dpn, alg_surface_type);
529 
530  //-----------------------------------------------------------------
531  // initialize matrices (LHS, K, M_{emi mesh}, M_{surface mesh})
532  //-----------------------------------------------------------------
533  SF::init_matrix(&lhs_emi);
534  SF::init_matrix(&stiffness_emi);
535  SF::init_matrix(&mass_emi);
536  SF::init_matrix(&mass_surf_emi);
537 
538  lhs_emi->init(emi_mesh, dpn, dpn, max_row_entries_emi*3);
539  stiffness_emi->init(emi_mesh, dpn, dpn, max_row_entries_emi*2);
540  mass_emi->init(emi_mesh, dpn, dpn, max_row_entries_emi*3);
541  mass_surf_emi->init(emi_mesh, dpn, dpn, max_row_entries_emi*3);
542  //---------------------------------------------------------------------
543  // initialize operator matrices B, Bs, Vm_myocyte_emi
544  //---------------------------------------------------------------------
545  mesh_int_t M = emi_surfmesh_w_counter_face.g_numelem;
546  mesh_int_t N = emi_mesh.g_numpts;
547  mesh_int_t m = emi_surfmesh_w_counter_face.l_numelem;
548  mesh_int_t m_one_side = emi_surfmesh_one_side.l_numelem;
549  mesh_int_t M_one_side = emi_surfmesh_one_side.g_numelem;
550  mesh_int_t m_unique_face = emi_surfmesh_unique_face.l_numelem;
551  mesh_int_t M_unique_face = emi_surfmesh_unique_face.g_numelem;
552  mesh_int_t n = ui->lsize();
553 
554  log_msg(NULL, 0, 0, "\n**********************************");
555  log_msg(NULL, 0, 0, "#elements of emi surfmesh unique face: %zu", emi_surfmesh_unique_face.g_numelem);
556  log_msg(NULL, 0, 0, "#elements of emi surfmesh one side: %zu", emi_surfmesh_one_side.g_numelem);
557  log_msg(NULL, 0, 0, "#elements of emi surfmesh: %zu", emi_surfmesh_w_counter_face.g_numelem);
558  log_msg(NULL, 0, 0, "#elements of emi mesh: %zu", emi_mesh.g_numelem);
559  log_msg(NULL, 0, 0, "#dofs for emi_mesh: %zu", emi_mesh.g_numpts);
560  log_msg(NULL, 0, 0, "#max_row_entries_emi: %zu", max_row_entries_emi);
561  log_msg(NULL, 0, 0, "**********************************\n");
562  log_mesh_local_element_ranges(emi_mesh, emi_surfmesh_w_counter_face, emi_surfmesh_unique_face);
563 
564  SF::vector<long int> layout;
565  SF::layout_from_count<long int>(emi_surfmesh_w_counter_face.l_numelem, layout, emi_surfmesh_w_counter_face.comm);
566  mesh_int_t m_l = layout[rank];
567  mesh_int_t n_l = emi_mesh.pl.algebraic_layout()[rank];
568  const SF::vector<mesh_int_t> & alg_nod_surface = emi_surfmesh_w_counter_face.pl.algebraic_nodes();
569 
570  SF::vector<long int> layout_one_side;
571  SF::layout_from_count<long int>(emi_surfmesh_one_side.l_numelem, layout_one_side, emi_surfmesh_one_side.comm);
572  mesh_int_t m_one_side_l = layout_one_side[rank];
573 
574  SF::vector<long int> layout_unique_face;
575  SF::layout_from_count<long int>(emi_surfmesh_unique_face.l_numelem, layout_unique_face, emi_surfmesh_unique_face.comm);
576  mesh_int_t m_unique_face_l = layout_unique_face[rank];
577 
578  // To support a higher number of possible intersections between faces (e.g., multiple contacts on the membrane),
579  // such as when one myocyte is surrounded by more than one extracellular region,
580  // we slightly overallocate the maximum number of row entries (max_row_entries)
581  // for the matrices B, Bi, BsM, and Vm_myocyte_emi to ensure sufficient storage for all coupling cases.
582  SF::init_matrix(&B);
583  B->init(M, N, m, n, m_l, max_row_entries_emi*2);
584  B->zero();
585  SF::init_matrix(&Bi);
586  Bi->init(M, N, m, n, m_l, max_row_entries_emi*2);
587  Bi->zero();
588  SF::init_matrix(&BsM);
589  BsM->init(N, M, n, m, n_l, max_row_entries_emi*2);
590  BsM->zero();
591  SF::init_matrix(&Vm_myocyte_emi);
592  Vm_myocyte_emi->init(N, N, n, n, n_l, max_row_entries_emi*2);
593  Vm_myocyte_emi->zero();
594 
595  // Build the direct unique-face <-> both-face transfer operators
596  SF::construct_direct_unique_both_operators(operator_unique_to_both_faces,
597  operator_both_to_unique_face,
598  map_elem_uniqueFace_to_elem_bothface,
599  map_elem_uniqueFace_to_elem_oneface,
600  vec_both_to_one_face,
601  emi_surfmesh_w_counter_face,
602  emi_surfmesh_unique_face,
603  max_row_entries_emi,
604  dpn,
605  alg_surface_type);
606 
607  // initialize Interpolation operator by B and Scatter operator by Bs
608  SF::assemble_restrict_operator(*B, *Bi, *BsM, elemTag_surface_w_counter_mesh, map_vertex_tag_to_dof_petsc, line_face, tri_face, quad_face, emi_surfmesh_w_counter_face, emi_mesh, UM2_to_CM2);
609 
610  // DEBUG: Check mesh sizes and mappings
611  #ifdef EMI_DEBUG_MESH
612  {
613  int local_rank;
614  MPI_Comm_rank(emi_mesh.comm, &local_rank);
615 
616  fprintf(stderr, "RANK %d MESH SIZES: one_side=%zu, counter=%zu, unique=%zu\n",
617  local_rank, emi_surfmesh_one_side.l_numelem,
618  emi_surfmesh_w_counter_face.l_numelem, emi_surfmesh_unique_face.l_numelem);
619 
620  fprintf(stderr, "RANK %d MAP SIZES: uniqueFace_to_oneface=%zu\n",
621  local_rank, map_elem_uniqueFace_to_elem_oneface.size());
622  fflush(stderr);
623  }
624  #endif
625 
626  //-----------------------------------------------------------------
627  // setup the ionic current
628  //-----------------------------------------------------------------
629  // get initial value of the vm and Iion on the face from the ionicOnFace class
630  sf_vec* vb_ptr = get_data(vm_emi_itf);
631  sf_vec* Ib_ptr = get_data(iion_emi_itf);
632 
633  if(!(vb_ptr != NULL && Ib_ptr != NULL)) {
634  log_msg(0,5,0, "%s error: global Vb and Ib vectors not properly set up! Ionics seem invalid! Aborting!",
635  __func__);
636  EXIT(1);
637  }
638 
639  SF::init_vector(&vb, vb_ptr);
640  vb->shallow_copy(*vb_ptr);
641  SF::init_vector(&Ib, Ib_ptr);
642  Ib->shallow_copy(*Ib_ptr);
643 
644  parab_tech = static_cast<parabolic_solver_emi::parabolic_t>(param_globals::parab_solve_emi);
645 
646  dur = timing(t1, t0);
647 }
648 
649 void parabolic_solver_emi::rebuild_matrices(MaterialType* mtype, limpet::MULTI_IF & miif, SF::vector<stimulus> & stimuli, FILE_SPEC logger)
650 {
652  double start, end, period;
653  get_time(start);
654  double t0, t1, dur;
655  mass_integrator mass_integ;
656  mass_integrator mass_integ_emi;
657  int dpn = 1;
658 
659  int log_flag = param_globals::output_level > 1 ? ECHO : 0;
660  MaterialType & mt = mtype[0];
661  const bool have_dbc = have_dbc_stims(stimuli);
662 
663  double Dt = user_globals::tm_manager->time_step;
664 
665  cond_t condType = intra_cond;
666  sf_mesh & mesh = get_mesh(emi_msh);
667  sf_mesh & emi_surfmesh_w_counter_face = get_mesh(emi_surface_counter_msh);
668 
669 #ifdef EMI_DEBUG_MESH
670  // to check LHS if it is symmetric
671  const auto log_lhs_symmetry = [&](const char* stage) {
672  if (param_globals::flavor != std::string("petsc")) return;
673 
674  auto* petsc_lhs = dynamic_cast<SF::petsc_matrix*>(lhs_emi);
675  if (petsc_lhs == nullptr) return;
676 
677  int rank = 0;
678  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
679 
680  for (PetscReal tol_sym : {PetscReal(1e-14), PetscReal(1e-12), PetscReal(1e-10), PetscReal(1e-8)}) {
681  PetscBool is_sym = PETSC_FALSE;
682  MatIsSymmetric(petsc_lhs->data, tol_sym, &is_sym);
683  if (rank == 0) {
684  PetscPrintf(PETSC_COMM_SELF, "%s: MatIsSymmetric(%g) = %s\n",
685  stage, double(tol_sym), is_sym ? "TRUE" : "FALSE");
686  }
687  }
688  };
689 #endif
690 
691  // set material and conductivity type
692  set_cond_type(mt, condType);
693 
694  // assemble EMI matrices
695  // fill the EMI system
696  {
697  // emi-> stiffness
698  log_msg(NULL, 0, 0, "assemble stiffness matrix");
699  get_time(t0);
700  elec_stiffness_integrator stfn_integ_emi(mt);
701  stiffness_emi->zero();
702  SF::assemble_matrix(*stiffness_emi, mesh, stfn_integ_emi);
703  dur = timing(t1,t0);
704  log_msg(logger,0,log_flag, "Computed parabolic stiffness matrix in %.5f seconds.", float(dur));
705 
706  log_msg(NULL, 0, 0, "assemble mass matrix on the volumetric mesh");
707  get_time(t0);
708  mass_integrator mass_integ;
709  mass_emi->zero();
710  SF::assemble_matrix(*mass_emi, mesh, mass_integ);
711  dur = timing(t1,t0);
712  log_msg(logger,0,log_flag, "Computed volumetric mass matrix in %.5f seconds.", float(dur));
713 
714  log_msg(NULL, 0, 0, "assemble LHS matrix and mass matrix on the surface mesh");
715  get_time(t0);
716  lhs_emi->zero();
717  // lhs_emi = (CmM/dt + K)
718  SF::assemble_lhs_emi(*lhs_emi, *mass_surf_emi, mesh, emi_surfmesh_w_counter_face, map_vertex_tag_to_dof_petsc, line_face, tri_face, quad_face, stfn_integ_emi, mass_integ_emi, -1., UM2_to_CM2 / Dt);
719  dur = timing(t1,t0);
720  log_msg(logger,0,log_flag, "Computed parabolic mass matrix in %.5f seconds.", float(dur));
721 #ifdef EMI_DEBUG_MESH
722  log_lhs_symmetry("lhs_emi after assemble_lhs_emi");
723  log_lhs_positive_definite_probe(lhs_emi, logger, "lhs_emi after assemble_lhs_emi");
724  log_lhs_operator_probe(lhs_emi, logger, "lhs_emi after assemble_lhs_emi");
725 #endif
726  }
727 
728 
729  bool same_nonzero = false;
730 
731  // set boundary conditions
732  if(have_dbc) {
733  log_msg(logger,0,log_flag, "lhs matrix enforcing Dirichlet boundaries.");
734  get_time(t0);
735 
736  if(dbc == nullptr)
737  dbc = new dbc_manager(*lhs_emi, stimuli);
738  else
739  dbc->recompute_dbcs();
740 
741  dbc->enforce_dbc_lhs();
742 #ifdef EMI_DEBUG_MESH
743  log_lhs_symmetry("lhs_emi after enforce_dbc_lhs");
744  log_lhs_positive_definite_probe(lhs_emi, logger, "lhs_emi after enforce_dbc_lhs");
745  log_lhs_operator_probe(lhs_emi, logger, "lhs_emi after enforce_dbc_lhs");
746 #endif
747  dur = timing(t1,t0);
748  log_msg(logger,0,log_flag, "lhs matrix Dirichlet enforcing done in %.5f seconds.", float(dur));
749  }
750  else {
751  log_msg(logger,0,ECHO, "without enforcing Dirichlet boundaries on the lhs matrix!");
752  // we are dealing with a singular system
753  phie_mat_has_nullspace = true;
754  }
755 
756  set_dir(INPUT);
757  get_time(t0);
758  setup_linear_solver(logger);
759  dur = timing(t1,t0);
760  log_msg(logger,0,log_flag, "Initializing parabolic solver in %.5f seconds.", float(dur));
761  set_dir(OUTPUT);
762 
763  period = timing(end, start);
764 }
765 
766 void parabolic_solver_emi::setup_linear_solver(FILE_SPEC logger)
767 {
768  tol = param_globals::cg_tol_parab;
769  max_it = param_globals::cg_maxit_parab;
770 
771  std::string default_opts;
772  std::string solver_file;
773  solver_file = param_globals::parab_options_file;
774  if (param_globals::flavor == std::string("ginkgo")) {
775  default_opts = std::string(
776  R"(
777 {
778  "type": "solver::Cg",
779  "preconditioner": {
780  "type": "solver::Multigrid",
781  "min_coarse_rows": 8,
782  "max_levels": 16,
783  "default_initial_guess": "zero",
784  "mg_level": [
785  {
786  "type": "multigrid::Pgm",
787  "deterministic": false
788  }
789  ],
790  "coarsest_solver": {
791  "type": "preconditioner::Schwarz",
792  "local_solver": {
793  "type": "preconditioner::Jacobi"
794  }
795  },
796  "criteria": [
797  {
798  "type": "Iteration",
799  "max_iters": 1
800  }
801  ]
802  },
803  "criteria": [
804  {
805  "type": "Iteration",
806  "max_iters": 100
807  },
808  {
809  "type": "ResidualNorm",
810  "reduction_factor": 1e-4
811  }
812  ]
813 }
814  )");
815  } else if (param_globals::flavor == std::string("petsc")) {
816  default_opts = std::string("-ksp_type cg -pc_type gamg -options_left");
817  }
818  lin_solver->setup_solver(*lhs_emi, tol, max_it * 100, param_globals::cg_norm_parab,
819  "parabolic PDE", phie_mat_has_nullspace, logger, solver_file.c_str(),
820  default_opts.c_str());
821 }
822 
823 void parabolic_solver_emi::solve()
824 {
825  switch (parab_tech) {
826  case SEMI_IMPLICIT: solve_semiImplicit(); break;
827  }
828 }
829 
830 void parabolic_solver_emi::solve_semiImplicit()
831 {
832  double t0,t1;
833  get_time(t0);
834 
835  if(dbc != nullptr){
836  CALI_MARK_BEGIN("apply_dbc_rhs");
837  dbc->enforce_dbc_rhs(*ui);
838  CALI_MARK_END("apply_dbc_rhs");
839  }
840 
841  *ui_pre = *ui;
842 
843  // -K*u (K is assembled as -K)
844  CALI_MARK_BEGIN("stiff_mat");
845  stiffness_emi->mult(*ui, *Iij_temp);
846  // rhs = Iij + K*u
847  *Irhs -= *Iij_temp;
848  CALI_MARK_END("stiff_mat");
849 
850  // add volumetric stimulus currents (I_ex, I_in)
851  // rhs = Iij + K*u + M*Iij_stim
852  if (Iij_stim->mag() > 0.0) {
853  CALI_MARK_BEGIN("stim_application");
854  mass_emi->mult(*Iij_stim, *Iij_temp);
855  *Irhs -= *Iij_temp;
856  CALI_MARK_BEGIN("stim_application");
857  }
858 
859  // rhs = -(Iij + K*u + M*Iij_stim)
860  CALI_MARK_BEGIN("rhs_update");
861  (*Irhs) *= (-1.0);
862  CALI_MARK_END("rhs_update");
863 
864 #ifdef EMI_DEBUG_MESH
865  log_rhs_probe(Irhs, "rhs before linear_solve");
866  log_linear_system_probe(lhs_emi, Irhs, "rhs before linear_solve");
867 #endif
868 
869  // compute step
870  CALI_MARK_BEGIN("linear_solve");
871  (*lin_solver)(*dui, *Irhs);
872  CALI_MARK_END("linear_solve");
873 
874  // logfile for solver
875  if(lin_solver->reason < 0) {
876  log_msg(0, 5, 0,"%s solver diverged. Reason: %s.", lin_solver->name.c_str(),
877  petsc_get_converged_reason_str(lin_solver->reason));
878  EXIT(1);
879  }
880 
881  // update solution:: ui_pre = ui_pre + dui
882  CALI_MARK_BEGIN("sol_update");
883  ui_pre->add_scaled(*dui, 1.0);
884  *ui *=0;
885  ui->add_scaled(*ui_pre, 1.0);
886  CALI_MARK_END("sol_update");
887 
888  // We need to enforce DBCs again after the solution vector was updated.
889  // Otherwise, matrix/solver tolerances allow small nonzero values in dui to accumulate in ui.
890  if(dbc != nullptr){
891  CALI_MARK_BEGIN("apply_dbc_rhs");
892  dbc->enforce_dbc_rhs(*ui);
893  CALI_MARK_END("apply_dbc_rhs");
894  }
895 
896  // treat solver statistics
897  stats.slvtime += timing(t1, t0);
898  stats.update_iter(lin_solver->niter);
899 }
900 
902  std::pair<SF::emi_face<mesh_int_t,SF::tuple<mesh_int_t>>,
903  SF::emi_face<mesh_int_t,SF::tuple<mesh_int_t>>>> & line_face,
905  std::pair<SF::emi_face<mesh_int_t,SF::triple<mesh_int_t>>,
906  SF::emi_face<mesh_int_t,SF::triple<mesh_int_t>>>> & tri_face,
908  std::pair<SF::emi_face<mesh_int_t,SF::quadruple<mesh_int_t>>,
909  SF::emi_face<mesh_int_t,SF::quadruple<mesh_int_t>>>> & quad_face,
910  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>, mesh_int_t> & map_vertex_tag_to_dof,
911  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>, std::pair<mesh_int_t,mesh_int_t>> & map_vertex_tag_to_dof_petsc,
912  std::vector<std::string> & tags_data)
913 {
914  sf_mesh & emi_surfmesh_w_counter_face = get_mesh(emi_surface_counter_msh);
915 
916  const SF::vector<mesh_int_t> & tags = emi_surfmesh_w_counter_face.tag;
917 
918  const SF::vector<mesh_int_t> & rnod = emi_surfmesh_w_counter_face.get_numbering(SF::NBR_REF);
920 
921  for(size_t i=0; i<rnod.size(); i++){
922  l2g[i] = rnod[i];
923  }
924 
925  for(size_t eidx=0; eidx<emi_surfmesh_w_counter_face.l_numelem; eidx++)
926  {
927  std::vector<int> elem_nodes;
928  mesh_int_t tag = emi_surfmesh_w_counter_face.tag[eidx];
929  for (int n = emi_surfmesh_w_counter_face.dsp[eidx]; n < emi_surfmesh_w_counter_face.dsp[eidx+1];n++)
930  {
931  mesh_int_t l_idx = emi_surfmesh_w_counter_face.con[n];
932 
933  std::pair <mesh_int_t,mesh_int_t> Index_tag_old;
934  Index_tag_old = std::make_pair(l2g[l_idx],tags[eidx]);
935  mesh_int_t dof = map_vertex_tag_to_dof[Index_tag_old];
936  elem_nodes.push_back(dof);
937  }
938 
939  mesh_int_t tag_first = 0;
940  mesh_int_t tag_second = 0;
941  std::string result_first;
942  std::string result_second;
943  std::sort(elem_nodes.begin(),elem_nodes.end()); // make the node tuple order-independent
944 
945  // Extract all surface face pairs separating regions with different material or boundary tags (ionicFaces)
946  if(elem_nodes.size()==2){
948 
949  key.v1 = elem_nodes[0];
950  key.v2 = elem_nodes[1];
951  std::pair<SF::emi_face<mesh_int_t,SF::tuple<mesh_int_t>>,
952  SF::emi_face<mesh_int_t,SF::tuple<mesh_int_t>>> value = line_face[key];
953 
954  tag_first = value.first.tag;
955  tag_second = value.second.tag;
956  result_first = std::to_string(tag_first) + ":" + std::to_string(tag_second);
957  result_second = std::to_string(tag_second) + ":" + std::to_string(tag_first);
958  }
959  else if(elem_nodes.size()==3){
961  key.v1 = elem_nodes[0];
962  key.v2 = elem_nodes[1];
963  key.v3 = elem_nodes[2];
964  std::pair<SF::emi_face<mesh_int_t,SF::triple<mesh_int_t>>,
965  SF::emi_face<mesh_int_t,SF::triple<mesh_int_t>>> value = tri_face[key];
966 
967  tag_first = value.first.tag;
968  tag_second = value.second.tag;
969  result_first = std::to_string(tag_first) + ":" + std::to_string(tag_second);
970  result_second = std::to_string(tag_second) + ":" + std::to_string(tag_first);
971  }
972  else if(elem_nodes.size()==4){
974  key.v1 = elem_nodes[0];
975  key.v2 = elem_nodes[1];
976  key.v3 = elem_nodes[2];
977  key.v4 = elem_nodes[3];
978  std::pair<SF::emi_face<mesh_int_t,SF::quadruple<mesh_int_t>>,
979  SF::emi_face<mesh_int_t,SF::quadruple<mesh_int_t>>> value = quad_face[key];
980 
981  tag_first = value.first.tag;
982  tag_second = value.second.tag;
983  result_first = std::to_string(tag_first) + ":" + std::to_string(tag_second);
984  result_second = std::to_string(tag_second) + ":" + std::to_string(tag_first);
985  }
986  // check if current element tag has a custom region ID
987  tags_data.push_back(result_first);
988  tags_data.push_back(result_second);
989  }
990 }
991 
992 void EMI::initialize()
993 {
994  double t1, t2;
995  get_time(t1);
996 
997  set_dir(OUTPUT);
998 
999  // open logger
1000  logger = f_open("emi.log", param_globals::experiment != 4 ? "w" : "r");
1001  const int verb = param_globals::output_level;
1002  // Mesh processing step: convert the input mesh into
1003  // - an EMI mesh (with discontinuities at gap junctions and membranes)
1004  // - and generate the corresponding surface mesh
1005  CALI_MARK_BEGIN("mesh_setup");
1006  setup_EMI_mesh();
1007 
1008  // setup mappings between extra and intra grids, algebraic and nodal,
1009  // and between PETSc and canonical orderings
1010  setup_mappings();
1011 
1012  // the ionicOnFace physics is currently triggered from inside the emi to have tighter
1013  // control over it
1014  ion.logger = logger;
1015 
1016  ion.set_surface_mesh_data(parab_solver.line_face,
1017  parab_solver.tri_face,
1018  parab_solver.quad_face,
1019  parab_solver.map_vertex_tag_to_dof);
1020 
1021  // builds per-face adjacency tags
1022  std::vector<std::string> tags_data;
1023  tags_onFace(parab_solver.line_face,
1024  parab_solver.tri_face,
1025  parab_solver.quad_face,
1026  parab_solver.map_vertex_tag_to_dof,
1027  parab_solver.map_vertex_tag_to_dof_petsc,
1028  tags_data);
1029  CALI_MARK_END("mesh_setup");
1030 
1031  ion.set_tags_onFace(tags_data);
1032 
1033  ion.set_face_region_data(parab_solver.intra_tags, parab_solver.map_elem_uniqueFace_to_tags);
1034  ion.initialize();
1035 
1036  // set up tissue properties on the extra and intracellular domains
1037  set_elec_tissue_properties_emi_volume(mtype_vol, parab_solver.extra_tags, parab_solver.intra_tags, logger);
1038  region_mask(emi_msh, mtype_vol[0].regions, mtype_vol[0].regionIDs, true, "gregion_vol");
1039 
1040  // add electrics timer for time stepping, add to time stepper tool (TS)
1041  double global_time = user_globals::tm_manager->time;
1042  timer_idx = user_globals::tm_manager->add_eq_timer(global_time, param_globals::tend, 0,
1043  param_globals::dt, 0, "elec::ref_dt", "TS");
1044 
1045  // EMI stimuli setup
1046  CALI_MARK_BEGIN("stimulus_setup");
1047  param_globals::operator_splitting = 0; // as long as we don't use operator splitting, we default to 0 so we can keep the stimuli scaling mint
1048  setup_stimuli();
1049  CALI_MARK_END("stimulus_setup");
1050 
1051  // set up the linear equation systems. this needs to happen after the stimuli have been
1052  // set up, since we need boundary condition info
1053  CALI_MARK_BEGIN("solver_setup");
1054  setup_solvers();
1055  CALI_MARK_END("solver_setup");
1056 
1057  // balance electrodes, we may need the extracellular mass matrix
1058  balance_electrodes();
1059  // total current scaling
1060  scale_total_stimulus_current(stimuli, *parab_solver.mass_emi, *parab_solver.mass_surf_emi, logger);
1061 
1062  sf_mesh & emi_mesh = get_mesh(emi_msh);
1063  sf_mesh & emi_surfmesh_w_counter_face = get_mesh(emi_surface_counter_msh);
1064 
1065  // - Extract the resting potential from v_b, which is defined only on the membrane faces.
1066  // (we do not care about gap junction in this case)
1067  // - Project this resting potential from the surface mesh onto the EMI mesh.
1068  // - assign the corresponding resting potential on the memebrane to all DOFs
1069  // myocyte that is associated with those membrane faces.
1070  // (dofs on the interfaces of myocytes and all inner dofs located inside of myocytes).
1071  // - any DOF with tag number belongs to the extra_tags set, is assigned to zero.
1072  // - We provide a mapping between the vertices of the surface mesh(vertices of the original input mesh)
1073  // and those of the EMI mesh (with updated indices).
1074  // Based on the face type (line_face, tri_face, or quad_face) and the tag numbers that
1075  // identify paired faces (face and counter-face),
1076  // we can determine whether a given face belongs to a membrane or a gap-junction region.
1077  // - Direct unique -> both mapping
1078  parab_solver.operator_unique_to_both_faces->mult(*parab_solver.vb, *parab_solver.vb_both_face);
1079  SF::assign_resting_potential_from_ionic_models_on_myocyte(*parab_solver.Vm_myocyte_emi,
1080  parab_solver.vb_both_face,
1081  parab_solver.elemTag_emi_mesh,
1082  parab_solver.map_vertex_tag_to_dof_petsc,
1083  parab_solver.line_face, parab_solver.tri_face, parab_solver.quad_face,
1084  emi_surfmesh_w_counter_face, emi_mesh);
1085 
1086  parab_solver.Vm_myocyte_emi->get_diagonal(*parab_solver.ui); // initialize all myocytes to the resting potential based on ionic models
1087  *parab_solver.vb_unique_face = *parab_solver.vb;
1088 
1089  CALI_MARK_BEGIN("output_setup");
1090  // prepare the electrics output. we skip it if we do post-processing
1091  if(param_globals::experiment != EXP_POSTPROCESS)
1092  setup_output();
1093  CALI_MARK_END("output_setup");
1094 
1095  this->initialize_time += timing(t2, t1);
1096 }
1097 
1098 void EMI::setup_mappings()
1099 {
1100  bool emi_exits = mesh_is_registered(emi_msh);
1101  assert(emi_exits);
1102  const int dpn = 1;
1103 
1104  if(get_scattering(emi_msh, ALG_TO_NODAL, dpn) == NULL)
1105  {
1106  log_msg(logger, 0, 0, "%s: Setting up intracellular algebraic-to-nodal scattering.", __func__);
1108  }
1109  if(get_permutation(emi_msh, PETSC_TO_CANONICAL, dpn) == NULL)
1110  {
1111  log_msg(logger, 0, 0, "%s: Setting up intracellular PETSc to canonical permutation.", __func__);
1113  }
1114 }
1115 
1116 void EMI::compute_step()
1117 {
1118  double t1, t2;
1119  get_time(t1);
1120 
1121  // activation checking
1122  const double time = user_globals::tm_manager->time,
1123  time_step = user_globals::tm_manager->time_step;
1124 
1125  const int verb = param_globals::output_level;
1126  // We are treating stimuli by type:
1127  // - Potential stimuli (Phi_ex, GND_ex, Phi_ex_ol) are managed by a dbc_manager and are applied
1128  // to the left and right hand side of the equation.
1129  // - Current stimuli (I_ex, I_in) are applied to the right hand side,
1130  // while (I_tm) is applied to the vector Ib directly.
1131  CALI_MARK_BEGIN("apply_dbc_lhs");
1132  apply_dbc_stimulus();
1133  CALI_MARK_END("apply_dbc_lhs");
1134 
1135  // compute ionics update
1136  CALI_MARK_BEGIN("ion_compute");
1137  ion.compute_step();
1138  CALI_MARK_END("ion_compute");
1139 
1140  CALI_MARK_BEGIN("apply_stim");
1141  apply_current_stimulus();
1142  CALI_MARK_END("apply_stim");
1143 
1144  // convert Ib -> Irhs
1145  parab_solver.operator_unique_to_both_faces->mult(*parab_solver.Ib, *parab_solver.Ib_both_face);
1146  parab_solver.BsM->mult(*parab_solver.Ib_both_face, *parab_solver.Irhs);
1147 
1148  // solver parabolic system
1149  CALI_MARK_BEGIN("parab_solve");
1150  parab_solver.solve();
1151  {
1152  // v_b = B_i * u
1153  parab_solver.B->mult(*parab_solver.ui, *parab_solver.vb_both_face);
1154  // Direct both -> unique mapping
1155  parab_solver.operator_both_to_unique_face->mult(*parab_solver.vb_both_face, *parab_solver.vb_unique_face);
1156  *parab_solver.vb = *parab_solver.vb_unique_face;
1157  }
1158  CALI_MARK_END("parab_solve");
1159 
1160  if(user_globals::tm_manager->trigger(iotm_console)) {
1161  // output lin solver stats
1162  parab_solver.stats.log_stats(user_globals::tm_manager->time, false);
1163  }
1164  this->compute_time += timing(t2, t1);
1165 
1166  // since the traces have their own timing, we check for trace dumps in the compute step loop
1167  if(user_globals::tm_manager->trigger(iotm_trace))
1169 }
1170 
1171 void EMI::output_step()
1172 {
1173  double t1, t2;
1174  get_time(t1);
1175 
1176  output_manager.write_data();
1177 
1178  double curtime = timing(t2, t1);
1179  this->output_time += curtime;
1180 
1181  IO_stats.calls++;
1182  IO_stats.tot_time += curtime;
1183 
1185  IO_stats.log_stats(user_globals::tm_manager->time, false);
1186 }
1187 
1191 void EMI::destroy()
1192 {
1194  // close logger
1195  f_close(logger);
1196 
1197  // close output files
1198  output_manager.close_files_and_cleanup();
1199 
1200  // destroy ionics
1201  ion.destroy();
1202 }
1203 
1204 void EMI::setup_stimuli()
1205 {
1206  // initialize basic stim info data (used units, supported types, etc)
1207  init_stim_info();
1208 
1209  stimuli.resize(param_globals::num_stim);
1210  for (int i = 0; i < param_globals::num_stim; i++) {
1211  // construct new stimulus
1212  stimulus & s = stimuli[i];
1213 
1215  s.translate(i);
1216 
1217  // we associate to the EMI mesh. this is needed for the stim_phys and stim_electrode setups.
1218  s.associated_intra_mesh = emi_msh, s.associated_extra_mesh = emi_msh;
1219 
1220  s.setup(i);
1221 
1222  if (s.phys.type == Illum) {
1223  log_msg(0, MAX_LOG_LEVEL, ECHO, "Stimulus of type Illum (=6) is not implemented in EMI. Abort.");
1224  EXIT(EXIT_FAILURE);
1225  }
1226 
1227  // Depending on the stimulus type, we make sure to only stimulate the correct regions of the mesh:
1228  // Extracellular stimuli restrict to all DOFs of the extracellular region (including the extracellular side of the membrane).
1229  // Equivalent for intracellular stimuli.
1230  // Stimuli that act directly on the membrane restrict to DOFs with ptsData > 0, i.e. membrane, gap junctions, and complex junctions.
1231  if (is_extra(s.phys.type)) {
1232  SF::vector<mesh_int_t> extra_vertices;
1233  const sf_mesh& mesh = get_mesh(s.associated_extra_mesh);
1234 
1235  // Gather vertices from emi_msh using extra_tags
1236  indices_from_region_tags(extra_vertices, mesh, parab_solver.extra_tags);
1237 
1238  // Restrict electrode vertices to extra region
1239  restrict_to_set(s.electrode.vertices, extra_vertices);
1240  } else if (s.phys.type == I_tm) {
1241  const sf_mesh& mesh = get_mesh(emi_msh);
1242  SF::restrict_to_membrane(s.electrode.vertices, dof2ptsData, mesh);
1243  } else {
1244  SF::vector<mesh_int_t> intra_vertices;
1245  const sf_mesh& mesh = get_mesh(s.associated_intra_mesh);
1246 
1247  // Gather vertices from emi_msh using intra_tags
1248  indices_from_region_tags(intra_vertices, mesh, parab_solver.intra_tags);
1249 
1250  // Restrict electrode vertices to intra region
1251  restrict_to_set(s.electrode.vertices, intra_vertices);
1252  }
1253 
1254  if (s.electrode.dump_vtx) {
1255  set_dir(OUTPUT);
1256  s.dump_vtx_file(i);
1257  }
1258 
1259  if(param_globals::stim[i].pulse.dumpTrace && get_rank() == 0) {
1260  set_dir(OUTPUT);
1261  s.pulse.wave.write_trace(s.name+".trc");
1262  }
1263 
1264  }
1265 }
1266 
1267 void EMI::apply_dbc_stimulus()
1268 {
1269  parabolic_solver_emi& ps = parab_solver;
1270 
1271  // we check if the DBC layout changed, if so we recompute the matrix and the dbc_manager
1272  bool dbcs_have_updated = ps.dbc != nullptr && ps.dbc->dbc_update();
1274 
1275  if (dbcs_have_updated && time_not_final) {
1276  parab_solver.rebuild_matrices(mtype_vol, *ion.miif, stimuli, logger);
1277  }
1278 }
1279 
1280 void EMI::apply_current_stimulus()
1281 {
1282  parabolic_solver_emi& ps = parab_solver;
1283  ps.Iij_stim->set(0.0);
1284 
1285  // iterate over stimuli
1286  for(stimulus & s : stimuli) {
1287  if(s.is_active()) {
1288  switch (s.phys.type) {
1289  case I_tm: {
1290  apply_stim_to_vector(s, *ps.Iij_temp, true);
1291  ps.Bi->mult(*ps.Iij_temp, *ps.Ib_both_face);
1292  ps.operator_both_to_unique_face->mult(*ps.Ib_both_face, *ps.Ib_unique_face);
1293  ps.Ib->add_scaled(*ps.Ib_unique_face, -0.5); // scale by 0.5 to revert 2x scaling in Bi*Iij
1294  } break;
1295 
1296  case I_ex:
1297  case I_in: {
1298  apply_stim_to_vector(s, *ps.Iij_stim, true);
1299  } break;
1300 
1301  default: break;
1302  }
1303  }
1304  }
1305 }
1306 
1307 void EMI::balance_electrodes()
1308 {
1309  for (int i = 0; i < param_globals::num_stim; i++) {
1310  if (param_globals::stim[i].crct.balance != -1) {
1311  int from = param_globals::stim[i].crct.balance;
1312  int to = i;
1313 
1314  log_msg(NULL, 0, 0, "Balancing stimulus %d with %d %s-wise.", from, to,
1315  is_current(stimuli[from].phys.type) ? "current" : "voltage");
1316 
1317  stimulus& s_from = stimuli[from];
1318  stimulus& s_to = stimuli[to];
1319 
1320  s_to.pulse = s_from.pulse;
1321  s_to.ptcl = s_from.ptcl;
1322  s_to.phys = s_from.phys;
1323  s_to.pulse.strength *= -1.0;
1324 
1325  if (s_from.phys.type == I_ex || s_from.phys.type == I_in) {
1326  // if from is total current, skip volume based adjustment of strength
1327  // otherwise, calling constant_total_stimulus_current() will undo the balanced scaling of to.pulse.strength
1328  // constant_total_stimulus_current() will do the scaling based on the volume
1329  if (!s_from.phys.total_current) {
1330  sf_mat& mass = *parab_solver.mass_emi;
1331  SF_real vol0 = get_volume_from_nodes(mass, s_from.electrode.vertices);
1332  SF_real vol1 = get_volume_from_nodes(mass, s_to.electrode.vertices);
1333 
1334  s_to.pulse.strength *= fabs(vol0 / vol1);
1335  }
1336  }
1337  }
1338  }
1339 }
1340 
1341 void EMI::scale_total_stimulus_current(SF::vector<stimulus>& stimuli,
1342  sf_mat& mass_vol,
1343  sf_mat& mass_surf,
1344  FILE_SPEC logger)
1345 {
1346  for (stimulus & s : stimuli){
1347  if(is_current(s.phys.type) && s.phys.total_current){
1348  switch (s.phys.type) {
1349  case I_in:
1350  case I_ex: {
1351  // compute affected volume in um^3
1352  SF_real vol = get_volume_from_nodes(mass_vol, s.electrode.vertices);
1353  // s->strength holds the total current in uA, compute current density in uA/cm^3
1354  // Theoretically, we don't need to scale the volume to cm^3 here since we later
1355  // multiply with the mass matrix and we get um^3 * uA/um^3 = uA.
1356  // However, for I_ex/I_in there is an additional um^3 to cm^3 scaling in phys.scale,
1357  // since I_ex/I_in is expected to be in uA/cm^3. Therefore, we need to compensate for that to arrive at uA later.
1358  assert(vol > 0);
1359  float scale = 1.e12 / vol;
1360 
1361  s.pulse.strength *= scale;
1362 
1363  log_msg(logger, 0, ECHO,
1364  "%s [Stimulus %d]: current density scaled to %.4g uA/cm^3\n",
1365  s.name.c_str(), s.idx, s.pulse.strength);
1366  } break;
1367 
1368  case I_tm: {
1369  // In the EMI model, I_tm only affects the membrane. Therefore, we compute
1370  // the affected membrane surface in um^2 using the membrane mass matrix.
1371  // The electrode vertices are resticted to the membrane during setup, hence this function returns a surface already.
1372  SF_real surf = get_volume_from_nodes(mass_surf, s.electrode.vertices);
1373 
1374  // convert to cm^2
1375  assert(surf > 0);
1376  surf /= 1.e8;
1377 
1378  // scale surface density now to result in correct total current
1379  s.pulse.strength /= surf;
1380  log_msg(logger, 0, ECHO,
1381  "%s [Stimulus %d]: current density scaled to %.4g uA/cm^2\n",
1382  s.name.c_str(), s.idx, s.pulse.strength);
1383  } break;
1384 
1385  default: break;
1386  }
1387  }
1388  }
1389 }
1390 
1391 // Assign a deterministic global element numbering for a surface mesh based on
1392 // element type, tag, and sorted node ids. This only affects output ordering.
1393 static void assign_deterministic_elem_numbering(sf_mesh & mesh)
1394 {
1395  const int KEY_SIZE = 6; // type, tag, n1, n2, n3, n4
1396  int rank = 0, size = 0;
1397  MPI_Comm_rank(mesh.comm, &rank);
1398  MPI_Comm_size(mesh.comm, &size);
1399 
1400  auto make_key = [&](size_t i) {
1401  std::array<long long, KEY_SIZE> k;
1402  k.fill(-1);
1403  k[0] = (long long)mesh.type[i];
1404  k[1] = (long long)mesh.tag[i];
1405 
1406  int nn = 0;
1407  if (mesh.type[i] == SF::Line) nn = 2;
1408  else if (mesh.type[i] == SF::Tri) nn = 3;
1409  else if (mesh.type[i] == SF::Quad) nn = 4;
1410  else nn = 0;
1411 
1412  std::vector<long long> nodes;
1413  nodes.reserve(nn);
1414  size_t off = mesh.dsp[i];
1415  for (int j = 0; j < nn; j++) {
1416  nodes.push_back((long long)mesh.con[off + j]);
1417  }
1418  std::sort(nodes.begin(), nodes.end());
1419  for (int j = 0; j < (int)nodes.size(); j++) {
1420  k[2 + j] = nodes[j];
1421  }
1422  return k;
1423  };
1424 
1425  // Pack local keys
1426  std::vector<long long> local_keys(mesh.l_numelem * KEY_SIZE, -1);
1427  for (size_t i = 0; i < mesh.l_numelem; i++) {
1428  auto k = make_key(i);
1429  for (int j = 0; j < KEY_SIZE; j++) local_keys[i * KEY_SIZE + j] = k[j];
1430  }
1431 
1432  // Gather sizes
1433  std::vector<int> counts(size, 0), displs(size, 0);
1434  int local_count = (int)local_keys.size();
1435  MPI_Allgather(&local_count, 1, MPI_INT, counts.data(), 1, MPI_INT, mesh.comm);
1436  int total = 0;
1437  for (int r = 0; r < size; r++) {
1438  displs[r] = total;
1439  total += counts[r];
1440  }
1441 
1442  std::vector<long long> all_keys;
1443  if (rank == 0) all_keys.resize(total, -1);
1444  MPI_Gatherv(local_keys.data(), local_count, MPI_LONG_LONG,
1445  rank == 0 ? all_keys.data() : nullptr, counts.data(), displs.data(), MPI_LONG_LONG,
1446  0, mesh.comm);
1447 
1448  // Build sorted global keys on rank 0
1449  std::vector<std::array<long long, KEY_SIZE>> sorted_keys;
1450  if (rank == 0) {
1451  const int nkeys = total / KEY_SIZE;
1452  sorted_keys.resize(nkeys);
1453  for (int i = 0; i < nkeys; i++) {
1454  std::array<long long, KEY_SIZE> k;
1455  for (int j = 0; j < KEY_SIZE; j++) k[j] = all_keys[i * KEY_SIZE + j];
1456  sorted_keys[i] = k;
1457  }
1458  std::sort(sorted_keys.begin(), sorted_keys.end());
1459  }
1460 
1461  // Broadcast sorted keys
1462  int nkeys = 0;
1463  if (rank == 0) nkeys = (int)sorted_keys.size();
1464  MPI_Bcast(&nkeys, 1, MPI_INT, 0, mesh.comm);
1465  std::vector<long long> flat_sorted(nkeys * KEY_SIZE, -1);
1466  if (rank == 0) {
1467  for (int i = 0; i < nkeys; i++) {
1468  for (int j = 0; j < KEY_SIZE; j++) flat_sorted[i * KEY_SIZE + j] = sorted_keys[i][j];
1469  }
1470  }
1471  MPI_Bcast(flat_sorted.data(), (int)flat_sorted.size(), MPI_LONG_LONG, 0, mesh.comm);
1472 
1473  // Reconstruct sorted_keys on all ranks
1474  if (rank != 0) {
1475  sorted_keys.resize(nkeys);
1476  for (int i = 0; i < nkeys; i++) {
1477  std::array<long long, KEY_SIZE> k;
1478  for (int j = 0; j < KEY_SIZE; j++) k[j] = flat_sorted[i * KEY_SIZE + j];
1479  sorted_keys[i] = k;
1480  }
1481  }
1482 
1483  // Assign deterministic element numbering (both REF and SUBMESH)
1484  SF::vector<mesh_int_t> & nbr_ref = mesh.register_numbering(SF::NBR_ELEM_REF);
1485  SF::vector<mesh_int_t> & nbr_sub = mesh.register_numbering(SF::NBR_ELEM_SUBMESH);
1486  nbr_ref.resize(mesh.l_numelem);
1487  nbr_sub.resize(mesh.l_numelem);
1488  for (size_t i = 0; i < mesh.l_numelem; i++) {
1489  auto k = make_key(i);
1490  auto it = std::lower_bound(sorted_keys.begin(), sorted_keys.end(), k);
1491  if (it == sorted_keys.end() || *it != k) {
1492  log_msg(0, 5, 0, "deterministic numbering failed to find key (rank %d, elem %zu)", rank, i);
1493  EXIT(1);
1494  }
1495  mesh_int_t gid = (mesh_int_t)(it - sorted_keys.begin());
1496  nbr_ref[i] = gid;
1497  nbr_sub[i] = gid;
1498  }
1499 }
1500 
1501 void EMI::setup_output()
1502 {
1503  bool write_binary = false;
1504  std::string output_base = get_basename(param_globals::meshname);
1505  set_dir(OUTPUT);
1506 
1507  // write entire mesh
1508  sf_mesh & mesh = get_mesh(emi_msh);
1509  std::string output_file = output_base + "_e";
1510  log_msg(0, 0, 0, "Writing \"%s\" mesh: %s", mesh.name.c_str(), output_file.c_str());
1511  write_mesh_parallel(mesh, write_binary, output_file.c_str());
1512  // register output for overall phi on the entire mesh
1513  output_manager.register_output(parab_solver.ui, emi_msh, 1, param_globals::phiefile, "mV", NULL);
1514 
1515  // write membrane mesh
1517  mesh_m.name = "Membrane";
1518  output_file = output_base + "_m";
1519  log_msg(0, 0, 0, "Writing \"%s\" mesh: %s", mesh_m.name.c_str(), output_file.c_str());
1520  write_mesh_parallel(mesh_m, write_binary, output_file.c_str());
1521  // register output for Vm on membrane interface
1522  // ensure deterministic element ordering in output
1525  }
1526  output_manager.register_output(parab_solver.vb_unique_face, emi_surface_unique_face_msh, 1, param_globals::vofile, "mV", NULL, true);
1527 
1528  if(param_globals::num_trace) {
1529  sf_mesh & imesh = get_mesh(emi_msh);
1530  open_trace(ion.miif, param_globals::num_trace, param_globals::trace_node, NULL, &imesh);
1531  }
1532 
1533  // initialize generic logger for IO timings per time_dt
1534  IO_stats.init_logger("IO_stats.dat");
1535 }
1536 
1537 void EMI::dump_matrices()
1538 {
1539  std::string bsname = param_globals::dump_basename;
1540  std::string fn;
1541 
1542  set_dir(OUTPUT);
1543 
1544  fn = bsname + "_lhs.bin";
1545  parab_solver.lhs_emi->write(fn.c_str());
1546 
1547  fn = bsname + "_K.bin";
1548  parab_solver.stiffness_emi->write(fn.c_str());
1549 
1550  fn = bsname + "_B.bin";
1551  parab_solver.B->write(fn.c_str());
1552 
1553  fn = bsname + "_Bi.bin";
1554  parab_solver.Bi->write(fn.c_str());
1555 
1556  fn = bsname + "_BsM.bin";
1557  parab_solver.BsM->write(fn.c_str());
1558 
1559  fn = bsname + "_M.bin";
1560  parab_solver.mass_emi->write(fn.c_str());
1561 
1562  fn = bsname + "_Ms.bin";
1563  parab_solver.mass_surf_emi->write(fn.c_str());
1564 
1565 }
1566 
1569 double EMI::timer_val(const int timer_id)
1570 {
1571  // determine
1572  int sidx = stimidx_from_timeridx(stimuli, timer_id);
1573  double val = 0.0;
1574  if(sidx != -1) {
1575  stimuli[sidx].value(val);
1576  }
1577  else
1578  val = std::nan("NaN");
1579 
1580  return val;
1581 }
1582 
1585 std::string EMI::timer_unit(const int timer_id)
1586 {
1587  int sidx = stimidx_from_timeridx(stimuli, timer_id);
1588  std::string s_unit;
1589 
1590  if(sidx != -1)
1591  // found a timer-linked stimulus
1592  s_unit = stimuli[sidx].pulse.wave.f_unit;
1593 
1594  return s_unit;
1595 }
1596 
1597 void EMI::setup_solvers()
1598 {
1599  set_dir(OUTPUT);
1600  parab_solver.init();
1601  parab_solver.rebuild_matrices(mtype_vol, *ion.miif, stimuli, logger);
1602 
1603  if(param_globals::dump2MatLab)
1604  dump_matrices();
1605 }
1606 
1607 void extract_unique_tag(SF::vector<mesh_int_t>& unique_tags)
1608 {
1609  MPI_Comm comm = SF_COMM;
1610  int size, rank;
1611  MPI_Comm_size(comm, &size);
1612  MPI_Comm_rank(comm, &rank);
1613 
1614  binary_sort(unique_tags);
1615  unique_resize(unique_tags); // unique_resize is done locally(for each rank)
1616  make_global(unique_tags, comm);
1617  binary_sort(unique_tags);
1618  unique_resize(unique_tags);
1619 }
1620 
1621 void compute_tags_per_rank(int num_tags, SF::vector<mesh_int_t>& num_tags_per_rank)
1622 {
1623  MPI_Comm comm = SF_COMM;
1624  int size, rank;
1625  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
1626  divide(num_tags, size, num_tags_per_rank);
1627 }
1628 
1629 void EMI::setup_EMI_mesh()
1630 {
1631  log_msg(0,0,0, "\n *** Processing EMI mesh ***\n");
1632 
1633  const std::string basename = param_globals::meshname;
1634  const int verb = param_globals::output_level;
1635  std::map<mesh_t, sf_mesh> & mesh_registry = user_globals::mesh_reg;
1636  assert(mesh_registry.count(emi_msh) == 1);
1637 
1638  set_dir(INPUT);
1639 
1640  sf_mesh & emi_mesh = mesh_registry[emi_msh];
1641  sf_mesh & emi_surfmesh_one_side = mesh_registry[emi_surface_msh];
1642  sf_mesh & emi_surfmesh_w_counter_face = mesh_registry[emi_surface_counter_msh];
1643  sf_mesh & emi_surfmesh_unique_face = mesh_registry[emi_surface_unique_face_msh];
1644 
1645  MPI_Comm comm = emi_mesh.comm;
1646 
1647  int size, rank;
1648  double t1, t2, s1, s2;
1649  const double total_setup_t0 = MPI_Wtime();
1650  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
1651 
1652  //-----------------------------------------------------------------
1653  // Validate mesh dimension: EMI only supports 3D volumetric meshes
1654  //-----------------------------------------------------------------
1655  if (emi_mesh.l_numelem > 0) {
1656  // Check first local element type (all elements should be same type in CARP)
1657  SF::elem_t first_elem = emi_mesh.type[0];
1658 
1659  if (first_elem == SF::Line || first_elem == SF::Tri || first_elem == SF::Quad) {
1660  const char* type_name = (first_elem == SF::Line) ? "1D (Line)" :
1661  (first_elem == SF::Tri) ? "2D (Tri)" :
1662  "2D (Quad)";
1663  if (rank == 0) {
1664  log_msg(0, 5, 0, "\n*** ERROR: EMI model requires a 3D volumetric mesh!");
1665  log_msg(0, 5, 0, "*** Current mesh element type: %s", type_name);
1666  log_msg(0, 5, 0, "*** EMI only supports 3D element types: Tetra, Pyramid, Prism, Hexa");
1667  log_msg(0, 5, 0, "*** Please provide a 3D mesh with volume elements.\n");
1668  }
1669  EXIT(EXIT_FAILURE);
1670  }
1671  }
1672 
1673  //-----------------------------------------------------------------
1674  // Step 1: READ *.intra and *.extra
1675  //-----------------------------------------------------------------
1676  int total_num_tags = 0;
1677  if(verb) log_msg(NULL, 0, 0,"\nReading tags for extra and intra regions from input files");
1678  t1 = MPI_Wtime();
1679  {
1680  SF::vector<mesh_int_t> unique_extra_tags;
1681  SF::vector<mesh_int_t> unique_intra_tags;
1682 
1683  if(verb) log_msg(NULL, 0, 0,"Read extracellular tags");
1684  read_indices_global(unique_extra_tags,basename+".extra", comm);
1685  for(mesh_int_t tag:unique_extra_tags){
1686  parab_solver.extra_tags.insert(tag);
1687  }
1688 
1689  if(verb) log_msg(NULL, 0, 0,"Read intracellular tags");
1690  read_indices_global(unique_intra_tags,basename+".intra", comm);
1691  for(mesh_int_t tag:unique_intra_tags){
1692  parab_solver.intra_tags.insert(tag);
1693  }
1694 
1695  total_num_tags = parab_solver.extra_tags.size() + parab_solver.intra_tags.size();
1696  if(total_num_tags < size){
1697  log_msg(0,5,0, "\nThe number of unique tags on EMI mesh is smaller than number of processors!");
1698  EXIT(EXIT_FAILURE);
1699  }
1700  if(verb) log_msg(NULL, 0, 0,"\nextra_tags=%lu, intra_tags=%lu", parab_solver.extra_tags.size(), parab_solver.intra_tags.size());
1701  }
1702  t2 = MPI_Wtime();
1703  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1704 
1705  //-----------------------------------------------------------------
1706  // Step 2: READ points from *.pts
1707  //-----------------------------------------------------------------
1709  SF::vector<mesh_int_t> ptsidx;
1710  SF::vector<mesh_int_t> ptsData;
1711  if(verb) log_msg(NULL, 0, 0,"\nReading points with data on each vertex");
1712  t1 = MPI_Wtime();
1713  SF::read_points(basename, comm, pts, ptsidx);
1714  ptsData.resize(ptsidx.size());
1715  assert(ptsidx.size()==ptsData.size());
1716  t2 = MPI_Wtime();
1717  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1718 
1719  std::list< sf_mesh* > meshlist;
1720  meshlist.push_back(&emi_mesh);
1721 
1722  //-----------------------------------------------------------------
1723  // Step 3: Distrubute mesh based on tag or *.part
1724  //-----------------------------------------------------------------
1725  if(verb) log_msg(NULL, 0, 0,"\nDistribute mesh based on tags");
1726  // should be replaced by scotch/pt-scotch for efficiency
1727  t1 = MPI_Wtime();
1728  distribute_elements_based_tags(emi_mesh, total_num_tags);
1729  t2 = MPI_Wtime();
1730  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1731 
1732  //-----------------------------------------------------------------
1733  // Step 4: insert coordinates to the mesh
1734  //-----------------------------------------------------------------
1735  // insert points
1736  if(verb) log_msg(NULL, 0, 0, "\nInserting points");
1737  t1 = MPI_Wtime();
1738  SF::insert_points(pts, ptsidx, meshlist);
1739  t2 = MPI_Wtime();
1740  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1741 
1742  //-----------------------------------------------------------------
1743  // Step 5: extract ptsdata based on the location of each vertex
1744  //-----------------------------------------------------------------
1745  if(verb) log_msg(NULL, 0, 0, "\nCompute location of all dofs on emi mesh(inner dofs, memebarne, gap junction) saved into ptsData from original mesh");
1746  t1 = MPI_Wtime();
1747  compute_ptsdata_from_original_mesh( emi_mesh,
1748  SF::NBR_REF,
1749  vertex2ptsdata,
1750  parab_solver.extra_tags,
1751  parab_solver.intra_tags);
1752  t2 = MPI_Wtime();
1753  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1754 
1755  //-----------------------------------------------------------------
1756  // Step 6: extract faces and counterfaces on the interfaces (mem/gap)
1757  // and mark the unique faces on surface mesh
1758  // generate required mapping between surface mesh and unique faces
1759  //-----------------------------------------------------------------
1760  t1 = MPI_Wtime();
1761  if(verb) log_msg(NULL, 0, 0, "\nExtract EMI surface mesh");
1762  hashmap::unordered_map<mesh_int_t, SF::emi_index_rank<mesh_int_t>> unused_map_elem_oneface_to_elem_uniqueFace;
1763  extract_face_based_tags(emi_mesh, SF::NBR_REF, vertex2ptsdata,
1764  parab_solver.line_face,
1765  parab_solver.tri_face,
1766  parab_solver.quad_face,
1767  parab_solver.extra_tags,
1768  parab_solver.intra_tags,
1769  emi_surfmesh_one_side, emi_surfmesh_w_counter_face, emi_surfmesh_unique_face,
1770  parab_solver.map_elem_uniqueFace_to_elem_oneface,
1771  unused_map_elem_oneface_to_elem_uniqueFace);
1772  meshlist.push_back(&emi_surfmesh_one_side);
1773  t2 = MPI_Wtime();
1774  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1775  //-----------------------------------------------------------------
1776  // Step 7: register surface meshes based on the {typeof}_face with counter face
1777  //-----------------------------------------------------------------
1778  compute_surface_mesh_with_counter_face(emi_surfmesh_w_counter_face, SF::NBR_REF,
1779  parab_solver.line_face,
1780  parab_solver.tri_face,
1781  parab_solver.quad_face);
1782 
1783  compute_surface_mesh_with_unique_face(emi_surfmesh_unique_face, SF::NBR_REF,
1784  parab_solver.line_face,
1785  parab_solver.tri_face,
1786  parab_solver.quad_face,
1787  parab_solver.map_elem_uniqueFace_to_tags);
1788 
1789  //-----------------------------------------------------------------
1790  // Step 8: create a map between emi_surfmesh_w_counter_face and emi_surfmesh (both -> one)
1791  //-----------------------------------------------------------------
1792  SF::create_reverse_elem_mapping_between_surface_meshes(parab_solver.line_face,
1793  parab_solver.tri_face,
1794  parab_solver.quad_face,
1795  parab_solver.vec_both_to_one_face,
1796  comm);
1797  t2 = MPI_Wtime();
1798  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1799 
1800  //-----------------------------------------------------------------
1801  // Step 9: global number of interfaces only on emi_surfmesh without counter face
1802  //-----------------------------------------------------------------
1803  if(verb) log_msg(NULL, 0, 0, "\ncompute global number of interface");
1804  size_t global_count_surf = 0;
1805  size_t numelem_surface = emi_surfmesh_one_side.l_numelem;
1806  unsigned long long local_count_surf = static_cast<unsigned long long>(numelem_surface);
1807  unsigned long long global_count_surf_ull = 0;
1808  MPI_Reduce(&local_count_surf, &global_count_surf_ull, 1, MPI_UNSIGNED_LONG_LONG, MPI_SUM, 0, MPI_COMM_WORLD);
1809  global_count_surf = static_cast<size_t>(global_count_surf_ull);
1810  if(verb && rank==0) fprintf(stdout, "global number of interfaces = %zu\n", global_count_surf);
1811 
1812  //-----------------------------------------------------------------
1813  // Step 10: submesh_numbering on the emi_mesh and generate parallel layout
1814  //-----------------------------------------------------------------
1815  {
1817  sub_numbering(emi_mesh);
1818  emi_mesh.generate_par_layout();
1819  }
1820 
1821  SF::meshdata<mesh_int_t, mesh_real_t> tmesh_backup_old = emi_mesh;
1822 
1823  //-----------------------------------------------------------------
1824  // Step 11: make a map bwteen key(vertex,tag)-> value(dof) after decoupling interfaces defined on emi mesh
1825  //-----------------------------------------------------------------
1826  // During interface decoupling and the introduction of new degrees of freedom (DoFs),
1827  // we assume that the mesh partitioning is performed at least on a per-tag basis.
1828  // This means that all elements sharing the same tag number belong to the same rank.
1829  t1 = MPI_Wtime();
1830  if(verb) log_msg(NULL, 0, 0, "\ndecouple emi interfaces");
1831  if(verb) log_msg(NULL, 0, 0, "\tcompute map oldIdx to dof");
1832  compute_map_vertex_to_dof(emi_mesh, SF::NBR_REF, vertex2ptsdata, parab_solver.extra_tags, parab_solver.map_vertex_tag_to_dof);
1833 
1834  //-----------------------------------------------------------------
1835  // Step 12: complete the map bwteen key(vertex,tag)-> value(dof) for counter faces defined on emi mesh
1836  //-----------------------------------------------------------------
1837  if(verb) log_msg(NULL, 0, 0, "\tcomplete map oldIdx to dof with counter interface");
1838  // add to the map the counter part of the interface
1839  complete_map_vertex_to_dof_with_counter_face(parab_solver.line_face, parab_solver.tri_face, parab_solver.quad_face, parab_solver.map_vertex_tag_to_dof);
1840 
1841  //-----------------------------------------------------------------
1842  // Step 13: update emi mesh with new dofs, so decoupling the interfaces applied on emi mesh
1843  //-----------------------------------------------------------------
1844  if(verb) log_msg(NULL, 0, 0, "\tupdate mesh with dof");
1845  update_emi_mesh_with_dofs(emi_mesh, SF::NBR_REF, parab_solver.map_vertex_tag_to_dof, parab_solver.dof2vertex);
1846  t2 = MPI_Wtime();
1847  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1848 
1849  //-----------------------------------------------------------------
1850  // Step 14: set map_vertex_tag_to_dof_petsc where the key(vertex,tag) and value(dof,petsc)
1851  //-----------------------------------------------------------------
1852  t1 = MPI_Wtime();
1853  if(verb) log_msg(NULL, 0, 0, "\nInitialize petsc =0 for map<oldIdx,tag> -><dof, petsc>");
1854  // Iterate over map to assign the (oldIndx, tag) -> (dof, petsc) atm petsc = 0
1855  for(const auto & key_value : parab_solver.map_vertex_tag_to_dof)
1856  {
1857  mesh_int_t gIndex_old = key_value.first.first;
1858  mesh_int_t tag_old = key_value.first.second;
1859  mesh_int_t dof = key_value.second;
1860 
1861  std::pair <mesh_int_t,mesh_int_t> dof_petsc = std::make_pair(dof,-1); // petsc numbering ...
1862  parab_solver.map_vertex_tag_to_dof_petsc.insert({key_value.first,dof_petsc});
1863  }
1864  t2 = MPI_Wtime();
1865  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1866 
1867  //-----------------------------------------------------------------
1868  // Step 16: insert the coordinates on new dofs to the emi mesh
1869  //-----------------------------------------------------------------
1870  t1 = MPI_Wtime();
1871  if(verb) log_msg(NULL, 0, 0, "Inserting points and ptsData of dofs to emi_mesh");
1872  insert_points_ptsData_to_dof(tmesh_backup_old, emi_mesh, SF::NBR_REF, parab_solver.dof2vertex, vertex2ptsdata, dof2ptsData);
1873  t2 = MPI_Wtime();
1874  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1875 
1876  //-----------------------------------------------------------------
1877  // Step 17: submesh_numbering on the emi_mesh and generate parallel layout after decoupling the interfaces on the emi mesh
1878  //-----------------------------------------------------------------
1879  t1 = MPI_Wtime();
1880  if(verb) log_msg(NULL, 0, 0, "Generating unique PETSc numberings");
1881  {
1883  sub_numbering(emi_mesh);
1884  emi_mesh.generate_par_layout();
1885  }
1886  t2 = MPI_Wtime();
1887  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1888 
1889  //-----------------------------------------------------------------
1890  // Step 18: register petsc_numbering on emi mesh
1891  //-----------------------------------------------------------------
1892  t1 = MPI_Wtime();
1893  if(verb) log_msg(NULL, 0, 0, "Generating unique PETSc numberings");
1894  {
1895  SF::petsc_numbering<mesh_int_t,mesh_real_t> petsc_numbering(emi_mesh.pl);
1896  petsc_numbering(emi_mesh);
1897  }
1898  t2 = MPI_Wtime();
1899  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1900 
1901  //-----------------------------------------------------------------
1902  // Step 19: complete map_vertex_tag_to_dof_petsc where the key(vertex,tag) and value(dof,petsc) for counter faces on emi mesh
1903  //-----------------------------------------------------------------
1904  if(verb) log_msg(NULL, 0, 0, "Updating the map between indices to PETSc numberings");
1905  t1 = MPI_Wtime();
1906  update_map_indices_to_petsc(emi_mesh, SF::NBR_REF, SF::NBR_PETSC, parab_solver.extra_tags, parab_solver.map_vertex_tag_to_dof_petsc, parab_solver.dof2vertex, parab_solver.elemTag_emi_mesh);
1907  t2 = MPI_Wtime();
1908  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1909 
1910  //-----------------------------------------------------------------
1911  // Step 20: update value of map {type_of}_face with new dofs after decoupling
1912  //-----------------------------------------------------------------
1913  t1 = MPI_Wtime();
1914  if(verb) log_msg(NULL, 0, 0, "Updating surface mesh with dof");
1915  update_faces_on_surface_mesh_after_decoupling_with_dofs(emi_surfmesh_one_side, parab_solver.map_vertex_tag_to_dof, parab_solver.line_face, parab_solver.tri_face, parab_solver.quad_face);
1916  t2 = MPI_Wtime();
1917  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1918 
1919  //-----------------------------------------------------------------
1920  // Step 21: register emi_surfmesh for SF::NBR_ELEM_REF, NBR_REF, SF::NBR_SUBMESH
1921  //-----------------------------------------------------------------
1922  t1 = MPI_Wtime();
1923  if(verb) log_msg(NULL, 0, 0, "Layout for element of EMI surfmesh");
1924  SF::vector<mesh_int_t> & emi_surfmesh_elem = emi_surfmesh_one_side.register_numbering(SF::NBR_ELEM_REF);
1925  SF::vector<long int> layout;
1926  SF::layout_from_count<long int>(emi_surfmesh_one_side.l_numelem, layout, emi_surfmesh_one_side.comm);
1927  size_t count = layout[rank+1] - layout[rank];
1928  emi_surfmesh_elem.resize(count);
1929  for (int i = 0; i < count; ++i){
1930  emi_surfmesh_elem[i] = layout[rank]+i;
1931  }
1932 
1933  emi_surfmesh_one_side.localize(SF::NBR_REF);
1934  emi_surfmesh_one_side.register_numbering(SF::NBR_SUBMESH);
1935  t2 = MPI_Wtime();
1936  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1937 
1938  //-----------------------------------------------------------------
1939  // Step 22: register emi_surfmesh_w_counter_face for SF::NBR_ELEM_REF, NBR_REF, SF::NBR_SUBMESH
1940  //-----------------------------------------------------------------
1941  t1 = MPI_Wtime();
1942  if(verb) log_msg(NULL, 0, 0, "Layout for element of EMI surfmesh w counter face");
1943  SF::vector<mesh_int_t> & emi_surfmesh_counter_elem = emi_surfmesh_w_counter_face.register_numbering(SF::NBR_ELEM_REF);
1944  SF::vector<long int> layout_counter;
1945  SF::layout_from_count<long int>(emi_surfmesh_w_counter_face.l_numelem, layout_counter, emi_surfmesh_w_counter_face.comm);
1946  size_t count_counter = layout_counter[rank+1] - layout_counter[rank];
1947  emi_surfmesh_counter_elem.resize(count_counter);
1948  for (int i = 0; i < count_counter; ++i){
1949  emi_surfmesh_counter_elem[i] = layout_counter[rank]+i;
1950  }
1951  emi_surfmesh_w_counter_face.localize(SF::NBR_REF);
1952  emi_surfmesh_w_counter_face.register_numbering(SF::NBR_SUBMESH);
1953  t2 = MPI_Wtime();
1954  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1955 
1956 
1957  //-----------------------------------------------------------------
1958  // Step 23: register emi_surfmesh_w_counter_face for SF::NBR_ELEM_REF, NBR_REF, SF::NBR_SUBMESH
1959  //-----------------------------------------------------------------
1960  t1 = MPI_Wtime();
1961  if(verb) log_msg(NULL, 0, 0, "Layout for element of EMI surfmesh w counter face");
1962  SF::vector<mesh_int_t> & emi_surfmesh_unique_elem = emi_surfmesh_unique_face.register_numbering(SF::NBR_ELEM_REF);
1963  SF::vector<long int> layout_unique;
1964  SF::layout_from_count<long int>(emi_surfmesh_unique_face.l_numelem, layout_unique, emi_surfmesh_unique_face.comm);
1965  size_t count_unique = layout_unique[rank+1] - layout_unique[rank];
1966  emi_surfmesh_unique_elem.resize(count_unique);
1967  for (int i = 0; i < count_unique; ++i){
1968  emi_surfmesh_unique_elem[i] = layout_unique[rank]+i;
1969  }
1970  emi_surfmesh_unique_face.localize(SF::NBR_REF);
1971  emi_surfmesh_unique_face.register_numbering(SF::NBR_SUBMESH);
1972  t2 = MPI_Wtime();
1973  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1974 
1975  //-----------------------------------------------------------------
1976  // Step 24: insert coorinates to emi_surfmesh_one_side and emi_surfmesh_w_counter_face
1977  //-----------------------------------------------------------------
1978  t1 = MPI_Wtime();
1979  if(verb) log_msg(NULL, 0, 0, "Inserting points to EMI surfmesh");
1980  insert_points_to_surface_mesh(tmesh_backup_old, emi_surfmesh_one_side, SF::NBR_REF, parab_solver.dof2vertex, parab_solver.extra_tags, parab_solver.elemTag_surface_mesh);
1981  insert_points_to_surface_mesh(tmesh_backup_old, emi_surfmesh_w_counter_face, SF::NBR_REF, parab_solver.dof2vertex, parab_solver.extra_tags, parab_solver.elemTag_surface_w_counter_mesh);
1982  insert_points_to_surface_mesh(tmesh_backup_old, emi_surfmesh_unique_face, SF::NBR_REF, parab_solver.dof2vertex);
1983  t2 = MPI_Wtime();
1984  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
1985 
1986  //-----------------------------------------------------------------
1987  // Step 25: submesh_numbering & petc numebring for both emi_surfmesh_one_side and emi_surfmesh_w_counter_face
1988  //-----------------------------------------------------------------
1989  t1 = MPI_Wtime();
1990  if(verb) log_msg(NULL, 0, 0, "Generating submesh_numbering and PETSc numberings for surface mesh");
1991  {
1993  sub_numbering(emi_surfmesh_one_side);
1994  emi_surfmesh_one_side.generate_par_layout();
1995 
1996  SF::petsc_numbering<mesh_int_t,mesh_real_t> petsc_numbering(emi_surfmesh_one_side.pl);
1997  petsc_numbering(emi_surfmesh_one_side);
1998  }
1999 
2000  {
2002  sub_numbering(emi_surfmesh_w_counter_face);
2003  emi_surfmesh_w_counter_face.generate_par_layout();
2004 
2005  SF::petsc_numbering<mesh_int_t,mesh_real_t> petsc_numbering(emi_surfmesh_w_counter_face.pl);
2006  petsc_numbering(emi_surfmesh_w_counter_face);
2007  }
2008 
2009  {
2011  sub_numbering(emi_surfmesh_unique_face);
2012  emi_surfmesh_unique_face.generate_par_layout();
2013 
2014  SF::petsc_numbering<mesh_int_t,mesh_real_t> petsc_numbering(emi_surfmesh_unique_face.pl);
2015  petsc_numbering(emi_surfmesh_unique_face);
2016  // Ensure deterministic element ordering across MPI for output
2017  assign_deterministic_elem_numbering(emi_surfmesh_unique_face);
2018  }
2019  t2 = MPI_Wtime();
2020  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
2021 
2022  //-----------------------------------------------------------------
2023  // Step 26: assign petsc numbering to map_vertex_tag_to_dof_petsc where the key(vertex,tag) and value(dof,petsc)
2024  //-----------------------------------------------------------------
2025  t1 = MPI_Wtime();
2026  if(verb) log_msg(NULL, 0, 0, "assign PETSc numbering for new faces");
2027  SF::assign_petsc_on_counter_face(parab_solver.map_vertex_tag_to_dof_petsc,comm);
2028  {
2029  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>, std::pair<mesh_int_t,mesh_int_t>>::iterator it;
2030  for (it = parab_solver.map_vertex_tag_to_dof_petsc.begin(); it != parab_solver.map_vertex_tag_to_dof_petsc.end(); it++)
2031  {
2032  std::pair <mesh_int_t,mesh_int_t> Index_tag_old = it->first;
2033  std::pair <mesh_int_t,mesh_int_t> dof_petsc = it->second;
2034  parab_solver.dof2petsc[dof_petsc.first] = dof_petsc.second;
2035  parab_solver.petsc2dof[dof_petsc.second] = dof_petsc.first;
2036  }
2037  }
2038 
2039  // DEBUG: Check for invalid PETSc indices after exchange
2040  #ifdef EMI_DEBUG_MESH
2041  {
2042  int invalid_count = 0;
2043  for (const auto& [key, val] : parab_solver.map_vertex_tag_to_dof_petsc) {
2044  if (val.second < 0) {
2045  invalid_count++;
2046  if (invalid_count <= 3) {
2047  fprintf(stderr, "RANK %d INVALID: vertex=%ld tag=%ld dof=%ld petsc=%ld\n",
2048  rank, (long)key.first, (long)key.second, (long)val.first, (long)val.second);
2049  }
2050  }
2051  }
2052  fprintf(stderr, "RANK %d: After Step 26: %d invalid PETSc indices out of %zu total\n",
2053  rank, invalid_count, parab_solver.map_vertex_tag_to_dof_petsc.size());
2054  fflush(stderr);
2055  }
2056  #endif
2057 
2058  t2 = MPI_Wtime();
2059  if(verb) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
2060 
2061  //-----------------------------------------------------------------
2062  // Step 27: update counter faces to the map{type_of}_face
2063  //-----------------------------------------------------------------
2064  added_counter_faces_to_map(parab_solver.line_face, parab_solver.tri_face, parab_solver.quad_face);
2065  const double total_setup = MPI_Wtime() - total_setup_t0;
2066  log_msg(0,0,0, "Total setup_EMI_mesh processing time: %.5f sec.", float(total_setup));
2067 
2068  log_msg(0,0,0, "\n *** EMI mesh processing Done ***\n");
2069 }
2070 
2071 void distribute_elements_based_tags(SF::meshdata<mesh_int_t, mesh_real_t>& mesh,
2072  int total_num_tags)
2073 {
2074  MPI_Comm comm = mesh.comm;
2075  int size, rank;
2076  double t1, t2, s1, s2;
2077  MPI_Comm_size(comm, &size); MPI_Comm_rank(comm, &rank);
2078  const int verb = param_globals::output_level;
2079 
2080  if(total_num_tags < size)
2081  {
2082  PetscPrintf(PETSC_COMM_WORLD,"\nThe number of processors should be less than the number of tags, size = %d & ntags = %d !!!\n",
2083  size, total_num_tags);
2084  cleanup_and_exit();
2085  }
2086 
2087  if(verb==10) log_msg(NULL, 0, 0,"\ncompute the number of tags which belongs to one rank");
2088  // compute the number of tags per rank
2089  SF::vector<mesh_int_t> ntags_per_rank;
2090  t1 = MPI_Wtime();
2091  compute_tags_per_rank(total_num_tags, ntags_per_rank);
2092  t2 = MPI_Wtime();
2093  if(verb==10) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
2094 
2095  // generate global destination w.r.t. tags
2096  // redistribute elements based on new partition
2097  t1 = MPI_Wtime();
2098  SF::vector<mesh_int_t> part_based_Tags(mesh.l_numelem);
2099  partition_based_tags(total_num_tags, mesh.tag, ntags_per_rank, part_based_Tags);
2100  SF::redistribute_elements(mesh,part_based_Tags);
2101  // permute elements locally first based on the tag then element index
2102  permute_mesh_locally_based_on_tag_elemIdx(mesh);
2103  t2 = MPI_Wtime();
2104  if(verb==10) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
2105 }
2106 
2107 // In the EMI model we cannot have partition boundaries crossing cells.
2108 // Therefore this function assigns a partition to each tag. All elements with this tag will go to the same partition.
2109 void partition_based_tags(int num_tags,
2111  SF::vector<mesh_int_t> num_tags_per_rank,
2112  SF::vector<mesh_int_t>& part)
2113 {
2114  const int verb = param_globals::output_level;
2115  MPI_Comm comm = SF_COMM;
2116  int size, rank;
2117  MPI_Comm_size(comm, &size);
2118  MPI_Comm_rank(comm, &rank);
2119 
2120  // we need to have tags_to_rank as a map globally.
2122  // Try to load the mapping from a .part file and fall back on the internal mapping if there is no such file.
2123  // If the mapping is read from file, the num_tags_per_rank argument is not used.
2124  if (!load_partitions_from_file(tags_to_rank_map, num_tags, comm)) {
2125  if(verb==10) log_msg(NULL, 0, 0,"\ncompute the number of unique tags");
2126  SF::vector<mesh_int_t> unique_tags = tag;
2127  double t1 = MPI_Wtime();
2128  extract_unique_tag(unique_tags);
2129  double t2 = MPI_Wtime();
2130  if(verb==10) log_msg(NULL, 0, 0, "Done in %f sec.", float(t2 - t1));
2131 
2132  if (unique_tags.size() != static_cast<size_t>(num_tags)) {
2133  log_msg(0,5,0, "\nerror: the number of tags in the EMI mesh does not match the total tags from *.extra and *.intra");
2134  EXIT(1);
2135  }
2136  map_tags_to_rank(size, unique_tags, num_tags_per_rank, tags_to_rank_map);
2137  }
2138 
2139  for (size_t i = 0; i < part.size(); ++i) {
2140  if(tags_to_rank_map.count(tag[i]))
2141  part[i] = tags_to_rank_map[tag[i]];
2142  }
2143 }
2144 
2145 void map_tags_to_rank(int size, const SF::vector<mesh_int_t> & unique_tags, const SF::vector<mesh_int_t> & num_tags_per_rank, hashmap::unordered_map<mesh_int_t, mesh_int_t> &tags_to_rank_map)
2146 {
2147  // the goal is to map the tag as a key -> rank, should be global one
2148  mesh_int_t t = 0; // tag index
2149  for (size_t r = 0; r < size; ++r) {
2150  for (size_t count = 0; count < num_tags_per_rank[r];) {
2151  mesh_int_t tag = unique_tags[t];
2152 
2153  // // here we consider each tag belongs to only one rank.
2154  tags_to_rank_map.insert({tag, r});
2155  count +=1;
2156  t+=1;
2157  }
2158  }
2159 }
2160 
2161 // This function creates a tag-to-rank mapping based only on the numerical order of the tags.
2162 bool load_partitions_from_file(hashmap::unordered_map<mesh_int_t, mesh_int_t>& tags_to_rank_map,
2163  int expected_num_tags,
2164  MPI_Comm comm)
2165 {
2166  FILE_SPEC tags_to_part = nullptr;
2167 
2168  SF::vector<int> parts;
2169  const std::string basename = param_globals::meshname;
2170  FILE* fd = fopen((basename + ".part").c_str(), "r");
2171  if (fd != NULL) { // check if the file exists otherwise read_indices_global() crash
2172  read_indices_global(parts, basename + ".part", comm);
2173  fclose(fd);
2174  }
2175 
2176  if (parts.size() == 0) return false; // file does not exist or is empty
2177 
2178  if (parts.size() % 2 != 0) {
2179  log_msg(0,5,0, "\nThe part file should contain 2 lines per tag, one for the tag number and the next for its associated partition number.!");
2180  EXIT(1);
2181  }
2182 
2183  for(int i = 0; i < parts.size(); i+=2) {
2184  tags_to_rank_map.insert({parts[i], parts[i + 1]});
2185  }
2186 
2187  if (tags_to_rank_map.size() != static_cast<size_t>(expected_num_tags)) {
2188  log_msg(0,5,0, "\nerror: the number of tags in the .part file does not match the total tags from *.extra and *.intra");
2189  EXIT(1);
2190  }
2191 
2192  return true;
2193 }
2194 
2195 void permute_mesh_locally_based_on_tag_elemIdx(SF::meshdata<mesh_int_t, mesh_real_t>& mesh)
2196 {
2197  mesh.globalize(SF::NBR_REF);
2200  interval(perm, 0, mesh.tag.size());
2201 
2202  SF::vector<mesh_int_t> tags = tmesh.tag;
2204  binary_sort_sort_copy(tags, elemIdx, perm);
2205  permute_mesh(tmesh, mesh, perm);
2206 
2207  mesh.localize(SF::NBR_REF);
2208 }
2209 
2210 } // namespace opencarp
2211 
2212 #endif
int mesh_int_t
Definition: SF_container.h:46
#define SF_COMM
the default SlimFem MPI communicator
Definition: SF_globals.h:27
double SF_real
Use the general double as real type.
Definition: SF_globals.h:38
#define MAX_LOG_LEVEL
Definition: basics.h:331
#define ECHO
Definition: basics.h:324
#define CALI_CXX_MARK_FUNCTION
Definition: caliper_hooks.h:5
#define CALI_MARK_BEGIN(_str)
Definition: caliper_hooks.h:3
#define CALI_MARK_END(_str)
Definition: caliper_hooks.h:4
void globalize(SF_nbr nbr_type)
Localize the connectivity data w.r.t. a given numbering.
Definition: SF_container.h:525
size_t l_numelem
local number of elements
Definition: SF_container.h:399
void localize(SF_nbr nbr_type)
Localize the connectivity data w.r.t. a given numbering.
Definition: SF_container.h:496
MPI_Comm comm
the parallel mesh is defined on a MPI world
Definition: SF_container.h:404
vector< T > & get_numbering(SF_nbr nbr_type)
Get the vector defining a certain numbering.
Definition: SF_container.h:464
vector< T > tag
element tag
Definition: SF_container.h:417
Functor class generating a numbering optimized for PETSc.
Definition: SF_numbering.h:231
Container for a PETSc VecScatter.
Functor class applying a submesh renumbering.
Definition: SF_numbering.h:68
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
hm_int count(const K &key) const
Check if key exists.
Definition: hashmap.hpp:623
void insert(InputIterator first, InputIterator last)
Insert Iterator range.
Definition: hashmap.hpp:583
size_t size() const
Definition: hashmap.hpp:731
size_t size() const
Definition: hashmap.hpp:1152
iterator find(const K &key)
Definition: hashmap.hpp:1092
hm_int erase(const K &key)
Definition: hashmap.hpp:1064
long d_time
current time instance index
Definition: timer_utils.h:77
double time_step
global reference time step
Definition: timer_utils.h:78
int add_eq_timer(double istart, double iend, int ntrig, double iintv, double idur, const char *iname, const char *poolname=nullptr)
Add a equidistant step timer to the array of timers.
Definition: timer_utils.cc:78
long d_end
final index in multiples of dt
Definition: timer_utils.h:82
double time
current time
Definition: timer_utils.h:76
EMI model based on computed current on the faces, main EMI physics class.
void init_solver(SF::abstract_linear_solver< T, S > **sol)
Definition: SF_init.h:220
void read_points(const std::string basename, const MPI_Comm comm, vector< S > &pts, vector< T > &ptsidx)
Read the points and insert them into a list of meshes.
Definition: SF_mesh_io.h:824
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 make_global(const vector< T > &vec, vector< T > &out, MPI_Comm comm)
make a parallel vector global
Definition: SF_network.h:225
void permute_mesh(const meshdata< T, S > &inmesh, meshdata< T, S > &outmesh, const vector< T > &perm)
Permute the element data of a mesh based on a given permutation.
Definition: SF_mesh_utils.h:56
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 divide(const size_t gsize, const size_t num_parts, vector< T > &loc_sizes)
divide gsize into num_parts local parts with even distribution of the remainder
Definition: SF_vector.h:358
void count(const vector< T > &data, vector< S > &cnt)
Count number of occurrences of indices.
Definition: SF_vector.h:332
void insert_points(const vector< S > &pts, const vector< T > &ptsidx, std::list< meshdata< T, S > * > &meshlist)
Insert the points from the read-in buffers into a list of distributed meshes.
Definition: SF_mesh_io.h:894
void assemble_matrix(abstract_matrix< T, S > &mat, meshdata< mesh_int_t, mesh_real_t > &domain, matrix_integrator< mesh_int_t, mesh_real_t > &integrator)
Generalized matrix assembly.
int max_nodal_edgecount(const meshdata< T, S > &mesh)
Compute the maximum number of node-to-node edges for a mesh.
Definition: SF_container.h:608
void redistribute_elements(meshdata< T, S > &mesh, meshdata< T, S > &sendbuff, vector< T > &part)
Redistribute the element data of a parallel mesh among the ranks based on a partitioning.
void restrict_to_set(vector< T > &v, const hashmap::unordered_set< T > &set)
size_t root_count_ascii_lines(std::string file, MPI_Comm comm)
count the lines in a ascii file
void init_vector(SF::abstract_vector< T, S > **vec)
Definition: SF_init.h:99
void binary_sort(vector< T > &_V)
Definition: SF_sort.h:284
void init_matrix(SF::abstract_matrix< T, S > **mat)
Definition: SF_init.h:199
void write_mesh_parallel(const meshdata< T, S > &mesh, bool binary, std::string basename)
Write a parallel mesh to harddisk without gathering it on one rank.
elem_t
element type enum
Definition: SF_container.h:53
@ Line
Definition: SF_container.h:61
@ Tri
Definition: SF_container.h:60
@ Quad
Definition: SF_container.h:59
void binary_sort_sort_copy(vector< T > &_V, vector< T > &_W, vector< S > &_A)
Definition: SF_sort.h:335
@ NBR_PETSC
PETSc numbering of nodes.
Definition: SF_container.h:203
@ 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
@ NBR_SUBMESH
Submesh nodal numbering: The globally ascending sorted reference indices are reindexed.
Definition: SF_container.h:202
@ NBR_ELEM_SUBMESH
Submesh element numbering: The globally ascending sorted reference indices are reindexed.
Definition: SF_container.h:205
constexpr T min(T a, T b)
Definition: ion_type.h:33
constexpr T max(T a, T b)
Definition: ion_type.h:31
void dump_trace(MULTI_IF *MIIF, limpet::Real time)
void open_trace(MULTI_IF *MIIF, int n_traceNodes, int *traceNodes, int *label, opencarp::sf_mesh *imesh)
Set up ionic model traces at some global node numbers.
timer_manager * tm_manager
a manager for the various physics timers
Definition: main.cc:58
bool using_legacy_stimuli
flag storing whether legacy stimuli are used
Definition: main.cc:64
std::map< mesh_t, sf_mesh > mesh_reg
Registry for the different meshes used in a multi-physics simulation.
Definition: main.cc:52
int stimidx_from_timeridx(const SF::vector< stimulus > &stimuli, const int timer_id)
determine link between timer and stimulus
Definition: electrics.cc:807
@ iotm_console
Definition: timer_utils.h:44
@ iotm_trace
Definition: timer_utils.h:44
sf_vec * get_data(datavec_t d)
Retrieve a petsc data vector from the data registry.
Definition: sim_utils.cc:925
SF::scattering * get_scattering(const int from, const int to, const SF::SF_nbr nbr, const int dpn)
Get a scattering from the global scatter registry.
void set_cond_type(MaterialType &m, cond_t type)
Definition: electrics.cc:859
SF_real get_volume_from_nodes(sf_mat &mass, SF::vector< mesh_int_t > &local_idx)
Definition: fem_utils.cc:217
sf_mesh & get_mesh(const mesh_t gt)
Get a mesh by specifying the gridID.
Definition: sf_interface.cc:33
SF::scattering * register_scattering(const int from, const int to, const SF::SF_nbr nbr, const int dpn)
Register a scattering between to grids, or between algebraic and nodal representation of data on the ...
Definition: sf_interface.cc:69
cond_t
description of electrical tissue properties
Definition: fem_types.h:42
@ intra_cond
Definition: fem_types.h:43
SF::scattering * get_permutation(const int mesh_id, const int perm_id, const int dpn)
Get the PETSC to canonical permutation scattering for a given mesh and number of dpn.
SF::meshdata< mesh_int_t, mesh_real_t > sf_mesh
Definition: sf_interface.h:48
void apply_stim_to_vector(const stimulus &s, sf_vec &vec, bool add)
Definition: electrics.cc:478
int set_dir(IO_t dest)
Definition: sim_utils.cc:495
void cleanup_and_exit()
Definition: sim_utils.cc:1312
void read_indices_global(SF::vector< T > &idx, const std::string filename, MPI_Comm comm)
Definition: fem_utils.h:53
int get_rank(MPI_Comm comm=PETSC_COMM_WORLD)
Definition: basics.h:292
void indices_from_region_tags(SF::vector< mesh_int_t > &idx, const sf_mesh &mesh, const hashmap::unordered_set< int > &tags)
Populate vertex data with the vertices of multiple tag regions.
Definition: fem_utils.cc:169
const char * get_mesh_type_name(mesh_t t)
get a char* to the name of a mesh type
Definition: sf_interface.cc:46
void region_mask(mesh_t meshspec, SF::vector< RegionSpecs > &regspec, SF::vector< int > &regionIDs, bool mask_elem, const char *reglist)
classify elements/points as belonging to a region
Definition: ionics.cc:363
void init_stim_info(void)
uses potential for stimulation
Definition: stimulate.cc:49
bool is_extra(stim_t type)
whether stimulus is on extra grid (or on intra)
Definition: stimulate.cc:82
FILE_SPEC f_open(const char *fname, const char *mode)
Open a FILE_SPEC.
Definition: basics.cc:138
bool have_dbc_stims(const SF::vector< stimulus > &stimuli)
return wheter any stimuli require dirichlet boundary conditions
Definition: electrics.cc:882
SF::scattering * register_permutation(const int mesh_id, const int perm_id, const int dpn)
Register a permutation between two orderings for a mesh.
bool is_current(stim_t type)
uses current as stimulation
Definition: stimulate.cc:72
@ OUTPUT
Definition: sim_utils.h:55
char * dupstr(const char *old_str)
Definition: basics.cc:44
void log_msg(FILE_SPEC out, int level, unsigned char flag, const char *fmt,...)
Definition: basics.cc:72
mesh_t
The enum identifying the different meshes we might want to load.
Definition: sf_interface.h:59
@ emi_surface_unique_face_msh
Definition: sf_interface.h:68
@ emi_surface_msh
Definition: sf_interface.h:66
@ emi_surface_counter_msh
Definition: sf_interface.h:67
void get_time(double &tm)
Definition: basics.h:452
bool mesh_is_registered(const mesh_t gt)
check wheter a SF mesh is set
Definition: sf_interface.cc:63
SF::abstract_vector< SF_int, SF_real > sf_vec
Definition: sf_interface.h:50
int get_size(MPI_Comm comm=PETSC_COMM_WORLD)
Definition: basics.h:306
SF::abstract_matrix< SF_int, SF_real > sf_mat
Definition: sf_interface.h:52
V timing(V &t2, const V &t1)
Definition: basics.h:464
std::string get_basename(const std::string &path)
Definition: basics.cc:61
void f_close(FILE_SPEC &f)
Close a FILE_SPEC.
Definition: basics.cc:165
@ ElecMat
Definition: fem_types.h:39
file_desc * FILE_SPEC
Definition: basics.h:138
Basic physics types.
#define UM2_to_CM2
convert um^2 to cm^2
Definition: physics_types.h:35
#define PETSC_TO_CANONICAL
Permute algebraic data from PETSC to canonical ordering.
Definition: sf_interface.h:79
#define ALG_TO_NODAL
Scatter algebraic to nodal.
Definition: sf_interface.h:77
#define ELEM_PETSC_TO_CANONICAL
Permute algebraic element data from PETSC to canonical ordering.
Definition: sf_interface.h:81
#define EXP_POSTPROCESS
Definition: sim_utils.h:172
Electrical stimulation functions.