openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
ionics.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 #include "ionics.h"
28 
29 #include "SF_init.h"
30 
31 namespace opencarp {
32 
33 void initialize_sv_dumps(limpet::MULTI_IF *pmiif, IMPregion* reg, int id, double t, double dump_dt);
34 
36 {
37  double t1, t2;
38  get_time(t1);
39 
41 
42  double comp_time = timing(t2, t1);
43  this->compute_time += comp_time;
44 
45  comp_stats.calls++;
46  comp_stats.tot_time += comp_time;
47 
50 }
51 
53 {
54  miif->free_MIIF();
55 }
56 
58 {}
59 
61 {
62  double t1, t2;
63  get_time(t1);
64 
65  set_dir(OUTPUT);
66 
67  // initialize generic logger for ODE timings per time_dt
68  comp_stats.init_logger("ODE_stats.dat");
69 
70  double tstart = 0.;
71  sf_mesh & mesh = get_mesh(ion_domain);
72  int loc_size = mesh.pl.num_algebraic_idx();
73 
74  // create ionic current vector and register it
75  sf_vec *IIon;
76  sf_vec *Vmv;
77  SF::init_vector(&IIon, mesh, 1, sf_vec::algebraic);
78  SF::init_vector(&Vmv, mesh, 1, sf_vec::algebraic);
79  register_data(IIon, iion_vec);
80  register_data(Vmv, vm_vec);
81 
82  // setup miif
83  miif = new limpet::MULTI_IF();
84  // store IIF_IDs and Plugins in arrays
85  miif->name = "myocardium";
86  miif->gdata[limpet::Vm] = Vmv;
87  miif->gdata[limpet::Iion] = IIon;
88 
89  int rank = get_rank();
90 
91  SF::vector<RegionSpecs> rs(param_globals::num_imp_regions);
92  for (size_t i = 0; i < rs.size(); i++ ) {
93  rs[i].nsubregs = param_globals::imp_region[i].num_IDs;
94  rs[i].subregtags = param_globals::imp_region[i].ID;
95  for (int j=0;j<rs[i].nsubregs;j++) {
96  if(rs[i].subregtags[j]==-1 && rank==0)
97  log_msg(NULL,3,ECHO, "Warning: not all %u IDs provided for imp_region[%u]!\n", rs[i].nsubregs, i);
98  }
99  }
100 
101  SF::vector<int> reg_mask;
102  region_mask(ion_domain, rs, reg_mask, false, "imp_region");
103 
104  int purkfLen = 1;
105  tstart = setup_MIIF(loc_size, param_globals::num_imp_regions, param_globals::imp_region,
106  reg_mask.data(), param_globals::start_statef, param_globals::num_adjustments,
107  param_globals::adjustment, param_globals::dt, purkfLen > 0);
108 
109  miif->extUpdateVm = !param_globals::operator_splitting;
110 
111  // if we start at a non-zero time (i.e. we have restored a state), we notify the
112  // timer manager
113  if(tstart > 0.) {
114  log_msg(logger, 0, 0, "Changing simulation start time to %.2lf", tstart);
115  user_globals::tm_manager->setup(param_globals::dt, tstart, param_globals::tend);
117  }
118 
119  set_dir(INPUT);
120 
121  this->initialize_time += timing(t2, t1);
122 }
123 
124 double Ionics::setup_MIIF(int nnodes, int nreg, IMPregion* impreg, int* mask,
125  const char *start_fn, int numadjust, IMPVariableAdjustment *adjust,
126  double time_step, bool close)
127 {
128  double tstart = 0;
129 
130  miif->N_IIF = nreg;
131  miif->numNode = nnodes;
132  miif->iontypes = {};
133  miif->numplugs = (int*)calloc( miif->N_IIF, sizeof(int));
134  miif->plugtypes = std::vector<limpet::IonTypeList>(miif->N_IIF);
135  miif->targets = std::vector<limpet::Target>(miif->N_IIF, limpet::Target::AUTO);
136 
137  log_msg(logger,0,ECHO, "\nSetting up ionic models and plugins\n" \
138  "-----------------------------------\n\n" \
139  "Assigning IMPS to tagged regions:" );
140 
141  for (int i=0;i<miif->N_IIF;i++) {
142  auto pT = limpet::get_ion_type(std::string(impreg[i].im));
143  if (pT != NULL)
144  {
145  miif->iontypes.push_back(*pT);
146  log_msg(logger, 0, ECHO|NONL, "\tIonic model: %s to tag region(s)", impreg[i].im);
147 
148  if(impreg[i].num_IDs > 0) {
149  for(int j = 0; j < impreg[i].num_IDs; j++)
150  log_msg(logger,0,ECHO|NONL, " [%d],", impreg[i].ID[j]);
151  log_msg(logger,0,ECHO,"\b.");
152  }
153  else {
154  log_msg(logger,0,ECHO, " [0] (implicitely)");
155  }
156  }
157  else {
158  log_msg(NULL,5,ECHO, "Illegal IM specified: %s\n", impreg[i].im );
159  log_msg(NULL,5,ECHO, "Run bench --list-imps for a list of all available models.\n" );
160  EXIT(1);
161  }
162  if (limpet::get_plug_flag( impreg[i].plugins, &miif->numplugs[i], miif->plugtypes[i]))
163  {
164  if(impreg[i].plugins[0] != '\0') {
165  log_msg(logger,0, ECHO|NONL, "\tPlug-in(s) : %s to tag region(s)", impreg[i].plugins);
166 
167  for(int j = 0; j < impreg[i].num_IDs; j++)
168  log_msg(logger,0,ECHO|NONL, " [%d],", impreg[i].ID[j]);
169  log_msg(logger,0,ECHO,"\b.");
170  }
171  }
172  else {
173  log_msg(NULL,5,ECHO,"Illegal plugin specified: %s\n", impreg[i].plugins);
174  log_msg(NULL,5,ECHO, "Run bench --list-imps for a list of all available plugins.\n" );
175  EXIT(1);
176  }
177  }
178 
180 
181  // The mask is a nodal vector with the region IDs of each node.
182  // It is already reduced during the region_mask call to guarantee unique values
183  // for overlapping interface nodes
184  if (mask) {
185  for (int i=0; i<miif->numNode; i++)
186  miif->IIFmask[i] = (limpet::IIF_Mask_t) mask[i];
187  }
188 
190 
191  for(int i=0; i<miif->N_IIF; i++)
192  if(!miif->IIF[i]->cgeom().SVratio)
193  miif->IIF[i]->cgeom().SVratio = param_globals::imp_region[i].cellSurfVolRatio;
194 
195  for (int i=0;i<miif->N_IIF;i++) {
196  // the IMP tuning does not handle spaces well, thus we remove them here
197  remove_char(impreg[i].im_param, strlen(impreg[i].im_param), ' ');
198  miif->IIF[i]->tune(impreg[i].im_param, impreg[i].plugins, impreg[i].plug_param);
199  }
200 
201  set_dir(INPUT);
202  miif->initialize_currents(time_step, param_globals::ode_fac);
203 
204  // overriding initial values goes here
205  // read in single cell state vector and spread it out over the entire region
206  set_dir(INPUT);
207  for (int i=0;i<miif->N_IIF;i++) {
208  if (impreg[i].im_sv_init && strlen(impreg[i].im_sv_init) > 0)
209  if (read_sv(miif, i, impreg[i].im_sv_init)) {
210  log_msg(NULL, 5, ECHO|FLUSH, "State vector initialization failed for %s.\n", impreg[i].name);
211  EXIT(-1);
212  }
213  }
214 
215  if( !start_fn || strlen(start_fn)>0 )
216  tstart = (double) miif->restore_state(start_fn, ion_domain, close);
217 
218  for (int i=0; i<numadjust; i++)
219  {
220  set_dir(INPUT);
221 
222  SF::vector<SF_int> indices;
223  SF::vector<SF_real> values;
224  bool restrict_to_algebraic = true;
225 
226  sf_mesh & imesh = get_mesh(ion_domain);
227  std::map<std::string,std::string> metadata;
228  read_metadata(adjust[i].file, metadata, PETSC_COMM_WORLD);
229 
230  SF::SF_nbr nbr = SF::NBR_REF;
231  if(metadata.count("grid") && metadata["grid"].compare("intra") == 0) {
232  nbr = SF::NBR_SUBMESH;
233  }
234  read_indices_with_data(indices, values, adjust[i].file, imesh, nbr, restrict_to_algebraic, 1, PETSC_COMM_WORLD);
235 
236  int rank = get_rank();
237 
238  for(size_t gi = 0; gi < indices.size(); gi++)
239  indices[gi] = SF::local_nodal_to_local_petsc<mesh_int_t,mesh_real_t>(imesh, rank, indices[gi]);
240 
241  // debug, output parameters on global intracellular vector
242  if(adjust[i].dump) {
243  sf_vec* adjPars;
244  SF::init_vector(&adjPars, imesh, 1, sf_vec::algebraic);
245  adjPars->set(indices, values);
246 
247  set_dir(OUTPUT);
248  char fname[2085];
249  snprintf(fname, sizeof fname, "adj_%s_perm.dat", adjust[i].variable);
250  adjPars->write_ascii(fname, false);
251 
252  // get the scattering to the canonical permutation
254  if(sc == NULL) {
255  log_msg(0,3,0, "%s warning: PETSC_TO_CANONICAL permutation needed registering!", __func__);
257  }
258 
259  (*sc)(*adjPars, true);
260  snprintf(fname, sizeof fname, "adj_%s_canonical.dat", adjust[i].variable);
261  adjPars->write_ascii(fname, false);
262  }
263 
264  int nc = miif->adjust_MIIF_variables(adjust[i].variable, indices, values);
265  log_msg(logger, 0, 0, "Adjusted %d values for %s", nc, adjust[i].variable);
266  }
267 
268  set_dir(OUTPUT);
269 
270  for (int i=0;i<miif->N_IIF;i++)
271  initialize_sv_dumps(miif, impreg+i, i, tstart, param_globals::spacedt);
272 
273  return tstart;
274 }
275 
285 void initialize_sv_dumps(limpet::MULTI_IF *pmiif, IMPregion* reg, int id, double t, double dump_dt)
286 {
287  char svs[1024], plgs[1024], plgsvs[1024], fname[1024];
288 
289  strcpy(svs, reg->im_sv_dumps ? reg->im_sv_dumps : "");
290  strcpy(plgs, reg->plugins ? reg->plugins : "");
291  strcpy(plgsvs, reg->plug_sv_dumps ? reg->plug_sv_dumps : "");
292 
293  if( !(strlen(svs)+strlen(plgsvs) ) )
294  return;
295 
296  /* The string passed to the "reg_name" argument (#4) of the sv_dump_add
297  * function is supposed to be "region name". It's only purpose is to
298  * provide the base name for the SV dump file, eg: Purkinje.Ca_i.bin.
299  * Thus, we pass: [vofile].[reg name], Otherwise, dumping SVs in
300  * batched runs would be extremely tedious.
301  */
302  strcpy(fname, param_globals::vofile); // [vofile].igb.gz
303 
304  if( !reg->name ) {
305  log_msg(NULL, 5, ECHO, "%s: a region name must be specified\n", __func__ );
306  exit(0);
307  }
308 
309  // We want to convert vofile.igb or vofile.igb.gz to vofile.regname
310  char* ext_start = NULL;
311  ext_start = strstr(fname, "igb.gz");
312  if(ext_start == NULL) ext_start = strstr(fname, "igb");
313  if(ext_start == NULL) ext_start = fname + strlen(fname);
314  strcpy(ext_start, reg->name);
315 
316  pmiif->sv_dump_add_by_name_list(id, reg->im, fname, svs, plgs, plgsvs, t, dump_dt);
317 }
318 
330  const char* gridname, const char* reglist)
331 {
332  bool AllTagsExist = true;
333 
335  tagset.insert(tags.begin(), tags.end());
336 
337  // cycle through all user-specified regions
338  for (size_t reg=0; reg<regspec.size(); reg++)
339  // cycle through all tags which belong to the region
340  for (int k=0; k<regspec[reg].nsubregs; k++) {
341  // check whether this tag exists in element list
342  int n = 0;
343  if(tagset.count(regspec[reg].subregtags[k])) n++;
344  // globalize n
345  int N = get_global(n, MPI_SUM);
346  if (N==0) {
347  if(strcmp(reglist, "gregion_vol"))
348  log_msg(NULL, 3, ECHO,
349  "Region tag %d in %s[%d] not found in element list for %s grid.\n",
350  regspec[reg].subregtags[k], reglist, reg, gridname);
351  AllTagsExist = false;
352  }
353  }
354 
355  if (!AllTagsExist) {
356  log_msg(NULL, 4, ECHO,"Missing element tags in element file!\n"
357  "Check region ID specs in input files!\n");
358  }
359 
360  return AllTagsExist;
361 }
362 
363 void region_mask(mesh_t meshspec, SF::vector<RegionSpecs> & regspec,
364  SF::vector<int> & regionIDs, bool mask_elem, const char* reglist)
365 {
366  if(regspec.size() == 1) return;
367 
368  sf_mesh & mesh = get_mesh(meshspec);
370 
371  // initialize the list with the default regionID, 0
372  size_t rIDsize = mask_elem ? mesh.l_numelem : mesh.l_numpts;
373 
374  size_t nelem = mesh.l_numelem;
375  const SF::vector<mesh_int_t> & tags = mesh.tag;
376 
377  // check whether all specified tags exist in the element list
378  check_tags_in_elems(tags, regspec, mesh.name.c_str(), reglist);
379 
380  regionIDs.assign(rIDsize, 0);
382 
383  // we generate a map from tags to region IDs. This has many benefits, mainly we can check
384  // whether a tag is assigned to multiple regions and simplify our regionIDs filling loop
385  for (size_t reg=1; reg < regspec.size(); reg++) {
386  int err = 0;
387  for (int k=0; k < regspec[reg].nsubregs; k++) {
388  int curtag = regspec[reg].subregtags[k];
389  if(tag_to_reg.count(curtag)) err++;
390 
391  if(get_global(err, MPI_SUM))
392  log_msg(0,4,0, "%s warning: Tag idx %d is assigned to multiple regions!\n"
393  "Its final assignment will be to the highest assigned region ID!",
394  __func__, curtag);
395 
396  tag_to_reg[curtag] = reg;
397  }
398  }
399 
400  SF::vector<int> mx_tag(regionIDs.size(), -1);
401 
402  // cycle through the element list
403  for(size_t i=0; i<nelem; i++) {
404  const mesh_int_t & cur_tag = tags[i];
405  // check if current element tag has a custom region ID
406  if (tag_to_reg.count(cur_tag))
407  {
408  int reg = tag_to_reg[cur_tag];
409 
410  if (mask_elem) regionIDs[i] = reg;
411  else {
412  eview.set_elem(i);
413 
414  for (int j=0; j < eview.num_nodes(); j++) {
415  mesh_int_t n = eview.node(j);
416  if (cur_tag > mx_tag[n]) mx_tag[n] = cur_tag;
417  }
418  }
419  }
420  }
421 
422  if(mask_elem) return;
423 
424  for(size_t i=0; i < regionIDs.size(); i++)
425  if(mx_tag[i] >= 0) regionIDs[i] = tag_to_reg[mx_tag[i]];
426 
427  // for POINT masks we need to do a nodal -> algebraic map
428  mesh.pl.reduce(regionIDs, "max");
429 
430  const SF::vector<mesh_int_t> & alg_nod = mesh.pl.algebraic_nodes();
432 
433  SF::vector<int> alg_reg(alg_nod.size());
434  SF::vector<int> gids (alg_nod.size());
435 
436  for(size_t i=0; i<alg_nod.size(); i++) {
437  alg_reg[i] = regionIDs[alg_nod[i]];
438  gids[i] = nbr[alg_nod[i]];
439  }
440 
441  regionIDs = alg_reg;
442 
443  if(param_globals::dump_imp_region) {
444  set_dir(OUTPUT);
445  write_data_ascii(PETSC_COMM_WORLD, gids, regionIDs, std::string(reglist)+".dat");
446  }
447 }
448 
451 double Ionics::timer_val(const int timer_id)
452 {
453  double val = std::nan("NaN");
454  return val;
455 }
456 
459 std::string Ionics::timer_unit(const int timer_id)
460 {
461  std::string s_unit;
462  return s_unit;
463 }
464 
466 {
467  if (impdata[limpet::Iion] != NULL)
468  impdata[limpet::Iion][n] = 0;
469 
470  pIF.for_each([&](limpet::IonIfBase& imp) {
471  update_ts(&imp.get_tstp());
472  imp.compute(n, n + 1, impdata);
473  });
474 }
475 
487 void* find_SV_in_IMP(limpet::MULTI_IF* miif, const int idx, const char *IMP, const char *SV,
488  int* offset, int* sz, int* plugin_idx)
489 {
490  *plugin_idx = -1;
491 
492  if(strcmp(IMP, miif->iontypes[idx].get().get_name().c_str()) == 0) {
493  return (void*) miif->iontypes[idx].get().get_sv_offset(SV, offset, sz);
494  }
495  else {
496  for(int k=0; k<miif->numplugs[idx]; k++)
497  if(strcmp(IMP, miif->plugtypes[idx][k].get().get_name().c_str()) == 0) {
498  *plugin_idx = k;
499  return (void*) miif->plugtypes[idx][k].get().get_sv_offset(SV, offset, sz);
500  }
501  }
502 
503  return NULL;
504 }
505 
506 
516 void alloc_gvec_data(const int nGVcs, const int nRegs,
517  GVecs *prmGVecs, gvec_data &glob_vecs)
518 {
519  glob_vecs.nRegs = nRegs;
520 
521  if (nGVcs) {
522  glob_vecs.vecs.resize(nGVcs);
523  glob_vecs.plugin_idx.resize(nGVcs);
524 
525  for (size_t i = 0; i < glob_vecs.vecs.size(); i++) {
526  sv_data &gvec = glob_vecs.vecs[i];
527 
528  gvec.name = dupstr(prmGVecs[i].name);
529  gvec.units = dupstr(prmGVecs[i].units);
530  gvec.bogus = prmGVecs[i].bogus;
531 
532  gvec.imps = (char**) calloc(nRegs, sizeof(char *));
533  gvec.svNames = (char**) calloc(nRegs, sizeof(char *));
534  gvec.svSizes = (int*) calloc(nRegs, sizeof(int));
535  gvec.svOff = (int*) calloc(nRegs, sizeof(int));
536  gvec.getsv = (void**) calloc(nRegs, sizeof(limpet::SVgetfcn));
537 
538  for (int j = 0; j < nRegs; j++) {
539  if (strlen(prmGVecs[i].imp)) gvec.imps[j] = dupstr(prmGVecs[i].imp);
540 #ifdef WITH_PURK
541  else if (j < param_globals::num_imp_regions)
542  gvec.imps[j] = dupstr(param_globals::imp_region[j].im);
543  else
544  gvec.imps[j] = dupstr(param_globals::PurkIon[j - param_globals::num_imp_regions].im);
545 #else
546  else if (j < param_globals::num_imp_regions)
547  gvec.imps[j] = dupstr(param_globals::imp_region[j].im);
548 #endif
549  gvec.svNames[j] = dupstr(prmGVecs[i].ID[j]);
550  }
551  }
552  }
553 }
554 
568  igb_output_manager & output_manager)
569 {
570  GVs.inclPS = false;
571  // int num_purk_regions = GVs->inclPS ? purk->ion.N_IIF : 0;
572  int num_purk_regions = 0;
573  int nRegs = param_globals::num_imp_regions + num_purk_regions;
574 
575  alloc_gvec_data(param_globals::num_gvecs, nRegs, param_globals::gvec, GVs);
576 
577 #ifdef WITH_PURK
578  if (GVs->inclPS) sample_PS_ionSVs(purk);
579 #endif
580 
581  for (unsigned int i = 0; i < GVs.vecs.size(); i++) {
582  sv_data & gv = GVs.vecs[i];
583  int noSV = 0;
584 
585  for (int j = 0; j < miif->N_IIF; j++) {
586  gv.getsv[j] = find_SV_in_IMP(miif, j, gv.imps[j], gv.svNames[j], gv.svOff + j, gv.svSizes + j, &GVs.plugin_idx[i]);
587 
588  if (gv.getsv[j] == NULL) {
589  log_msg(NULL, 3, ECHO, "\tWarning: SV(%s) not found in region %d\n", gv.svNames[j], j);
590  noSV++;
591  }
592  }
593 
594  SF::init_vector(&gv.ordered, &tmpl);
595  output_manager.register_output(gv.ordered, intra_elec_msh, 1, gv.name, gv.units);
596 
597 #ifdef WITH_PURK
598  // same procedure for Purkinje
599  if (GVs->inclPS) {
600  IF_PURK_PROC(purk) {
601  MULTI_IF* pmiif = &purk->ion;
602  for (int j = miif->N_IIF; j < nRegs; j++) {
603  gv.getsv[j] = find_SV_in_IMP(pmiif, j - miif->N_IIF, gv.imps[j],
604  gv.svNames[j], gv.svOff + j, gv.svSizes + j);
605 
606  if (gv.getsv[j] == NULL) {
607  LOG_MSG(NULL, 3, ECHO, "\tWarning: state variable \"%s\" not found in region %d\n", gv.svNames[j], j);
608  noSV++;
609  }
610  }
611  RVector_dup(purk->vm_pt, &gv.orderedPS_PS);
612  MYO_COMM(purk);
613  }
614  RVector_dup(purk->vm_pt_over, &gv.orderedPS);
615  initialize_grid_output(grid, NULL, tmo, intra_elec_msh, GRID_WRITE, 0., 1., gv.units, gv.orderedPS,
616  1, gv.GVcName, -purk->npt, param_globals::output_level);
617  }
618 #endif
619 
620  if (noSV == nRegs) {
621  log_msg(NULL, 5, ECHO, "\tError: no state variables found for global vector %d\n", i);
622  log_msg(NULL, 5, ECHO, "Run bench --imp=YourModel --imp-info to get a list of all parameters.\\n");
623  exit(1);
624  }
625  }
626 }
627 
639 {
640  for(size_t i=0; i<gvecs.vecs.size(); i++ ) {
641  sv_data & gv = gvecs.vecs[i];
642  int plugin_idx = gvecs.plugin_idx[i];
643 
644  // set to the defalt value
645  gv.ordered->set(gv.bogus);
646  SF_int start, stop;
647  gv.ordered->get_ownership_range(start, stop);
648 
649  for( int n = 0; n<miif->N_IIF; n++ ) {
650  if( !gv.getsv[n] ) continue;
651 
652  SF::vector<SF_real> data (miif->N_Nodes[n]);
653  SF::vector<SF_int> indices(miif->N_Nodes[n]);
654 
655  limpet::IonIfBase* base_im = miif->IIF[n];
656  limpet::IonIfBase* imp = plugin_idx > -1 ? base_im->plugins()[plugin_idx] : base_im;
657 
658  for( int j=0; j<miif->N_Nodes[n]; j++ ) {
659  indices[j] = miif->NodeLists[n][j] + start;
660  data[j] = ((limpet::SVgetfcn)(gv.getsv[n]))( *miif->IIF[n], j, gv.svOff[n]);
661  }
662 
663  bool add = false;
664  gv.ordered->set(indices, data, add);
665  }
666  }
667 }
668 
669 
670 } // namespace opencarp
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
virtual void get_ownership_range(T &start, T &stop) const =0
virtual void set(const vector< T > &idx, const vector< S > &vals, const bool additive=false, const bool local=false)=0
Comfort class. Provides getter functions to access the mesh member variables more comfortably.
Definition: SF_fem_utils.h:704
const T & node(short nidx) const
Access the connectivity information.
Definition: SF_fem_utils.h:793
void set_elem(size_t eidx)
Set the view to a new element.
Definition: SF_fem_utils.h:731
T num_nodes() const
Getter function for the number of nodes.
Definition: SF_fem_utils.h:761
overlapping_layout< T > pl
nodal parallel layout
Definition: SF_container.h:429
size_t l_numelem
local number of elements
Definition: SF_container.h:399
std::string name
the mesh name
Definition: SF_container.h:407
size_t l_numpts
local number of points
Definition: SF_container.h:401
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
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 resize(size_t n)
Resize a vector.
Definition: SF_vector.h:209
const T * end() const
Pointer to the vector's end.
Definition: SF_vector.h:128
void assign(InputIterator s, InputIterator e)
Assign a memory range.
Definition: SF_vector.h:161
const T * begin() const
Pointer to the vector's start.
Definition: SF_vector.h:116
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
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
std::vector< IonIfBase * > & plugins()
Returns a vector containing the plugins of this IMP.
Definition: ION_IF.cc:184
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
int numNode
local number of nodes
Definition: MULTI_ION_IF.h:210
bool extUpdateVm
flag indicating update function for Vm
Definition: MULTI_ION_IF.h:208
std::vector< IonIfBase * > IIF
array of IIF's
Definition: MULTI_ION_IF.h:202
opencarp::sf_vec * gdata[NUM_IMP_DATA_TYPES]
data used by all IMPs
Definition: MULTI_ION_IF.h:216
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< Target > targets
target for each region
Definition: MULTI_ION_IF.h:213
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
void initialize_currents(double, int)
float restore_state(const char *, opencarp::mesh_t gid, bool)
int N_IIF
how many different IIF's
Definition: MULTI_ION_IF.h:211
void compute_ionic_current(bool flag_send=1, bool flag_receive=1)
GPU kernel to emulate the add_scaled call made to adjust the Vm values when the update to Vm is not m...
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
IIF_Mask_t * IIFmask
region for each node
Definition: MULTI_ION_IF.h:214
int adjust_MIIF_variables(const char *variable, const SF::vector< SF_int > &indices, const SF::vector< SF_real > &values)
std::string name
name for MIIF region
Definition: MULTI_ION_IF.h:198
FILE_SPEC logger
The logger of the physic, each physic should have one.
Definition: physics_types.h:64
const char * name
The name of the physic, each physic should have one.
Definition: physics_types.h:62
void output_step()
Definition: ionics.cc:57
mesh_t ion_domain
Definition: ionics.h:67
generic_timing_stats comp_stats
Definition: ionics.h:69
limpet::MULTI_IF * miif
Definition: ionics.h:66
void compute_step()
Definition: ionics.cc:35
void initialize()
Definition: ionics.cc:60
void destroy()
Definition: ionics.cc:52
std::string timer_unit(const int timer_id)
figure out units of a signal linked to a given timer
Definition: ionics.cc:459
double timer_val(const int timer_id)
figure out current value of a signal linked to a given timer
Definition: ionics.cc:451
void register_output(sf_vec *inp_data, const mesh_t inp_meshid, const int dpn, const char *name, const char *units, const SF::vector< mesh_int_t > *idx=NULL, bool elem_data=false)
Register a data vector for output.
Definition: sim_utils.cc:1461
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 and LIMPET wrappers.
void write_data_ascii(const MPI_Comm comm, const vector< T > &idx, const vector< S > &data, std::string file, short dpn=1)
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
@ 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
void * find_SV_in_IMP(limpet::MULTI_IF *miif, const int idx, const char *IMP, const char *SV, int *offset, int *sz, int *plugin_idx)
Definition: ionics.cc:487
@ 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
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
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
void initialize_sv_dumps(limpet::MULTI_IF *pmiif, IMPregion *reg, int id, double t, double dump_dt)
Definition: ionics.cc:285
void init_sv_gvec(gvec_data &GVs, limpet::MULTI_IF *miif, sf_vec &tmpl, igb_output_manager &output_manager)
Definition: ionics.cc:567
void assemble_sv_gvec(gvec_data &gvecs, limpet::MULTI_IF *miif)
Definition: ionics.cc:638
char * dupstr(const char *old_str)
Definition: basics.cc:44
bool check_tags_in_elems(const SF::vector< int > &tags, SF::vector< RegionSpecs > &regspec, const char *gridname, const char *reglist)
Check whether the tags in the region spec struct matches with an array of tags.
Definition: ionics.cc:329
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 alloc_gvec_data(const int nGVcs, const int nRegs, GVecs *prmGVecs, gvec_data &glob_vecs)
Definition: ionics.cc:516
void get_time(double &tm)
Definition: basics.h:436
SF::abstract_vector< SF_int, SF_real > sf_vec
Definition: sf_interface.h:50
void compute_IIF(limpet::IonIfBase &pIF, limpet::GlobalData_t **impdata, int n)
Definition: ionics.cc:465
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
void log_stats(double tm, bool cflg)
Definition: timers.cc:93
void init_logger(const char *filename)
Definition: timers.cc:77
int calls
# calls for this interval, this is incremented externally
Definition: timers.h:70
double tot_time
total time, this is incremented externally
Definition: timers.h:72
SF::vector< int > plugin_idx
if we use a plugin, its index in the plugins list of the IMP will be stored here, else -1.
Definition: ionics.h:149
SF::vector< sv_data > vecs
store sv dump indices for global vectors
Definition: ionics.h:148
unsigned int nRegs
number of imp regions
Definition: ionics.h:146
bool inclPS
include PS if exists
Definition: ionics.h:147
float bogus
value indicating sv not in region
Definition: ionics.h:142
void ** getsv
functions to retrieve sv
Definition: ionics.h:139
char * name
Name of global composite sv vector.
Definition: ionics.h:132
int * svOff
sv size in bytes
Definition: ionics.h:137
char ** svNames
sv names of components forming global vector
Definition: ionics.h:134
int * svSizes
sv size in bytes
Definition: ionics.h:136
char ** imps
Name of imp to which sv belongs.
Definition: ionics.h:133
char * units
units of sv
Definition: ionics.h:141
sf_vec * ordered
vector in which to place ordered data
Definition: ionics.h:140