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  SF::vector<RegionSpecs> rs(param_globals::num_imp_regions);
90  for (size_t i = 0; i < rs.size(); i++ ) {
91  rs[i].nsubregs = param_globals::imp_region[i].num_IDs;
92  rs[i].subregtags = param_globals::imp_region[i].ID;
93  }
94 
95  SF::vector<int> reg_mask;
96  region_mask(ion_domain, rs, reg_mask, false, "imp_region");
97 
98  int purkfLen = 1;
99  tstart = setup_MIIF(loc_size, param_globals::num_imp_regions, param_globals::imp_region,
100  reg_mask.data(), param_globals::start_statef, param_globals::num_adjustments,
101  param_globals::adjustment, param_globals::dt, purkfLen > 0);
102 
103  miif->extUpdateVm = !param_globals::operator_splitting;
104 
105  // if we start at a non-zero time (i.e. we have restored a state), we notify the
106  // timer manager
107  if(tstart > 0.) {
108  log_msg(logger, 0, 0, "Changing simulation start time to %.2lf", tstart);
109  user_globals::tm_manager->setup(param_globals::dt, tstart, param_globals::tend);
111  }
112 
113  set_dir(INPUT);
114 
115  this->initialize_time += timing(t2, t1);
116 }
117 
118 double Ionics::setup_MIIF(int nnodes, int nreg, IMPregion* impreg, int* mask,
119  const char *start_fn, int numadjust, IMPVariableAdjustment *adjust,
120  double time_step, bool close)
121 {
122  double tstart = 0;
123 
124  miif->N_IIF = nreg;
125  miif->numNode = nnodes;
126  miif->iontypes = {};
127  miif->numplugs = (int*)calloc( miif->N_IIF, sizeof(int));
128  miif->plugtypes = std::vector<limpet::IonTypeList>(miif->N_IIF);
129  miif->targets = std::vector<limpet::Target>(miif->N_IIF, limpet::Target::AUTO);
130 
131  log_msg(logger,0,ECHO, "\nSetting up ionic models and plugins\n" \
132  "-----------------------------------\n\n" \
133  "Assigning IMPS to tagged regions:" );
134 
135  for (int i=0;i<miif->N_IIF;i++) {
136  auto pT = limpet::get_ion_type(std::string(impreg[i].im));
137  if (pT != NULL)
138  {
139  miif->iontypes.push_back(*pT);
140  log_msg(logger, 0, ECHO|NONL, "\tIonic model: %s to tag region(s)", impreg[i].im);
141 
142  if(impreg[i].num_IDs > 0) {
143  for(int j = 0; j < impreg[i].num_IDs; j++)
144  log_msg(logger,0,ECHO|NONL, " [%d],", impreg[i].ID[j]);
145  log_msg(logger,0,ECHO,"\b.");
146  }
147  else {
148  log_msg(logger,0,ECHO, " [0] (implicitely)");
149  }
150  }
151  else {
152  log_msg(NULL,5,ECHO, "Illegal IM specified: %s\n", impreg[i].im );
153  log_msg(NULL,5,ECHO, "Run bench --list-imps for a list of all available models.\n" );
154  EXIT(1);
155  }
156  if (limpet::get_plug_flag( impreg[i].plugins, &miif->numplugs[i], miif->plugtypes[i]))
157  {
158  if(impreg[i].plugins[0] != '\0') {
159  log_msg(logger,0, ECHO|NONL, "\tPlug-in(s) : %s to tag region(s)", impreg[i].plugins);
160 
161  for(int j = 0; j < impreg[i].num_IDs; j++)
162  log_msg(logger,0,ECHO|NONL, " [%d],", impreg[i].ID[j]);
163  log_msg(logger,0,ECHO,"\b.");
164  }
165  }
166  else {
167  log_msg(NULL,5,ECHO,"Illegal plugin specified: %s\n", impreg[i].plugins);
168  log_msg(NULL,5,ECHO, "Run bench --list-imps for a list of all available plugins.\n" );
169  EXIT(1);
170  }
171  }
172 
174 
175  // The mask is a nodal vector with the region IDs of each node.
176  // It is already reduced during the region_mask call to guarantee unique values
177  // for overlapping interface nodes
178  if (mask) {
179  for (int i=0; i<miif->numNode; i++)
180  miif->IIFmask[i] = (limpet::IIF_Mask_t) mask[i];
181  }
182 
184 
185  for(int i=0; i<miif->N_IIF; i++)
186  if(!miif->IIF[i]->cgeom().SVratio)
187  miif->IIF[i]->cgeom().SVratio = param_globals::imp_region[i].cellSurfVolRatio;
188 
189  for (int i=0;i<miif->N_IIF;i++) {
190  // the IMP tuning does not handle spaces well, thus we remove them here
191  remove_char(impreg[i].im_param, strlen(impreg[i].im_param), ' ');
192  miif->IIF[i]->tune(impreg[i].im_param, impreg[i].plugins, impreg[i].plug_param);
193  }
194 
195  set_dir(INPUT);
196  miif->initialize_currents(time_step, param_globals::ode_fac);
197 
198  // overriding initial values goes here
199  // read in single cell state vector and spread it out over the entire region
200  set_dir(INPUT);
201  for (int i=0;i<miif->N_IIF;i++) {
202  if (impreg[i].im_sv_init && strlen(impreg[i].im_sv_init) > 0)
203  if (read_sv(miif, i, impreg[i].im_sv_init)) {
204  log_msg(NULL, 5, ECHO|FLUSH, "State vector initialization failed for %s.\n", impreg[i].name);
205  EXIT(-1);
206  }
207  }
208 
209  if( !start_fn || strlen(start_fn)>0 )
210  tstart = (double) miif->restore_state(start_fn, ion_domain, close);
211 
212  for (int i=0; i<numadjust; i++)
213  {
214  set_dir(INPUT);
215 
216  SF::vector<int> indices;
217  SF::vector<double> values;
218  bool restrict_to_algebraic = true;
219 
220  sf_mesh & imesh = get_mesh(ion_domain);
221  std::map<std::string,std::string> metadata;
222  read_metadata(adjust[i].file, metadata, PETSC_COMM_WORLD);
223 
224  SF::SF_nbr nbr = SF::NBR_REF;
225  if(metadata.count("grid") && metadata["grid"].compare("intra") == 0) {
226  nbr = SF::NBR_SUBMESH;
227  }
228  read_indices_with_data(indices, values, adjust[i].file, imesh, nbr, restrict_to_algebraic, 1, PETSC_COMM_WORLD);
229 
230  int rank = get_rank();
231 
232  for(size_t gi = 0; gi < indices.size(); gi++)
233  indices[gi] = SF::local_nodal_to_local_petsc(imesh, rank, indices[gi]);
234 
235  // debug, output parameters on global intracellular vector
236  if(adjust[i].dump) {
237  sf_vec* adjPars;
238  SF::init_vector(&adjPars, imesh, 1, sf_vec::algebraic);
239  adjPars->set(indices, values);
240 
241  set_dir(OUTPUT);
242  char fname[2085];
243  snprintf(fname, sizeof fname, "adj_%s_perm.dat", adjust[i].variable);
244  adjPars->write_ascii(fname, false);
245 
246  // get the scattering to the canonical permutation
248  if(sc == NULL) {
249  log_msg(0,3,0, "%s warning: PETSC_TO_CANONICAL permutation needed registering!", __func__);
251  }
252 
253  (*sc)(*adjPars, true);
254  snprintf(fname, sizeof fname, "adj_%s_canonical.dat", adjust[i].variable);
255  adjPars->write_ascii(fname, false);
256  }
257 
258  int nc = miif->adjust_MIIF_variables(adjust[i].variable, indices, values);
259  log_msg(logger, 0, 0, "Adjusted %d values for %s", nc, adjust[i].variable);
260  }
261 
262  set_dir(OUTPUT);
263 
264  for (int i=0;i<miif->N_IIF;i++)
265  initialize_sv_dumps(miif, impreg+i, i, tstart, param_globals::spacedt);
266 
267  return tstart;
268 }
269 
279 void initialize_sv_dumps(limpet::MULTI_IF *pmiif, IMPregion* reg, int id, double t, double dump_dt)
280 {
281  char svs[1024], plgs[1024], plgsvs[1024], fname[1024];
282 
283  strcpy(svs, reg->im_sv_dumps ? reg->im_sv_dumps : "");
284  strcpy(plgs, reg->plugins ? reg->plugins : "");
285  strcpy(plgsvs, reg->plug_sv_dumps ? reg->plug_sv_dumps : "");
286 
287  if( !(strlen(svs)+strlen(plgsvs) ) )
288  return;
289 
290  /* The string passed to the "reg_name" argument (#4) of the sv_dump_add
291  * function is supposed to be "region name". It's only purpose is to
292  * provide the base name for the SV dump file, eg: Purkinje.Ca_i.bin.
293  * Thus, we pass: [vofile].[reg name], Otherwise, dumping SVs in
294  * batched runs would be extremely tedious.
295  */
296  strcpy(fname, param_globals::vofile); // [vofile].igb.gz
297 
298  if( !reg->name ) {
299  log_msg(NULL, 5, ECHO, "%s: a region name must be specified\n", __func__ );
300  exit(0);
301  }
302 
303  // We want to convert vofile.igb or vofile.igb.gz to vofile.regname
304  char* ext_start = NULL;
305  ext_start = strstr(fname, "igb.gz");
306  if(ext_start == NULL) ext_start = strstr(fname, "igb");
307  if(ext_start == NULL) ext_start = fname + strlen(fname);
308  strcpy(ext_start, reg->name);
309 
310  pmiif->sv_dump_add_by_name_list(id, reg->im, fname, svs, plgs, plgsvs, t, dump_dt);
311 }
312 
324  const char* gridname, const char* reglist)
325 {
326  bool AllTagsExist = true;
327 
329  tagset.insert(tags.begin(), tags.end());
330 
331  // cycle through all user-specified regions
332  for (size_t reg=0; reg<regspec.size(); reg++)
333  // cycle through all tags which belong to the region
334  for (int k=0; k<regspec[reg].nsubregs; k++) {
335  // check whether this tag exists in element list
336  int n = 0;
337  if(tagset.count(regspec[reg].subregtags[k])) n++;
338  // globalize n
339  int N = get_global(n, MPI_SUM);
340  if (N==0) {
341  log_msg(NULL, 3, ECHO,
342  "Region tag %d in %s[%d] not found in element list for %s grid.\n",
343  regspec[reg].subregtags[k], reglist, reg, gridname);
344  AllTagsExist = false;
345  }
346  }
347 
348  if (!AllTagsExist) {
349  log_msg(NULL, 4, ECHO,"Missing element tags in element file!\n"
350  "Check region ID specs in input files!\n");
351  }
352 
353  return AllTagsExist;
354 }
355 
356 void region_mask(mesh_t meshspec, SF::vector<RegionSpecs> & regspec,
357  SF::vector<int> & regionIDs, bool mask_elem, const char* reglist)
358 {
359  if(regspec.size() == 1) return;
360 
361  sf_mesh & mesh = get_mesh(meshspec);
363 
364  // initialize the list with the default regionID, 0
365  size_t rIDsize = mask_elem ? mesh.l_numelem : mesh.l_numpts;
366 
367  size_t nelem = mesh.l_numelem;
368  const SF::vector<mesh_int_t> & tags = mesh.tag;
369 
370  // check whether all specified tags exist in the element list
371  check_tags_in_elems(tags, regspec, mesh.name.c_str(), reglist);
372 
373  regionIDs.assign(rIDsize, 0);
375 
376  // we generate a map from tags to region IDs. This has many benefits, mainly we can check
377  // whether a tag is assigned to multiple regions and simplify our regionIDs filling loop
378  for (size_t reg=1; reg < regspec.size(); reg++) {
379  int err = 0;
380  for (int k=0; k < regspec[reg].nsubregs; k++) {
381  int curtag = regspec[reg].subregtags[k];
382  if(tag_to_reg.count(curtag)) err++;
383 
384  if(get_global(err, MPI_SUM))
385  log_msg(0,4,0, "%s warning: Tag idx %d is assigned to multiple regions!\n"
386  "Its final assignment will be to the highest assigned region ID!",
387  __func__, curtag);
388 
389  tag_to_reg[curtag] = reg;
390  }
391  }
392 
393  SF::vector<int> mx_tag(regionIDs.size(), -1);
394 
395  // cycle through the element list
396  for(size_t i=0; i<nelem; i++) {
397  const mesh_int_t & cur_tag = tags[i];
398  // check if current element tag has a custom region ID
399  if (tag_to_reg.count(cur_tag))
400  {
401  int reg = tag_to_reg[cur_tag];
402 
403  if (mask_elem) regionIDs[i] = reg;
404  else {
405  eview.set_elem(i);
406 
407  for (int j=0; j < eview.num_nodes(); j++) {
408  mesh_int_t n = eview.node(j);
409  if (cur_tag > mx_tag[n]) mx_tag[n] = cur_tag;
410  }
411  }
412  }
413  }
414 
415  if(mask_elem) return;
416 
417  for(size_t i=0; i < regionIDs.size(); i++)
418  if(mx_tag[i] >= 0) regionIDs[i] = tag_to_reg[mx_tag[i]];
419 
420  // for POINT masks we need to do a nodal -> algebraic map
421  mesh.pl.reduce(regionIDs, "max");
422 
423  const SF::vector<mesh_int_t> & alg_nod = mesh.pl.algebraic_nodes();
425 
426  SF::vector<int> alg_reg(alg_nod.size());
427  SF::vector<int> gids (alg_nod.size());
428 
429  for(size_t i=0; i<alg_nod.size(); i++) {
430  alg_reg[i] = regionIDs[alg_nod[i]];
431  gids[i] = nbr[alg_nod[i]];
432  }
433 
434  regionIDs = alg_reg;
435 
436  if(param_globals::dump_imp_region) {
437  set_dir(OUTPUT);
438  write_data_ascii(PETSC_COMM_WORLD, gids, regionIDs, std::string(reglist)+".dat");
439  }
440 }
441 
444 double Ionics::timer_val(const int timer_id)
445 {
446  double val = std::nan("NaN");
447  return val;
448 }
449 
452 std::string Ionics::timer_unit(const int timer_id)
453 {
454  std::string s_unit;
455  return s_unit;
456 }
457 
459 {
460  if (impdata[limpet::Iion] != NULL)
461  impdata[limpet::Iion][n] = 0;
462 
463  pIF.for_each([&](limpet::IonIfBase& imp) {
464  update_ts(&imp.get_tstp());
465  imp.compute(n, n + 1, impdata);
466  });
467 }
468 
480 void* find_SV_in_IMP(limpet::MULTI_IF* miif, const int idx, const char *IMP, const char *SV,
481  int* offset, int* sz)
482 {
483  if(strcmp(IMP, miif->iontypes[idx].get().get_name().c_str()) == 0) {
484  return (void*) miif->iontypes[idx].get().get_sv_offset(SV, offset, sz);
485  }
486  else {
487  for( int k=0; k<miif->numplugs[idx]; k++ )
488  if(strcmp(IMP, miif->plugtypes[idx][k].get().get_name().c_str()) == 0) {
489  return (void*) miif->plugtypes[idx][k].get().get_sv_offset(SV, offset, sz);
490  }
491  }
492 
493  return NULL;
494 }
495 
496 
506 void alloc_gvec_data(const int nGVcs, const int nRegs,
507  GVecs *prmGVecs, gvec_data &glob_vecs)
508 {
509  glob_vecs.nRegs = nRegs;
510 
511  if (nGVcs) {
512  glob_vecs.vecs.resize(nGVcs);
513 
514  for (size_t i = 0; i < glob_vecs.vecs.size(); i++) {
515  sv_data &gvec = glob_vecs.vecs[i];
516 
517  gvec.name = dupstr(prmGVecs[i].name);
518  gvec.units = dupstr(prmGVecs[i].units);
519  gvec.bogus = prmGVecs[i].bogus;
520 
521  gvec.imps = (char**) calloc(nRegs, sizeof(char *));
522  gvec.svNames = (char**) calloc(nRegs, sizeof(char *));
523  gvec.svSizes = (int*) calloc(nRegs, sizeof(int));
524  gvec.svOff = (int*) calloc(nRegs, sizeof(int));
525  gvec.getsv = (void**) calloc(nRegs, sizeof(limpet::SVgetfcn));
526 
527  for (int j = 0; j < nRegs; j++) {
528  if (strlen(prmGVecs[i].imp)) gvec.imps[j] = dupstr(prmGVecs[i].imp);
529 #ifdef WITH_PURK
530  else if (j < param_globals::num_imp_regions)
531  gvec.imps[j] = dupstr(param_globals::imp_region[j].im);
532  else
533  gvec.imps[j] = dupstr(param_globals::PurkIon[j - param_globals::num_imp_regions].im);
534 #else
535  else if (j < param_globals::num_imp_regions)
536  gvec.imps[j] = dupstr(param_globals::imp_region[j].im);
537 #endif
538  gvec.svNames[j] = dupstr(prmGVecs[i].ID[j]);
539  }
540  }
541  }
542 }
543 
557  igb_output_manager & output_manager)
558 {
559  GVs.inclPS = false;
560  // int num_purk_regions = GVs->inclPS ? purk->ion.N_IIF : 0;
561  int num_purk_regions = 0;
562  int nRegs = param_globals::num_imp_regions + num_purk_regions;
563 
564  alloc_gvec_data(param_globals::num_gvecs, nRegs, param_globals::gvec, GVs);
565 
566 #ifdef WITH_PURK
567  if (GVs->inclPS) sample_PS_ionSVs(purk);
568 #endif
569 
570  for (unsigned int i = 0; i < GVs.vecs.size(); i++) {
571  sv_data & gv = GVs.vecs[i];
572  int noSV = 0;
573 
574  for (int j = 0; j < miif->N_IIF; j++) {
575  gv.getsv[j] = find_SV_in_IMP(miif, j, gv.imps[j], gv.svNames[j],
576  gv.svOff + j, gv.svSizes + j);
577 
578  if (gv.getsv[j] == NULL) {
579  log_msg(NULL, 3, ECHO, "\tWarning: SV(%s) not found in region %d\n", gv.svNames[j], j);
580  noSV++;
581  }
582  }
583 
584  SF::init_vector(&gv.ordered, &tmpl);
585  output_manager.register_output(gv.ordered, intra_elec_msh, 1, gv.name, gv.units);
586 
587 #ifdef WITH_PURK
588  // same procedure for Purkinje
589  if (GVs->inclPS) {
590  IF_PURK_PROC(purk) {
591  MULTI_IF* pmiif = &purk->ion;
592  for (int j = miif->N_IIF; j < nRegs; j++) {
593  gv.getsv[j] = find_SV_in_IMP(pmiif, j - miif->N_IIF, gv.imps[j],
594  gv.svNames[j], gv.svOff + j, gv.svSizes + j);
595 
596  if (gv.getsv[j] == NULL) {
597  LOG_MSG(NULL, 3, ECHO, "\tWarning: state variable \"%s\" not found in region %d\n", gv.svNames[j], j);
598  noSV++;
599  }
600  }
601  RVector_dup(purk->vm_pt, &gv.orderedPS_PS);
602  MYO_COMM(purk);
603  }
604  RVector_dup(purk->vm_pt_over, &gv.orderedPS);
605  initialize_grid_output(grid, NULL, tmo, intra_elec_msh, GRID_WRITE, 0., 1., gv.units, gv.orderedPS,
606  1, gv.GVcName, -purk->npt, param_globals::output_level);
607  }
608 #endif
609 
610  if (noSV == nRegs) {
611  log_msg(NULL, 5, ECHO, "\tError: no state variables found for global vector %d\n", i);
612  log_msg(NULL, 5, ECHO, "Run bench --imp=YourModel --imp-info to get a list of all parameters.\\n");
613  exit(1);
614  }
615  }
616 }
617 
629 {
630  for(size_t i=0; i<gvecs.vecs.size(); i++ ) {
631  sv_data & gv = gvecs.vecs[i];
632 
633  // set to the defalt value
634  gv.ordered->set(gv.bogus);
635  SF_int start, stop;
636  gv.ordered->get_ownership_range(start, stop);
637 
638  for( int n = 0; n<miif->N_IIF; n++ ) {
639  if( !gv.getsv[n] ) continue;
640 
641  SF::vector<SF_real> data (miif->N_Nodes[n]);
642  SF::vector<SF_int> indices(miif->N_Nodes[n]);
643 
644  for( int j=0; j<miif->N_Nodes[n]; j++ ) {
645  indices[j] = miif->NodeLists[n][j] + start;
646  data[j] = ((limpet::SVgetfcn)(gv.getsv[n]))( *miif->IIF[n], j, gv.svOff[n]);
647  }
648 
649  bool add = false;
650  gv.ordered->set(indices, data, add);
651  }
652 
653 #ifdef WITH_PURK
654  // include purkinje in sv dump, if exists
655  if(gvecs->inclPS) {
656  RVector_set( gv.orderedPS, gv.bogus ); // set this to the default value
657 
658  if(IS_PURK_PROC(purk))
659  {
660  MULTI_IF *pmiif = &purk->ion;
661  Real *data = new Real[purk->gvec_cab.nitems];
662 
663  int ci=0;
664  for( int n = 0; n<pmiif->N_IIF; n++ ) {
665  int gvidx = n+miif->N_IIF;
666  if( !gv.getsv[gvidx] ) continue;
667  for( int j=0; j<purk->gvec_ion[n].nitems; j++ )
668  data[ci++] = ((SVgetfcn)(gv.getsv[gvidx]))( pmiif->IIF+n, ((int*)(purk->gvec_ion[n].data))[j],
669  gv.svOff[gvidx]);
670  }
671  RVector_setvals(gv.orderedPS, purk->gvec_cab.nitems, (int*)purk->gvec_cab.data, data, true);
672  delete [] data;
673  }
674  RVector_sync( gv.orderedPS );
675  }
676 #endif
677  }
678 }
679 
680 
681 } // 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 set(const vector< T > &idx, const vector< S > &vals, const bool additive=false)=0
virtual void get_ownership_range(T &start, T &stop) const =0
Comfort class. Provides getter functions to access the mesh member variables more comfortably.
Definition: SF_fem_utils.h:705
const T & node(short nidx) const
Access the connectivity information.
Definition: SF_fem_utils.h:795
void set_elem(size_t eidx)
Set the view to a new element.
Definition: SF_fem_utils.h:733
T num_nodes() const
Getter function for the number of nodes.
Definition: SF_fem_utils.h:763
overlapping_layout< T > pl
nodal parallel layout
Definition: SF_container.h:419
size_t l_numelem
local number of elements
Definition: SF_container.h:388
std::string name
the mesh name
Definition: SF_container.h:396
size_t l_numpts
local number of points
Definition: SF_container.h:390
vector< T > & get_numbering(SF_nbr nbr_type)
Get the vector defining a certain numbering.
Definition: SF_container.h:454
vector< T > tag
element tag
Definition: SF_container.h:406
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
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:138
ts & get_tstp()
Gets the time stepper.
Definition: ION_IF.cc:202
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:411
void compute(int start, int end, GlobalData_t **data)
Perform ionic model computation for 1 time step.
Definition: ION_IF.cc:264
int numNode
local number of nodes
Definition: MULTI_ION_IF.h:209
bool extUpdateVm
flag indicating update function for Vm
Definition: MULTI_ION_IF.h:207
std::vector< IonIfBase * > IIF
array of IIF's
Definition: MULTI_ION_IF.h:201
opencarp::sf_vec * gdata[NUM_IMP_DATA_TYPES]
data used by all IMPs
Definition: MULTI_ION_IF.h:215
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:208
std::vector< Target > targets
target for each region
Definition: MULTI_ION_IF.h:212
std::vector< IonTypeList > plugtypes
plugins types for each region
Definition: MULTI_ION_IF.h:214
IonTypeList iontypes
type for each region
Definition: MULTI_ION_IF.h:211
void initialize_currents(double, int)
int adjust_MIIF_variables(const char *variable, const SF::vector< int > &indices, const SF::vector< double > &values)
float restore_state(const char *, opencarp::mesh_t gid, bool)
int N_IIF
how many different IIF's
Definition: MULTI_ION_IF.h:210
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:199
int ** NodeLists
local partitioned node lists for each IMP stored
Definition: MULTI_ION_IF.h:200
IIF_Mask_t * IIFmask
region for each node
Definition: MULTI_ION_IF.h:213
std::string name
name for MIIF region
Definition: MULTI_ION_IF.h:197
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:452
double timer_val(const int timer_id)
figure out current value of a signal linked to a given timer
Definition: ionics.cc:444
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:1409
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)
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:122
SF_nbr
Enumeration encoding the different supported numberings.
Definition: SF_container.h:189
@ NBR_PETSC
PETSc numbering of nodes.
Definition: SF_container.h:192
@ NBR_REF
The nodal numbering of the reference mesh (the one stored on HD).
Definition: SF_container.h:190
@ NBR_SUBMESH
Submesh nodal numbering: The globally ascending sorted reference indices are reindexed.
Definition: SF_container.h:191
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:460
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:44
timer_manager * tm_manager
a manager for the various physics timers
Definition: main.cc:52
@ iotm_console
Definition: timer_utils.h:44
SF::abstract_vector< mesh_int_t, double > sf_vec
Definition: sf_interface.h:49
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:47
int set_dir(IO_t dest)
Definition: sim_utils.cc:452
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
void * find_SV_in_IMP(limpet::MULTI_IF *miif, const int idx, const char *IMP, const char *SV, int *offset, int *sz)
Definition: ionics.cc:480
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:356
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:850
@ OUTPUT
Definition: sim_utils.h:52
void initialize_sv_dumps(limpet::MULTI_IF *pmiif, IMPregion *reg, int id, double t, double dump_dt)
Definition: ionics.cc:279
void init_sv_gvec(gvec_data &GVs, limpet::MULTI_IF *miif, sf_vec &tmpl, igb_output_manager &output_manager)
Definition: ionics.cc:556
void assemble_sv_gvec(gvec_data &gvecs, limpet::MULTI_IF *miif)
Definition: ionics.cc:628
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:323
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:60
@ intra_elec_msh
Definition: sf_interface.h:61
void alloc_gvec_data(const int nGVcs, const int nRegs, GVecs *prmGVecs, gvec_data &glob_vecs)
Definition: ionics.cc:506
void get_time(double &tm)
Definition: basics.h:436
void compute_IIF(limpet::IonIfBase &pIF, limpet::GlobalData_t **impdata, int n)
Definition: ionics.cc:458
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:76
void log_stats(double tm, bool cflg)
Definition: timers.cc:91
void init_logger(const char *filename)
Definition: timers.cc:75
int calls
# calls for this interval, this is incremented externally
Definition: timers.h:35
double tot_time
total time, this is incremented externally
Definition: timers.h:37
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