openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
ionicsOnFace.cc
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 // openCARP is an open cardiac electrophysiology simulator.
3 //
4 // Copyright (C) 2020 openCARP project
5 //
6 // This program is licensed under the openCARP Academic Public License (APL)
7 // v1.0: You can use and redistribute it and/or modify it in non-commercial
8 // academic environments under the terms of APL as published by the openCARP
9 // project v1.0, or (at your option) any later version. Commercial use requires
10 // a commercial license (info@opencarp.org).
11 //
12 // This program is distributed without any warranty; see the openCARP APL for
13 // more details.
14 //
15 // You should have received a copy of the openCARP APL along with this program
16 // and can find it online: http://www.opencarp.org/license
17 // ----------------------------------------------------------------------------
18 
27 #if WITH_EMI_MODEL
28 #include "ionicsOnFace.h"
29 
30 #include "SF_init.h"
31 
32 namespace opencarp {
33 
34 void initialize_sv_dumps_onFace(limpet::MULTI_IF *pmiif, IMPregion_EMI* reg, int id, double t, double dump_dt);
35 
36 void IonicsOnFace::compute_step()
37 {
38  double t1, t2;
39  get_time(t1);
40 
41  miif->compute_ionic_current();
42 
43  double comp_time = timing(t2, t1);
44  this->compute_time += comp_time;
45 
46  comp_stats.calls++;
47  comp_stats.tot_time += comp_time;
48 
50  comp_stats.log_stats(user_globals::tm_manager->time, false);
51 }
52 
53 void IonicsOnFace::destroy()
54 {
55  miif->free_MIIF();
56 }
57 
58 void IonicsOnFace::output_step()
59 {}
60 
61 void IonicsOnFace::initialize()
62 {
63  double t1, t2;
64  get_time(t1);
65 
66  set_dir(OUTPUT);
67 
68  // initialize generic logger for ODE timings per time_dt
69  comp_stats.init_logger("ODE_stats.dat");
70 
71  double tstart = 0.;
72  sf_mesh & mesh = get_mesh(ion_domain);
73  int loc_size = mesh.l_numelem;
74 
75  // create ionic current vector and register it
76  sf_vec *IIon;
77  sf_vec *Vmv;
78  SF::init_vector(&IIon, mesh, 1, sf_vec::elemwise);
79  SF::init_vector(&Vmv, mesh, 1, sf_vec::elemwise);
82 
83  // setup miif
84  miif = new limpet::MULTI_IF();
85  // store IIF_IDs and Plugins in arrays
86  miif->name = "myocardium";
87  miif->gdata[limpet::Vm] = Vmv;
88  miif->gdata[limpet::Iion] = IIon;
89 
90  SF::vector<RegionSpecs_EMI> rs(param_globals::num_imp_regions);
91 
92  // remove user‑specified pairs from the default sets
93  for (size_t i = 2; i < rs.size(); i++ ) {
94  for (int j = 0; j < param_globals::imp_region_emi[i].num_IDs; ++j)
95  {
96  std::string t1t2 = param_globals::imp_region_emi[i].ID[j];
97 
98  // Find the position of the ':' character
99  size_t colon_pos = t1t2.find(':');
100  if (colon_pos != std::string::npos) {
101  // Extract the two numbers using substr and convert to int
102  int tag1 = std::stoi(t1t2.substr(0, colon_pos));
103  int tag2 = std::stoi(t1t2.substr(colon_pos + 1));
104 
105  membraneFace_default.erase(std::make_pair(tag1, tag2));
106  gapjunctionFace_default.erase(std::make_pair(tag1, tag2));
107  // removed the other combinations of tags
108  membraneFace_default.erase(std::make_pair(tag2, tag1));
109  gapjunctionFace_default.erase(std::make_pair(tag2, tag1));
110  } else {
111  std::cerr << "Error: ':' the face tags are not defined properly in input file" << std::endl;
112  }
113  }
114  }
115 
116  if (rs.size() <2) {
117  log_msg(NULL, 5, ECHO, "\tError: num_imp_regions must be at least 2: the first region is reserved for the ionic model, and the second is designated for the gap junction.\n");
118  log_msg(NULL, 5, ECHO, "set a default ionic model and gapjuntion model in parameter file\\n");
119  exit(1);
120  }
121 
122  // default membrane faces
123  // imp_region_emi[0] selects the model for default membrane faces
124  rs[0].nsubregs = membraneFace_default.size();
125  rs[0].subregtags = new std::string[rs[0].nsubregs];
126  int j = 0;
127  for (const auto& v : membraneFace_default) {
128  int t1 = std::get<0>(v);
129  int t2 = std::get<1>(v);
130  rs[0].subregtags[j] = std::to_string(t1)+":"+std::to_string(t2); // Convert char* to std::string
131  j++;
132  }
133 
134  // default gap‑junction faces
135  // imp_region_emi[1] selects the model for default gap‑junction faces
136  j = 0;
137  rs[1].nsubregs = gapjunctionFace_default.size();
138  rs[1].subregtags = new std::string[rs[1].nsubregs];
139  for (const auto& v : gapjunctionFace_default) {
140  int t1 = std::get<0>(v);
141  int t2 = std::get<1>(v);
142  rs[1].subregtags[j] = std::to_string(t1)+":"+std::to_string(t2); // Convert char* to std::string
143  j++;
144  }
145 
146  // imp_region_emi[2+] assign models to specific face tag pairs like "t1:t2"
147  for (size_t i = 2; i < rs.size(); i++ ) {
148  rs[i].nsubregs = param_globals::imp_region_emi[i].num_IDs;
149  rs[i].subregtags = new std::string[rs[i].nsubregs]; // Allocate memory for the string array
150  for (int j = 0; j < rs[i].nsubregs;j++) {
151  std::string t1t2 = param_globals::imp_region_emi[i].ID[j];
152  size_t colon_pos = t1t2.find(':');
153  if (colon_pos == std::string::npos) {
154  std::cerr << "Error: ':' the face tags are not defined properly in input file, it should be tag1:tag2 as a string" << std::endl;
155  exit(1);
156  }
157 
158  rs[i].subregtags[j] = param_globals::imp_region_emi[i].ID[j]; // Convert char* to std::string
159  }
160  }
161 
162  SF::vector<int> reg_mask;
163  region_mask_onFace(ion_domain, tags_data, line_face, tri_face, quad_face, map_vertex_tag_to_dof, map_elem_uniqueFace_to_tags, rs, reg_mask, true, "imp_region_emi");
164 
165  for (size_t i = 0; i < rs.size(); i++ ) {
166  delete[] rs[i].subregtags;
167  rs[i].subregtags = nullptr;
168  }
169 
170  int purkfLen = 1;
171  tstart = setup_MIIF(loc_size, param_globals::num_imp_regions, param_globals::imp_region_emi,
172  reg_mask.data(), param_globals::start_statef, param_globals::num_adjustments,
173  param_globals::adjustment, param_globals::dt, purkfLen > 0);
174 
175  miif->extUpdateVm = !param_globals::operator_splitting;
176 
177  // if we start at a non-zero time (i.e. we have restored a state), we notify the
178  // timer manager
179  if(tstart > 0.) {
180  log_msg(logger, 0, 0, "Changing simulation start time to %.2lf", tstart);
181  user_globals::tm_manager->setup(param_globals::dt, tstart, param_globals::tend);
183  }
184 
185  set_dir(INPUT);
186 
187  this->initialize_time += timing(t2, t1);
188 }
189 
190 double IonicsOnFace::setup_MIIF(int nnodes, int nreg, IMPregion_EMI* impreg, int* mask,
191  const char *start_fn, int numadjust, IMPVariableAdjustment *adjust,
192  double time_step, bool close)
193 {
194  double tstart = 0;
195 
196  miif->N_IIF = nreg;
197  miif->numNode = nnodes;
198  miif->iontypes = {};
199  miif->numplugs = (int*)calloc( miif->N_IIF, sizeof(int));
200  miif->plugtypes = std::vector<limpet::IonTypeList>(miif->N_IIF);
201  miif->targets = std::vector<limpet::Target>(miif->N_IIF, limpet::Target::AUTO);
202 
203  log_msg(logger,0,ECHO, "\nSetting up ionic models on EMI Face and plugins\n" \
204  "-----------------------------------\n\n" \
205  "Assigning IMPS to tagged regions:" );
206 
207  for (int i=0;i<miif->N_IIF;i++) {
208  auto pT = limpet::get_ion_type(std::string(impreg[i].im));
209  if (pT != NULL)
210  {
211  miif->iontypes.push_back(*pT);
212  log_msg(logger, 0, ECHO|NONL, "\tIonic model: %s to tag region(s)", impreg[i].im);
213 
214  if(impreg[i].num_IDs > 0) {
215  for(int j = 0; j < impreg[i].num_IDs; j++)
216  log_msg(logger,0,ECHO|NONL, " [%s],", impreg[i].ID[j]);
217  log_msg(logger,0,ECHO,"\b.");
218  }
219  else {
220  log_msg(logger,0,ECHO, " [0] (implicitely)");
221  }
222  }
223  else {
224  log_msg(NULL,5,ECHO, "Illegal IM specified: %s\n", impreg[i].im );
225  log_msg(NULL,5,ECHO, "Run bench --list-imps for a list of all available models.\n" );
226  EXIT(1);
227  }
228  if (limpet::get_plug_flag( impreg[i].plugins, &miif->numplugs[i], miif->plugtypes[i]))
229  {
230  if(impreg[i].plugins[0] != '\0') {
231  log_msg(logger,0, ECHO|NONL, "\tPlug-in(s) : %s to tag region(s)", impreg[i].plugins);
232 
233  for(int j = 0; j < impreg[i].num_IDs; j++)
234  log_msg(logger,0,ECHO|NONL, " [%d],", impreg[i].ID[j]);
235  log_msg(logger,0,ECHO,"\b.");
236  }
237  }
238  else {
239  log_msg(NULL,5,ECHO,"Illegal plugin specified: %s\n", impreg[i].plugins);
240  log_msg(NULL,5,ECHO, "Run bench --list-imps for a list of all available plugins.\n" );
241  EXIT(1);
242  }
243  }
244 
245  miif->IIFmask = (limpet::IIF_Mask_t*)calloc(miif->numNode, sizeof(limpet::IIF_Mask_t));
246 
247  // The mask is a nodal vector with the region IDs of each node.
248  // It is already reduced during the region_mask call to guarantee unique values
249  // for overlapping interface nodes
250  if (mask) {
251  for (int i=0; i<miif->numNode; i++)
252  miif->IIFmask[i] = (limpet::IIF_Mask_t) mask[i];
253  }
254 
255  miif->initialize_MIIF();
256 
257  for (int i=0;i<miif->N_IIF;i++) {
258  // the IMP tuning does not handle spaces well, thus we remove them here
259  remove_char(impreg[i].im_param, strlen(impreg[i].im_param), ' ');
260  miif->IIF[i]->tune(impreg[i].im_param, impreg[i].plugins, impreg[i].plug_param);
261  }
262 
263  set_dir(INPUT);
264  miif->initialize_currents(time_step, param_globals::ode_fac);
265 
266  // overriding initial values goes here
267  // read in single cell state vector and spread it out over the entire region
268  set_dir(INPUT);
269  for (int i=0;i<miif->N_IIF;i++) {
270  if (impreg[i].im_sv_init && strlen(impreg[i].im_sv_init) > 0)
271  if (read_sv(miif, i, impreg[i].im_sv_init)) {
272  log_msg(NULL, 5, ECHO|FLUSH, "State vector initialization failed for %s.\n", impreg[i].name);
273  EXIT(-1);
274  }
275  }
276 
277  if( !start_fn || strlen(start_fn)>0 )
278  tstart = (double) miif->restore_state(start_fn, ion_domain, close);
279 
280  for (int i=0; i<numadjust; i++)
281  {
282  set_dir(INPUT);
283 
284  SF::vector<int> indices;
285  SF::vector<double> values;
286  bool restrict_to_algebraic = true;
287 
288  sf_mesh & imesh = get_mesh(ion_domain);
289  std::map<std::string,std::string> metadata;
290  read_metadata(adjust[i].file, metadata, PETSC_COMM_WORLD);
291 
292  SF::SF_nbr nbr = SF::NBR_REF;
293  if(metadata.count("grid") && metadata["grid"].compare("intra") == 0) {
294  nbr = SF::NBR_SUBMESH;
295  }
296  read_indices_with_data(indices, values, adjust[i].file, imesh, nbr, restrict_to_algebraic, 1, PETSC_COMM_WORLD);
297 
298  int rank = get_rank();
299 
300  for(size_t gi = 0; gi < indices.size(); gi++)
301  indices[gi] = SF::local_nodal_to_local_petsc(imesh, rank, indices[gi]);
302 
303  // debug, output parameters on global intracellular vector
304  if(adjust[i].dump) {
305  sf_vec* adjPars;
306  SF::init_vector(&adjPars, imesh, 1, sf_vec::algebraic);
307  adjPars->set(indices, values);
308 
309  set_dir(OUTPUT);
310  char fname[2085];
311  snprintf(fname, sizeof fname, "adj_%s_perm.dat", adjust[i].variable);
312  adjPars->write_ascii(fname, false);
313 
314  // get the scattering to the canonical permutation
315  SF::scattering* sc = get_permutation(ion_domain, PETSC_TO_CANONICAL, 1);
316  if(sc == NULL) {
317  log_msg(0,3,0, "%s warning: PETSC_TO_CANONICAL permutation needed registering!", __func__);
318  sc = register_permutation(ion_domain, PETSC_TO_CANONICAL, 1);
319  }
320 
321  (*sc)(*adjPars, true);
322  snprintf(fname, sizeof fname, "adj_%s_canonical.dat", adjust[i].variable);
323  adjPars->write_ascii(fname, false);
324  }
325 
326  int nc = miif->adjust_MIIF_variables(adjust[i].variable, indices, values);
327  log_msg(logger, 0, 0, "Adjusted %d values for %s", nc, adjust[i].variable);
328  }
329 
330  set_dir(OUTPUT);
331 
332  for (int i=0;i<miif->N_IIF;i++)
333  initialize_sv_dumps_onFace(miif, impreg+i, i, tstart, param_globals::spacedt);
334 
335  return tstart;
336 }
337 
347 void initialize_sv_dumps_onFace(limpet::MULTI_IF *pmiif, IMPregion_EMI* reg, int id, double t, double dump_dt)
348 {
349  char svs[1024], plgs[1024], plgsvs[1024], fname[1024];
350 
351  strcpy(svs, reg->im_sv_dumps ? reg->im_sv_dumps : "");
352  strcpy(plgs, reg->plugins ? reg->plugins : "");
353  strcpy(plgsvs, reg->plug_sv_dumps ? reg->plug_sv_dumps : "");
354 
355  if( !(strlen(svs)+strlen(plgsvs) ) )
356  return;
357 
358  /* The string passed to the "reg_name" argument (#4) of the sv_dump_add
359  * function is supposed to be "region name". It's only purpose is to
360  * provide the base name for the SV dump file, eg: Purkinje.Ca_i.bin.
361  * Thus, we pass: [vofile].[reg name], Otherwise, dumping SVs in
362  * batched runs would be extremely tedious.
363  */
364  strcpy(fname, param_globals::vofile); // [vofile].igb.gz
365 
366  if( !reg->name ) {
367  log_msg(NULL, 5, ECHO, "%s: a region name must be specified\n", __func__ );
368  exit(0);
369  }
370 
371  // We want to convert vofile.igb or vofile.igb.gz to vofile.regname
372  char* ext_start = NULL;
373  ext_start = strstr(fname, "igb.gz");
374  if(ext_start == NULL) ext_start = strstr(fname, "igb");
375  if(ext_start == NULL) ext_start = fname + strlen(fname);
376  strcpy(ext_start, reg->name);
377 
378  pmiif->sv_dump_add_by_name_list(id, reg->im, fname, svs, plgs, plgsvs, t, dump_dt);
379 }
380 
391 bool check_tags_in_elems_onFace(std::vector<std::string> & tags_data, SF::vector<RegionSpecs_EMI> & regspec,
392  const char* gridname, const char* reglist)
393 {
394  bool AllTagsExist = true;
396 
397  tagset.insert(tags_data.begin(), tags_data.end());
398 
399  // cycle through all user-specified regions
400  for (size_t reg=0; reg<regspec.size(); reg++)
401  // cycle through all tags which belong to the region
402  for (int k=0; k<regspec[reg].nsubregs; k++) {
403  // check whether this tag exists in element list
404  int n = 0;
405  size_t colon_pos = regspec[reg].subregtags[k].find(':');
406  int tag1 = std::stoi(regspec[reg].subregtags[k].substr(0, colon_pos));
407  int tag2 = std::stoi(regspec[reg].subregtags[k].substr(colon_pos + 1));
408 
409  // considered the other combinations of tags
410  std::string inverse_pair;
411  inverse_pair = std::to_string(tag2) + ":" + std::to_string(tag1);
412  if(tagset.count(regspec[reg].subregtags[k]) || tagset.count(inverse_pair)) {
413  n++;
414  }
415 
416  // globalize n
417  int N = get_global(n, MPI_SUM);
418 
419  if (N==0) {
420  log_msg(NULL, 3, ECHO,
421  "on face Region tag %s in %s[%d] not found in element list for %s grid.\n",
422  regspec[reg].subregtags[k].c_str(), reglist, reg, gridname);
423  AllTagsExist = false;
424  }
425  }
426 
427  if (!AllTagsExist) {
428  log_msg(NULL, 4, ECHO,"Assigned wrong pair of tags on the face of EMI surface mesh!\n"
429  "Check region ID specs in input files!\n");
430  }
431 
432  return AllTagsExist;
433 }
434 
435 void region_mask_onFace(mesh_t meshspec,
436  std::vector<std::string> & tags_data,
438  std::pair<SF::emi_face<mesh_int_t,SF::tuple<mesh_int_t>>,
439  SF::emi_face<mesh_int_t,SF::tuple<mesh_int_t>>>> & line_face,
441  std::pair<SF::emi_face<mesh_int_t,SF::triple<mesh_int_t>>,
442  SF::emi_face<mesh_int_t,SF::triple<mesh_int_t>>>> & tri_face,
444  std::pair<SF::emi_face<mesh_int_t,SF::quadruple<mesh_int_t>>,
445  SF::emi_face<mesh_int_t,SF::quadruple<mesh_int_t>>>> & quad_face,
446  hashmap::unordered_map<std::pair<mesh_int_t,mesh_int_t>, mesh_int_t> & map_vertex_tag_to_dof,
447  hashmap::unordered_map<mesh_int_t, std::pair<mesh_int_t, mesh_int_t>> & map_elem_uniqueFace_to_tags,
448  SF::vector<RegionSpecs_EMI> & regspec,
449  SF::vector<int> & regionIDs, bool mask_elem, const char* reglist)
450 {
451  if(regspec.size() == 1) return;
452 
453  sf_mesh & mesh = get_mesh(meshspec);
454  SF::element_view<mesh_int_t,mesh_real_t> eview(mesh, SF::NBR_PETSC); // error is here!!!
455 
456  // initialize the list with the default regionID, 0
457  size_t rIDsize = mask_elem ? mesh.l_numelem : mesh.l_numpts;
458 
459  size_t nelem = mesh.l_numelem;
460 
461  // check whether all specified tags exist in the element list
462  check_tags_in_elems_onFace(tags_data, regspec, mesh.name.c_str(), reglist);
463 
464  regionIDs.assign(rIDsize, 0);
465 
466  int* rid = regionIDs.data();
468 
469  // we generate a map from tags to region IDs. This has many benefits, mainly we can check
470  // whether a tag is assigned to multiple regions and simplify our regionIDs filling loop
471  for (size_t reg=0; reg < regspec.size(); reg++) {
472  int err = 0;
473  for (int k=0; k < regspec[reg].nsubregs; k++) {
474  std::string curtag = regspec[reg].subregtags[k];
475  if(tag_to_reg.count(curtag)) err++;
476 
477  if(get_global(err, MPI_SUM))
478  log_msg(0,4,0, "%s warning: Tag idx %s is assigned to multiple regions!\n"
479  "Its final assignment will be to the highest assigned region ID!",
480  __func__, curtag.c_str());
481 
482  tag_to_reg[curtag] = reg;
483  }
484  }
485 
486  std::vector<int> mx_tag( rIDsize, -1 );
487 
488  // cycle through the element list
489  for(size_t eidx=0; eidx<nelem; eidx++)
490  {
491  std::vector<int> elem_nodes;
492  mesh_int_t tag = mesh.tag[eidx];
493  std::string result_pair_orginal;
494  std::string result_pair_reverse;
495  std::pair<mesh_int_t,mesh_int_t> value = map_elem_uniqueFace_to_tags[eidx];
496  result_pair_orginal = std::to_string(value.first) + ":" + std::to_string(value.second);
497  result_pair_reverse = std::to_string(value.second) + ":" + std::to_string(value.first);
498  if (tag_to_reg.count(result_pair_orginal))
499  {
500  rid[eidx] = tag_to_reg[result_pair_orginal];
501  }else if(tag_to_reg.count(result_pair_reverse)){
502  rid[eidx] = tag_to_reg[result_pair_reverse];
503  }
504  if(tag!=value.first){
505  log_msg(NULL, 5, ECHO, "error in tags on surface mesh with the unique faces!!!.\\n");
506  exit(1);
507  }
508  }
509 }
510 
513 double IonicsOnFace::timer_val(const int timer_id)
514 {
515  double val = std::nan("NaN");
516  return val;
517 }
518 
521 std::string IonicsOnFace::timer_unit(const int timer_id)
522 {
523  std::string s_unit;
524  return s_unit;
525 }
526 
527 void compute_IIF_OnFace(limpet::IonIfBase& pIF, limpet::GlobalData_t** impdata, int n)
528 {
529  if (impdata[limpet::Iion] != NULL)
530  impdata[limpet::Iion][n] = 0;
531 
532  pIF.for_each([&](limpet::IonIfBase& imp) {
533  update_ts(&imp.get_tstp());
534  imp.compute(n, n + 1, impdata);
535  });
536 }
537 
549 void* find_SV_in_IMP_onFace(limpet::MULTI_IF* miif, const int idx, const char *IMP, const char *SV,
550  int* offset, int* sz)
551 {
552  if(strcmp(IMP, miif->iontypes[idx].get().get_name().c_str()) == 0) {
553  return (void*) miif->iontypes[idx].get().get_sv_offset(SV, offset, sz);
554  }
555  else {
556  for( int k=0; k<miif->numplugs[idx]; k++ )
557  if(strcmp(IMP, miif->plugtypes[idx][k].get().get_name().c_str()) == 0) {
558  return (void*) miif->plugtypes[idx][k].get().get_sv_offset(SV, offset, sz);
559  }
560  }
561 
562  return NULL;
563 }
564 
565 
575 void alloc_gvec_data_onFace(const int nGVcs, const int nRegs,
576  GVecs *prmGVecs, gvec_data_OnFace &glob_vecs)
577 {
578  glob_vecs.nRegs = nRegs;
579 
580  if (nGVcs) {
581  glob_vecs.vecs.resize(nGVcs);
582 
583  for (size_t i = 0; i < glob_vecs.vecs.size(); i++) {
584  sv_data_onFace &gvec = glob_vecs.vecs[i];
585 
586  gvec.name = dupstr(prmGVecs[i].name);
587  gvec.units = dupstr(prmGVecs[i].units);
588  gvec.bogus = prmGVecs[i].bogus;
589 
590  gvec.imps = (char**) calloc(nRegs, sizeof(char *));
591  gvec.svNames = (char**) calloc(nRegs, sizeof(char *));
592  gvec.svSizes = (int*) calloc(nRegs, sizeof(int));
593  gvec.svOff = (int*) calloc(nRegs, sizeof(int));
594  gvec.getsv = (void**) calloc(nRegs, sizeof(limpet::SVgetfcn));
595 
596  for (int j = 0; j < nRegs; j++) {
597  if (strlen(prmGVecs[i].imp)) gvec.imps[j] = dupstr(prmGVecs[i].imp);
598 #ifdef WITH_PURK
599  else if (j < param_globals::num_imp_regions)
600  gvec.imps[j] = dupstr(param_globals::imp_region_emi[j].im);
601  else
602  gvec.imps[j] = dupstr(param_globals::PurkIon[j - param_globals::num_imp_regions].im);
603 #else
604  else if (j < param_globals::num_imp_regions)
605  gvec.imps[j] = dupstr(param_globals::imp_region_emi[j].im);
606 #endif
607  gvec.svNames[j] = dupstr(prmGVecs[i].ID[j]);
608  }
609  }
610  }
611 }
612 
625 void init_sv_gvec_onFace(gvec_data_OnFace& GVs, limpet::MULTI_IF* miif, sf_vec & tmpl,
626  igb_output_manager & output_manager)
627 {
628  GVs.inclPS = false;
629  // int num_purk_regions = GVs->inclPS ? purk->ion.N_IIF : 0;
630  int num_purk_regions = 0;
631  int nRegs = param_globals::num_imp_regions + num_purk_regions;
632 
633  alloc_gvec_data_onFace(param_globals::num_gvecs, nRegs, param_globals::gvec, GVs);
634 
635 #ifdef WITH_PURK
636  if (GVs->inclPS) sample_PS_ionSVs(purk);
637 #endif
638 
639  for (unsigned int i = 0; i < GVs.vecs.size(); i++) {
640  sv_data_onFace & gv = GVs.vecs[i];
641  int noSV = 0;
642 
643  for (int j = 0; j < miif->N_IIF; j++) {
644  gv.getsv[j] = find_SV_in_IMP_onFace(miif, j, gv.imps[j], gv.svNames[j],
645  gv.svOff + j, gv.svSizes + j);
646 
647  if (gv.getsv[j] == NULL) {
648  log_msg(NULL, 3, ECHO, "\tWarning: SV(%s) not found in region %d\n", gv.svNames[j], j);
649  noSV++;
650  }
651  }
652 
653  SF::init_vector(&gv.ordered, &tmpl);
654  output_manager.register_output(gv.ordered, intra_elec_msh, 1, gv.name, gv.units);
655 
656 #ifdef WITH_PURK
657  // same procedure for Purkinje
658  if (GVs->inclPS) {
659  IF_PURK_PROC(purk) {
660  MULTI_IF* pmiif = &purk->ion;
661  for (int j = miif->N_IIF; j < nRegs; j++) {
662  gv.getsv[j] = find_SV_in_IMP_onFace(pmiif, j - miif->N_IIF, gv.imps[j],
663  gv.svNames[j], gv.svOff + j, gv.svSizes + j);
664 
665  if (gv.getsv[j] == NULL) {
666  LOG_MSG(NULL, 3, ECHO, "\tWarning: state variable \"%s\" not found in region %d\n", gv.svNames[j], j);
667  noSV++;
668  }
669  }
670  RVector_dup(purk->vm_pt, &gv.orderedPS_PS);
671  MYO_COMM(purk);
672  }
673  RVector_dup(purk->vm_pt_over, &gv.orderedPS);
674  initialize_grid_output(grid, NULL, tmo, intra_elec_msh, GRID_WRITE, 0., 1., gv.units, gv.orderedPS,
675  1, gv.GVcName, -purk->npt, param_globals::output_level);
676  }
677 #endif
678 
679  if (noSV == nRegs) {
680  log_msg(NULL, 5, ECHO, "\tError: no state variables found for global vector %d\n", i);
681  log_msg(NULL, 5, ECHO, "Run bench --imp=YourModel --imp-info to get a list of all parameters.\\n");
682  exit(1);
683  }
684  }
685 }
686 
697 void assemble_sv_gvec_onFace(gvec_data_OnFace & gvecs, limpet::MULTI_IF *miif)
698 {
699  for(size_t i=0; i<gvecs.vecs.size(); i++ ) {
700  sv_data_onFace & gv = gvecs.vecs[i];
701 
702  // set to the defalt value
703  gv.ordered->set(gv.bogus);
704  SF_int start, stop;
705  gv.ordered->get_ownership_range(start, stop);
706 
707  for( int n = 0; n<miif->N_IIF; n++ ) {
708  if( !gv.getsv[n] ) continue;
709 
710  SF::vector<SF_real> data (miif->N_Nodes[n]);
711  SF::vector<SF_int> indices(miif->N_Nodes[n]);
712 
713  for( int j=0; j<miif->N_Nodes[n]; j++ ) {
714  indices[j] = miif->NodeLists[n][j] + start;
715  data[j] = ((limpet::SVgetfcn)(gv.getsv[n]))( *miif->IIF[n], j, gv.svOff[n]);
716  }
717 
718  bool add = false;
719  gv.ordered->set(indices, data, add);
720  }
721 
722 #ifdef WITH_PURK
723  // include purkinje in sv dump, if exists
724  if(gvecs->inclPS) {
725  RVector_set( gv.orderedPS, gv.bogus ); // set this to the default value
726 
727  if(IS_PURK_PROC(purk))
728  {
729  MULTI_IF *pmiif = &purk->ion;
730  Real *data = new Real[purk->gvec_cab.nitems];
731 
732  int ci=0;
733  for( int n = 0; n<pmiif->N_IIF; n++ ) {
734  int gvidx = n+miif->N_IIF;
735  if( !gv.getsv[gvidx] ) continue;
736  for( int j=0; j<purk->gvec_ion[n].nitems; j++ )
737  data[ci++] = ((SVgetfcn)(gv.getsv[gvidx]))( pmiif->IIF+n, ((int*)(purk->gvec_ion[n].data))[j],
738  gv.svOff[gvidx]);
739  }
740  RVector_setvals(gv.orderedPS, purk->gvec_cab.nitems, (int*)purk->gvec_cab.data, data, true);
741  delete [] data;
742  }
743  RVector_sync( gv.orderedPS );
744  }
745 #endif
746  }
747 }
748 
749 
750 } // namespace opencarp
751 #endif
int mesh_int_t
Definition: SF_container.h:46
std::int32_t SF_int
Use the general std::int32_t as int type.
Definition: SF_globals.h:37
#define FLUSH
Definition: basics.h:311
#define ECHO
Definition: basics.h:308
#define NONL
Definition: basics.h:312
Comfort class. Provides getter functions to access the mesh member variables more comfortably.
Definition: SF_fem_utils.h:704
Container for a PETSc VecScatter.
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 assign(InputIterator s, InputIterator e)
Assign a memory range.
Definition: SF_vector.h:161
T * data()
Pointer to the vector's start.
Definition: SF_vector.h:91
hm_int count(const K &key) const
Check if key exists.
Definition: hashmap.hpp:579
Custom unordered_set implementation.
Definition: hashmap.hpp:706
hm_int count(const K &key) const
Definition: hashmap.hpp:992
void insert(InputIterator first, InputIterator last)
Definition: hashmap.hpp:962
Represents the ionic model and plug-in (IMP) data structure.
Definition: ION_IF.h:139
ts & get_tstp()
Gets the time stepper.
Definition: ION_IF.cc:208
void for_each(const std::function< void(IonIfBase &)> &consumer)
Executes the consumer functions on this IMP and each of its plugins.
Definition: ION_IF.cc:417
void compute(int start, int end, GlobalData_t **data)
Perform ionic model computation for 1 time step.
Definition: ION_IF.cc:270
std::vector< IonIfBase * > IIF
array of IIF's
Definition: MULTI_ION_IF.h:202
void sv_dump_add_by_name_list(int, char *, char *, char *, char *, char *, double, double)
int * numplugs
number of plugins for each region
Definition: MULTI_ION_IF.h:209
std::vector< IonTypeList > plugtypes
plugins types for each region
Definition: MULTI_ION_IF.h:215
IonTypeList iontypes
type for each region
Definition: MULTI_ION_IF.h:212
int N_IIF
how many different IIF's
Definition: MULTI_ION_IF.h:211
int * N_Nodes
#nodes for each IMP
Definition: MULTI_ION_IF.h:200
int ** NodeLists
local partitioned node lists for each IMP stored
Definition: MULTI_ION_IF.h:201
void setup(double inp_dt, double inp_start, double inp_end)
Initialize the timer_manager.
Definition: timer_utils.cc:36
void reset_timers()
Reset time in timer_manager and then reset registered timers.
Definition: timer_utils.h:115
Electrical ionics functions (gap junction + ionic current) on the face of interface mesh for EMI mesh...
T local_nodal_to_local_petsc(const meshdata< T, S > &mesh, int rank, T local_nodal)
void init_vector(SF::abstract_vector< T, S > **vec)
Definition: SF_init.h:99
SF_nbr
Enumeration encoding the different supported numberings.
Definition: SF_container.h:200
@ NBR_PETSC
PETSc numbering of nodes.
Definition: SF_container.h:203
@ 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
int get_plug_flag(char *plgstr, int *out_num_plugins, IonTypeList &out_plugins)
GlobalData_t(* SVgetfcn)(IonIfBase &, int, int)
Definition: ion_type.h:48
double Real
Definition: MULTI_ION_IF.h:151
@ AUTO
Definition: target.h:46
IonType * get_ion_type(const std::string &name)
SF_real GlobalData_t
Definition: limpet_types.h:27
void update_ts(ts *ptstp)
Definition: ION_IF.cc:466
int read_sv(MULTI_IF *, int, const char *)
char IIF_Mask_t
Definition: ion_type.h:50
std::map< int, std::string > units
Definition: stimulate.cc:41
timer_manager * tm_manager
a manager for the various physics timers
Definition: main.cc:58
@ iotm_console
Definition: timer_utils.h:44
sf_mesh & get_mesh(const mesh_t gt)
Get a mesh by specifying the gridID.
Definition: sf_interface.cc:33
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
int set_dir(IO_t dest)
Definition: sim_utils.cc:490
void read_metadata(const std::string filename, std::map< std::string, std::string > &metadata, MPI_Comm comm)
Read metadata from the header.
Definition: fem_utils.cc:72
int get_rank(MPI_Comm comm=PETSC_COMM_WORLD)
Definition: basics.h:276
T get_global(T in, MPI_Op OP, MPI_Comm comm=PETSC_COMM_WORLD)
Do a global reduction on a variable.
Definition: basics.h:233
SF::scattering * register_permutation(const int mesh_id, const int perm_id, const int dpn)
Register a permutation between two orderings for a mesh.
void register_data(sf_vec *dat, datavec_t d)
Register a data vector in the global registry.
Definition: sim_utils.cc:897
@ OUTPUT
Definition: sim_utils.h:53
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
@ intra_elec_msh
Definition: sf_interface.h:60
void get_time(double &tm)
Definition: basics.h:436
SF::abstract_vector< SF_int, SF_real > sf_vec
Definition: sf_interface.h:50
void remove_char(char *buff, const int buffsize, const char c)
Definition: basics.h:356
void read_indices_with_data(SF::vector< T > &idx, SF::vector< S > &dat, const std::string filename, const hashmap::unordered_map< mesh_int_t, mesh_int_t > &dd_map, const int dpn, MPI_Comm comm)
like read_indices, but with associated data for each index
Definition: fem_utils.h:269
V timing(V &t2, const V &t1)
Definition: basics.h:448
#define PETSC_TO_CANONICAL
Permute algebraic data from PETSC to canonical ordering.
Definition: sf_interface.h:79